QwGeant4
QweakSimTriggerScintillator_DetectorSD.cc
Go to the documentation of this file.
1 // QweakSimTriggerScintillator_DetectorSD.cc
2 // Klaus Hans Grimm
3 // 2005-12-27
4 
5 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
6 
8 
9 // user classes
13 
14 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
16 :G4VSensitiveDetector(name)
17 {
18  //G4cout << G4endl << "###### Calling QweakSimTriggerScintillator_DetectorSD::QweakSimTriggerScintillator_DetectorSD() " << G4endl << G4endl;
19 
20  collectionName.insert("TriggerScintillatorCollection");
22 
23  //G4cout << G4endl << "###### Leaving QweakSimTriggerScintillator_DetectorSD::QweakSimTriggerScintillator_DetectorSD() " << G4endl << G4endl;
24 }
25 
26 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
28 {
29  //G4cout << G4endl << "###### Calling QweakSimTriggerScintillator_DetectorSD::~QweakSimTriggerScintillator_DetectorSD() " << G4endl << G4endl;
30 
31  //delete TriggerScintillator_HitsCollection;
32 
33 }
34 
35 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
37 {
38  //G4cout << G4endl << "###### Calling QweakSimTriggerScintillator_DetectorSD::Initialize() " << G4endl << G4endl;
39 
41 
42  //if(TriggerScintillator_CollectionID<0)
43 
44  TriggerScintillator_CollectionID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
45 
47 
48  //G4cout << G4endl << "###### Leaving QweakSimTriggerScintillator_DetectorSD::Initialize() " << G4endl << G4endl;
49 }
50 
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52 //G4bool QweakSimTriggerScintillatorSD::ProcessHits(G4Step* aStep,G4TouchableHistory* /*ROhist*/)
53 G4bool QweakSimTriggerScintillator_DetectorSD::ProcessHits(G4Step* aStep, G4TouchableHistory* )
54 {
55 // G4cout << G4endl << "###### Calling QweakSimTriggerScintillator_DetectorSD::ProcessHits() " << G4endl << G4endl;
56 
57  G4Track *track = aStep->GetTrack();
58  G4String particlename = track->GetDefinition()->GetParticleName();
59  G4int particletype = track->GetDefinition()->GetPDGEncoding();
60 
61  G4ThreeVector worldPos;
62  G4ThreeVector localPos;
63  G4ThreeVector worldMomentum;
64 
65  G4TouchableHandle theTouchable = aStep->GetPreStepPoint()->GetTouchableHandle();
66 
67  if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ) {
68  // Entering Geometry
69  worldPos = aStep->GetPreStepPoint()->GetPosition();
70  localPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
71  worldMomentum = aStep->GetPreStepPoint()->GetMomentum();
72  } else {
73  return false;
74  }
75 
76  G4ThreeVector localMomentum = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldMomentum);
77 
78  // Vertex: where this track was generated
79  G4ThreeVector originVertexPosition = aStep->GetTrack()->GetVertexPosition();
80  G4ThreeVector originVertexMomentumDirection = aStep->GetTrack()->GetVertexMomentumDirection();
81  G4double originVertexKineticEnergy = aStep->GetTrack()->GetVertexKineticEnergy();
82 
83  G4double currentKineticEnergy = aStep->GetTrack()->GetKineticEnergy();
84  G4double currentTotalEnergy = aStep->GetTrack()->GetTotalEnergy();
85  G4ThreeVector currentMomentumDirection = aStep->GetTrack()->GetMomentumDirection();
86 
87  G4int trackID = aStep->GetTrack()->GetTrackID();
88 
89  G4int MotherCopyNo = theTouchable->GetVolume(1)->GetCopyNo(); // one Mother Volume
90 
92 
93  // Energy deposited in TS
94  G4double depositedEnergy = aStep->GetTotalEnergyDeposit();
95  aHit->StoreDepositedEnergy(depositedEnergy);
96 
97  aHit->StoreTrackID(trackID);
98 
99  aHit->StoreGlobalTime(aStep->GetPreStepPoint()->GetGlobalTime());
100 
101  aHit->StoreParticleName(particlename);
102  aHit->StoreParticleType(particletype);
103 
104  aHit->StoreWorldPosition(worldPos);
105  aHit->StoreLocalPosition(localPos);
106 
107  aHit->StoreWorldMomentum(worldMomentum);
108  aHit->StoreLocalMomentum(localMomentum);
109 
110  aHit->StoreOriginVertexPosition(originVertexPosition);
111  aHit->StoreOriginVertexKineticEnergy(originVertexKineticEnergy);
112  aHit->StoreOriginVertexTotalEnergy(originVertexKineticEnergy); /// \todo beware: total.neq.kinetic (testing only)
113  aHit->StoreOriginVertexMomentumDirection(originVertexMomentumDirection);
114 
115  aHit->StoreMomentumDirection(currentMomentumDirection);
116  aHit->StoreKineticEnergy(currentKineticEnergy);
117  aHit->StoreTotalEnergy(currentTotalEnergy);
118 
119  // check if it is first touch
120  if(!(aHit->GetLogVolume()))
121  {
122  // store translation and rotation matrix of the drift cell
123  // for the sake of drawing the hit
124  aHit->StoreLogVolume(theTouchable->GetVolume()->GetLogicalVolume());
125  G4AffineTransform aTrans = theTouchable->GetHistory()->GetTopTransform();
126  aTrans.Invert();
127  aHit->StoreCellRotation(aTrans.NetRotation());
128  aHit->StoreCellPosition(aTrans.NetTranslation());
129  }
130 
131 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% before inserting the hit" << G4endl;
132 
134 
135 // G4cout << G4endl << "###### Leaving QweakSimTriggerScintillatorSD::ProcessHits() " << G4endl << G4endl;
136 
137  return true;
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141 //void QweakSimTriggerScintillatorSD::EndOfEvent(G4HCofThisEvent* HCE)
143 {
144  //G4cout << G4endl << "###### Calling QweakSimTriggerScintillator_DetectorSD::EndOfEvent() " << G4endl << G4endl;
145 
146 // G4int NbDCHits = DC_hitsCollection->entries();
147 
148 // G4cout << "\n-------->Hits Collection: in this event they are " << NbDCHits
149 // << " hits in the Drift Cells : " << G4endl;
150 // for (G4int i=0;i<NbDCHits;i++) (*DC_hitsCollection)[i]->Print();
151 
152 
153 
154  //G4cout << G4endl << "###### Leaving QweakSimTriggerScintillator_DetectorSD::EndOfEvent() " << G4endl << G4endl;
155 
156 }
157 
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159 
160 
G4THitsCollection< QweakSimTriggerScintillator_DetectorHit > QweakSimTriggerScintillator_DetectorHitsCollection
G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist)
QweakSimTriggerScintillator_DetectorHitsCollection * TriggerScintillator_DetectorHitsCollection