QwGeant4
QweakSimVDC_WirePlaneSD.cc
Go to the documentation of this file.
1 //=============================================================================
2 //
3 // ---------------------------
4 // | Doxygen File Information |
5 // ---------------------------
6 //
7 /**
8 
9  \file QweakSimVDC_WirePlaneSD.cc
10 
11  $Revision: 1.5 $
12  $Date: 2006/05/05 21:51:07 $
13 
14  \author Klaus Hans Grimm
15 
16 */
17 //=============================================================================
18 
19 //=============================================================================
20 // -----------------------
21 // | CVS File Information |
22 // -----------------------
23 //
24 // Last Update: $Author: grimm $
25 // Update Date: $Date: 2006/05/05 21:51:07 $
26 // CVS/RCS Revision: $Revision: 1.5 $
27 // Status: $State: Exp $
28 //
29 // ===================================
30 // CVS Revision Log at end of file !!
31 // ===================================
32 //
33 //============================================================================
34 
35 
36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
37 
39 
40 // user includes
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
48 :G4VSensitiveDetector(name)
49 {
50 // G4cout << G4endl << "###### Calling QweakSimVDC_WirePlaneSD::QweakSimVDC_WirePlaneSD() " << G4endl << G4endl;
51 
52  collectionName.insert("VDCWirePlaneCollection"); //collectionName[0]
53  HCID = -1;
54 
55 // G4cout << G4endl << "###### Leaving QweakSimVDC_WirePlaneSD::QweakSimVDC_WirePlaneSD() " << G4endl << G4endl;
56 }
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60 {
61  //G4cout << G4endl << "###### Calling QweakSimVDC_WirePlaneSD::~QweakSimVDC_WirePlaneSD() " << G4endl << G4endl;
62  //delete hitsCollection;
63  //G4cout << G4endl << "###### Leaving QweakSimVDC_WirePlaneSD::~QweakSimVDC_WirePlaneSD() " << G4endl << G4endl;
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 void QweakSimVDC_WirePlaneSD::Initialize(G4HCofThisEvent* HCE)
68 {
69  //G4cout << G4endl << "###### Calling QweakSimVDC_WirePlaneSD::Initialize() " << G4endl << G4endl;
70 
71  hitsCollection = new QweakSimVDC_WirePlane_HitsCollection(SensitiveDetectorName,collectionName[0]);
72 
73  // if(HCID<0)
74  HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
75 
76  HCE->AddHitsCollection( HCID, hitsCollection );
77 
78  //G4cout << G4endl << "###### Leaving QweakSimVDC_WirePlaneSD::Initialize() " << G4endl << G4endl;
79 }
80 
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82 G4bool QweakSimVDC_WirePlaneSD::ProcessHits(G4Step* aStep,G4TouchableHistory* /*ROhist*/)
83 {
84  //G4cout << G4endl << "###### Calling QweakSimVDC_WirePlaneSD::ProcessHits() " << G4endl << G4endl;
85 
86 
87  G4double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
88  // reject non-charged particle hits
89  if (fabs(charge)<0.1)
90  return false;
91 
92  G4String particlename = aStep->GetTrack()->GetDefinition()->GetParticleName();
93  G4int particletype = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
94 
95 
96 
97  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
98  G4StepPoint* postStepPoint = aStep->GetPreStepPoint();
99  G4TouchableHistory* theTouchable = (G4TouchableHistory*)(preStepPoint->GetTouchable());
100  G4VPhysicalVolume* physVol = theTouchable->GetVolume();
101 
102  // we are only interested if the particle is crossing a boundary (glass->vaccuum)
103  if (postStepPoint->GetStepStatus() != fGeomBoundary) {
104 // G4cout << "=======================================================================" << G4endl;
105 // G4cout << ">>>> QweakSimVDC_WirePlaneSD: We are NOT crossing a boundary <<<<<<" << G4endl;
106 // G4cout << ">>>> Aborting QweakSimVDC_WirePlaneSD::ProcessHits <<<<<<" << G4endl;
107 // G4cout << "=======================================================================" << G4endl;
108 
109  return false;
110  }
111 
112 // G4cout << "================== We are in volume :" << physVol->GetName() << G4endl << G4endl;
113 
114 
115  if( (strcmp(physVol->GetName(),"VDC_UPlane_Physical")==0) || (strcmp(physVol->GetName(),"VDC_VPlane_Physical")==0) ) {
116 
117 
118 // if( (strcmp(physVol->GetName(),"VDC_UPlane_Physical")==0) )
119 // {
120 // G4cout << "=============================================================" << G4endl;
121 // G4cout << ">>>>>>> Particle crossing : VDC_UPlane_Physical <<<<<<<<<" << G4endl;
122 // G4cout << "=============================================================" << G4endl;
123 // }
124 
125 
126 // if( (strcmp(physVol->GetName(),"VDC_VPlane_Physical")==0) )
127 // {
128 // G4cout << "=============================================================" << G4endl;
129 // G4cout << ">>>>>>> Particle crossing : VDC_VPlane_Physical <<<<<<<<<" << G4endl;
130 // G4cout << "=============================================================" << G4endl;
131 // }
132 
133 
134  //G4double currentKineticEnergy = aStep->GetTrack()->GetKineticEnergy();
135  G4double currentKineticEnergy = preStepPoint->GetKineticEnergy();
136 
137  //G4double currentTotalEnergy = aStep->GetTrack()->GetTotalEnergy();
138  G4double currentTotalEnergy = preStepPoint->GetTotalEnergy();
139 
140 
141  //G4ThreeVector currentMomentumDirection = aStep->GetTrack()->GetMomentumDirection();
142  G4ThreeVector currentMomentumDirection = preStepPoint->GetMomentumDirection();
143 
144 // G4int trackID = aStep->GetTrack()->GetTrackID();
145 // G4cout << "====> Track ID : " << trackID << G4endl;
146 
147 // G4int parentID = aStep->GetTrack()->GetParentID();
148 // G4cout << "====> Parent ID : " << parentID << G4endl;
149 
150 
151 
152  //G4int MotherCopyNo = theTouchable->GetVolume(1)->GetCopyNo(); // Several MotherVolumes
153  G4int WirePlaneCopyNo = theTouchable->GetVolume()->GetCopyNo(); // but only one WirePlane per MV
154  //G4int MotherReplicaNo = theTouchable->GetReplicaNumber(1); // Several MotherVolumes
155 
156  //G4int MotherCopyNo2 = theTouchable->GetVolume(2)->GetCopyNo(); // Several MotherVolumes
157  G4int MotherReplicaNo2 = theTouchable->GetReplicaNumber(2); // Several MotherVolumes
158 
159 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane MV CopyNumber :" << MotherCopyNo << G4endl;
160 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane CopyNumber :" << WirePlaneCopyNo << G4endl;
161 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane ReplicaNumber :" << MotherReplicaNo << G4endl;
162 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane MV CopyNumber2 :" << MotherCopyNo2 << G4endl;
163 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane ReplicaNumber2 :" << MotherReplicaNo2 << G4endl;
164 
165 
166  G4ThreeVector worldPos = preStepPoint->GetPosition();
167 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldPos.X [cm]:" << worldPos.x()/cm << G4endl;
168 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldPos.Y [cm]:" << worldPos.y()/cm << G4endl;
169 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldPos.Z [cm]:" << worldPos.z()/cm << G4endl;
170 
171  G4ThreeVector localPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
172 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalPos.X [cm]:" << localPos.x()/cm << G4endl;
173 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalPos.Y [cm]:" << localPos.y()/cm << G4endl;
174 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalPos.Z [cm]:" << localPos.z()/cm << G4endl;
175 
176  G4ThreeVector worldMomentum = preStepPoint->GetMomentum();
177 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldMomentum.X [GeV]:" << worldMomentum.x()/GeV << G4endl;
178 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldMomentum.Y [GeV]:" << worldMomentum.y()/GeV << G4endl;
179 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldMomentum.Z [GeV]:" << worldMomentum.z()/GeV << G4endl;
180 
181  G4ThreeVector localMomentum = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldMomentum);
182 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalMomentum.X [GeV]:" << localMomentum.x()/GeV << G4endl;
183 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalMomentum.Y [GeV]:" << localMomentum.y()/GeV << G4endl;
184 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalMomentum.Z [GeV]:" << localMomentum.z()/GeV << G4endl;
185 
186 
187 
188  // Vertex: where this track was generated
189  G4ThreeVector originVertexPosition = aStep->GetTrack()->GetVertexPosition();
190  G4ThreeVector originVertexMomentumDirection = aStep->GetTrack()->GetVertexMomentumDirection();
191  G4double originVertexKineticEnergy = aStep->GetTrack()->GetVertexKineticEnergy();
192 
193  // G4cout << "====> originVertexPosition_X [cm] : " << originVertexPosition.x()/cm << G4endl;
194 // G4cout << "====> originVertexPosition_Y [cm] : " << originVertexPosition.y()/cm << G4endl;
195 // G4cout << "====> originVertexPosition_Z [cm] : " << originVertexPosition.z()/cm << G4endl;
196 
197 // G4cout << "====> originVertexMomentumDirection_X : " << originVertexMomentumDirection.x() << G4endl;
198 // G4cout << "====> originVertexMomentumDirection_Y : " << originVertexMomentumDirection.y() << G4endl;
199 // G4cout << "====> originVertexMomentumDirection_Z : " << originVertexMomentumDirection.z() << G4endl;
200 
201 // G4cout << "====> originVertexKineticEnergy [Mev] : " << originVertexKineticEnergy/MeV << G4endl;
202 
203  // get User Track Info
204  //QweakSimTrackInformation* info = (QweakSimTrackInformation*) (aStep->GetTrack()->GetUserInformation());
205  //info->Print();
206 
207 // G4double primaryQ2 = info->GetPrimaryQ2();
208 // G4double crossSection = info->GetCrossSection();
209 // G4double crossSectionWeight = info->GetCrossSectionWeight();
210 // G4int primaryEventNumber = info->GetPrimaryEventNumber();
211 
212 // G4cout << " @@@@@@@@@@@@@@@@@@@@@@@ Original Track Q2 : " << primaryQ2 << G4endl;
213 // G4cout << " @@@@@@@@@@@@@@@@@@@@@@@ Original Track Event Number : " << primaryEventNumber << G4endl;
214 // G4cout << " @@@@@@@@@@@@@@@@@@@@@@@ Original CrossSection Weight (ub) : " << crossSectionWeight << G4endl;
215 
216 
217 
218 
219  // QweakSimVDCHit* aHit = new QweakSimVDCHit(MotherCopyNo); // there is only one plane per motherVolume
220  QweakSimVDC_WirePlaneHit* aHit = new QweakSimVDC_WirePlaneHit(); // there is only one plane per motherVolume
221 
222  aHit->StoreVDCID(MotherReplicaNo2); // 0: Front VDC , 1: Back VDC
223 
224  // old VDC version
225  //aHit->StoreWirePlaneID(MotherReplicaNo); // 0: Front WirePlane , 1: Back WirePlane
226 
227  // new VDC version with foils etc
228  aHit->StoreWirePlaneID(WirePlaneCopyNo); // 0: Front WirePlane , 1: Back WirePlane
229 
230  aHit->StoreWorldPosition(worldPos);
231  aHit->StoreLocalPosition(localPos);
232  aHit->StoreGlobalTime(preStepPoint->GetGlobalTime());
233 
234  aHit->StoreOriginVertexPosition(originVertexPosition);
235  aHit->StoreOriginVertexMomentumDirection(originVertexMomentumDirection);
236  aHit->StoreOriginVertexKineticEnergy(originVertexKineticEnergy);
237 
238 // aHit->StorePrimaryQ2(primaryQ2);
239 // aHit->StoreCrossSection(crossSection);
240 // aHit->StoreCrossSectionWeight(crossSectionWeight);
241 // aHit->StorePrimaryEventNumber(primaryEventNumber);
242 
243  aHit->StoreWorldMomentum(worldMomentum);
244  aHit->StoreLocalMomentum(localMomentum);
245 
246  aHit->StoreParticleName(particlename);
247  aHit->StoreParticleType(particletype);
248 
249  //aHit->StoreMomentumDirection(currentMomentumDirection);
250  aHit->StoreKineticEnergy(currentKineticEnergy);
251  aHit->StoreTotalEnergy(currentTotalEnergy);
252 
253  // now store all the hit information into the hit collection
254  hitsCollection->insert(aHit);
255 
256  } // end of if(strcmp(physVol->GetName(),"VDC_WirePlane_Physical")==0)
257 
258 
259  //G4cout << G4endl << "###### Leaving QweakSimVDC_WirePlaneSD::ProcessHits() " << G4endl << G4endl;
260 
261  return true;
262 }
263 
264 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
265 //void QweakSimVDC_WirePlaneSD::EndOfEvent(G4HCofThisEvent* HCE)
266 void QweakSimVDC_WirePlaneSD::EndOfEvent(G4HCofThisEvent* )
267 {
268  //G4cout << G4endl << "###### Calling QweakSimVDC_WirePlaneSD::EndOfEvent() " << G4endl << G4endl;
269 
270 // G4int NbHits = hitsCollection->entries();
271 // G4cout << "\n-------->Hits Collection: in this event they are " << NbHits
272 // << " hits in the Wire Planes: " << G4endl;
273 // for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Print();
274 
275  //G4cout << G4endl << "###### Leaving QweakSimVDC_WirePlaneSD::EndOfEvent() " << G4endl << G4endl;
276 
277 }
278 
279 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
280 
281 //=======================================================
282 // -----------------------
283 // | CVS File Information |
284 // -----------------------
285 //
286 // $Revisions$
287 // $Log: QweakSimVDC_WirePlaneSD.cc,v $
288 // Revision 1.5 2006/05/05 21:51:07 grimm
289 // Records now the kinetic and total energy of track/particle
290 //
291 // Revision 1.4 2006/01/06 21:32:50 grimm
292 // Storing ParticleName and ParticleType
293 //
294 // Revision 1.3 2006/01/13 23:22:00 grimm
295 // Due to the updated the VDC geometry code, the variable for the wireplane index has changed.
296 // Commented the target origin cut: now all charged particles hits,regardless of origin, will be saved.
297 //
298 // Revision 1.2 2005/12/27 19:19:55 grimm
299 // - Redesign of Doxygen header containing CVS info like revision and date
300 // - Added CVS revision log at the end of file
301 //
302 //
void StoreOriginVertexKineticEnergy(G4double ekin)
void StoreKineticEnergy(G4double ekin)
QweakSimVDC_WirePlane_HitsCollection * hitsCollection
void StoreLocalMomentum(G4ThreeVector lpxyz)
void EndOfEvent(G4HCofThisEvent *HCE)
void StoreLocalPosition(G4ThreeVector lxyz)
void StoreWirePlaneID(G4int wireplane_ID)
G4THitsCollection< QweakSimVDC_WirePlaneHit > QweakSimVDC_WirePlane_HitsCollection
void StoreWorldMomentum(G4ThreeVector gpxyz)
void StoreWorldPosition(G4ThreeVector gxyz)
void StoreOriginVertexMomentumDirection(G4ThreeVector pxyz)
Handling of a U-WirePlane and/or V-WirePlane Hit of the VDC.
void StoreOriginVertexPosition(G4ThreeVector xyz)
void Initialize(G4HCofThisEvent *HCE)
G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist)