ACTARSim
ACTAR TPC Simulation Reference Guide
ActarSimSilRingSD.cc
Go to the documentation of this file.
1 // - AUTHOR: Hector Alvarez-Pol 04/2008
2 /******************************************************************
3  * Copyright (C) 2005-2016, Hector Alvarez-Pol *
4  * All rights reserved. *
5  * *
6  * License according to GNU LESSER GPL (see lgpl-3.0.txt). *
7  * For the list of contributors see CREDITS. *
8  ******************************************************************/
9 //////////////////////////////////////////////////////////////////
10 /// \class ActarSimSilRingSD
11 /// SD for the Silicon Ring
12 /////////////////////////////////////////////////////////////////
13 
14 #include "ActarSimSilRingSD.hh"
15 
16 #include "G4HCofThisEvent.hh"
17 #include "G4Step.hh"
18 #include "G4ThreeVector.hh"
19 #include "G4SDManager.hh"
20 #include "G4ios.hh"
21 #include "G4UnitsTable.hh"
22 
23 #include "G4VProcess.hh"
24 #include "G4VPhysicalVolume.hh"
25 #include "G4TouchableHistory.hh"
26 #include "G4VTouchable.hh"
27 
28 #include "G4PhysicalConstants.hh"
29 #include "G4SystemOfUnits.hh"
30 
31 //////////////////////////////////////////////////////////////////
32 /// Constructor, just naming the Hit collection
34  :G4VSensitiveDetector(name){
35  G4String HCname;
36  collectionName.insert(HCname="SilRingCollection");
37 }
38 
39 //////////////////////////////////////////////////////////////////
40 /// Destructor, nothing to do
42 }
43 
44 //////////////////////////////////////////////////////////////////
45 /// Initializing the ActarSimSilRingGeantHitsCollection object
46 /// Invoked automatically at the beggining of each event
47 void ActarSimSilRingSD::Initialize(G4HCofThisEvent* HCE){
49  (SensitiveDetectorName,collectionName[0]);
50  static G4int HCID = -1;
51  if(HCID<0)
52  HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
53 
54  HCE->AddHitsCollection( HCID, hitsCollection );
55 }
56 
57 //////////////////////////////////////////////////////////////////
58 /// Filling the ActarSimSilRingGeantHit information with the step info
59 /// Invoked by G4SteppingManager for each step
60 G4bool ActarSimSilRingSD::ProcessHits(G4Step* aStep,G4TouchableHistory*){
61  G4double edep = aStep->GetTotalEnergyDeposit();
62 
63  if(edep==0.) return false;
64 
66 
67  newHit->SetEdep(edep/MeV);
68  newHit->SetEBeforeSil(aStep->GetPreStepPoint()->GetKineticEnergy()/MeV);
69  newHit->SetEAfterSil(aStep->GetPostStepPoint()->GetKineticEnergy()/MeV);
70 
71  newHit->SetPos(aStep->GetPostStepPoint()->GetPosition());
72  //The elements should be taken of the PostStep, the real parameters of the
73  //interaction. The PreStep keeps record of where the particle comes from
74  newHit->SetPrePos(aStep->GetPreStepPoint()->GetPosition());
75 
76  //for the local translation, x = R-1 (x_m - T)
77  G4ThreeVector theLocalPos =
78  (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetObjectRotationValue().inverse()) *
79  (aStep->GetPostStepPoint()->GetPosition()-aStep->GetPostStepPoint()->GetPhysicalVolume()->GetObjectTranslation());
80  newHit->SetLocalPos(theLocalPos);
81 
82  theLocalPos = //just using the same variable...
83  (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetObjectRotationValue().inverse())*
84  (aStep->GetPreStepPoint()->GetPosition()-aStep->GetPreStepPoint()->GetPhysicalVolume()->GetObjectTranslation());
85  newHit->SetLocalPrePos(theLocalPos);
86 
87  newHit->SetDetName(aStep->GetTrack()->GetVolume()->GetName());
88  newHit->SetPreDetName(aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName());
89  newHit->SetPostDetName(aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName());
90  newHit->SetDetID(aStep->GetTrack()->GetVolume()->GetCopyNo());
91 
92  newHit->SetToF(aStep->GetPostStepPoint()->GetGlobalTime());
93 
94  newHit->SetTrackID(aStep->GetTrack()->GetTrackID());
95  newHit->SetParentID(aStep->GetTrack()->GetParentID());
96  newHit->SetParticleID(aStep->GetTrack()->GetDefinition()->GetPDGEncoding());
97  newHit->SetParticleCharge(aStep->GetTrack()->GetDefinition()->GetPDGCharge());
98  newHit->SetParticleMass(aStep->GetTrack()->GetDefinition()->GetPDGMass());
99 
100  hitsCollection->insert(newHit);
101 
102  // newHit cannot be deleted here !
103  // It should be deleted after the end of the event
104 
105  return true;
106 }
107 
108 //////////////////////////////////////////////////////////////////
109 /// Just prints and draws the event hits (class ActarSimSilRingGeantHit)
110 /// The recollection of the hits energy deposition in the plastic
111 /// is done in the ActarSimROOTAnalysis::EndOfEventAction()
112 void ActarSimSilRingSD::EndOfEvent(G4HCofThisEvent*){
113  G4int NbHits = hitsCollection->entries();
114  if (verboseLevel>0) {
115  G4cout << "Hits Collection: in this event they are " << NbHits
116  << " (GEANT-like) hits in the Sil: " << G4endl;
117  for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Print();
118  }
119  //for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Draw();
120 }
G4THitsCollection< ActarSimSilRingGeantHit > ActarSimSilRingGeantHitsCollection
void SetPrePos(G4ThreeVector xyz)
void Initialize(G4HCofThisEvent *)
void SetPos(G4ThreeVector xyz)
G4bool ProcessHits(G4Step *, G4TouchableHistory *)
void SetLocalPrePos(G4ThreeVector xyz)
void SetParticleCharge(G4double charge)
~ActarSimSilRingSD()
Destructor, nothing to do.
ActarSimSilRingGeantHitsCollection * hitsCollection
Geant step-like hits collect.
void SetLocalPos(G4ThreeVector xyz)
ActarSimSilRingSD(G4String)
Constructor, just naming the Hit collection.
void SetPreDetName(G4String Name)
void SetParticleMass(G4double mass)
void SetPostDetName(G4String Name)
void EndOfEvent(G4HCofThisEvent *)