QwGeant4
QweakSimTrackingAction.cc
Go to the documentation of this file.
1 //=============================================================================
2 //
3 // ---------------------------
4 // | Doxygen File Information |
5 // ---------------------------
6 //
7 /**
8 
9  \file QweakSimTrackingAction.cc
10 
11  $Revision: 1.5 $
12  $Date: 2006/01/18 20:27:32 $
13 
14  \author Klaus Hans Grimm
15 
16 */
17 //=============================================================================
18 
19 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
20 
22 
23 // geant4 includes
24 #include "G4TrackingManager.hh"
25 #include "G4OpticalPhoton.hh"
26 
27 // user includes
31 #include "QweakSimTrajectory.hh"
32 
33 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
34 
36  : myUserInfo(myUI)
37 {
38  G4cout << G4endl << "###### Calling QweakSimTrackingAction::QweakSimTrackingAction() " << G4endl << G4endl;
39 
41 
42  G4cout << G4endl << "###### Leaving QweakSimTrackingAction::QweakSimTrackingAction() " << G4endl << G4endl;
43 }
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47 {
48  G4cout << G4endl << "###### Calling QweakSimTrackingAction::~QweakSimTrackingAction() " << G4endl << G4endl;
49 
51 
52  G4cout << G4endl << "###### Leaving QweakSimTrackingAction::~QweakSimTrackingAction() " << G4endl << G4endl;
53 }
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 {
58 
59  //G4cout << G4endl << "###### Calling QweakSimTrackingAction::PreUserTrackingAction() " << G4endl << G4endl;
60 
61  // store only trajectory of primary
62  // BUT allow tracking of secondaries with an origin of z>568*cm
63 
64 // if( (aTrack->GetParentID()==0) || (aTrack->GetVertexPosition().z() > 568.0*cm) )
65 // {
66 // fpTrackingManager->SetStoreTrajectory(true);
67 // //fpTrackingManager->SetTrajectory(new G4Trajectory(aTrack));
68 // }
69 // else
70 // {
71 // fpTrackingManager->SetStoreTrajectory(false);
72 // }
73 
74  //-----------------------------------------------------------
75  //
76 
77  /*
78  if (TrackingFlag==1)
79  {
80  // store only trajectory of primary particles (== electrons from event generator) for tracking
81  //if( (aTrack->GetParentID()==0) || (aTrack->GetVertexPosition().z() > 568.0*cm) )
82 
83  if( (aTrack->GetParentID()==0) || (aTrack->GetDefinition() == G4OpticalPhoton::OpticalPhotonDefinition()) )
84  {
85  fpTrackingManager->SetStoreTrajectory(true);
86 
87  // store track
88  //fpTrackingManager->SetTrajectory(new G4Trajectory(aTrack));
89  fpTrackingManager->SetTrajectory(new QweakSimTrajectory(aTrack));
90  }
91  else
92  {
93  // G4cout << "===> Not a primary event, Trajectory not stored !!!" << G4endl;
94  fpTrackingManager->SetStoreTrajectory(false);
95  }
96  }
97  else
98  {
99  // store all trajectories (primary+secondaries)
100  fpTrackingManager->SetStoreTrajectory(true);
101  fpTrackingManager->SetTrajectory(new QweakSimTrajectory(aTrack));
102 
103  }*/
104 
105 
106 
107  if (TrackingFlag==0) //track primary electron only
108  {
109  // store only trajectory of primary particles (== electrons from event generator) for tracking
110  //if( (aTrack->GetParentID()==0) || (aTrack->GetVertexPosition().z() > 568.0*cm) )
111 
112  if( aTrack->GetParentID()==0 )
113  {
114  fpTrackingManager->SetStoreTrajectory(true);
115 
116  // store track
117  fpTrackingManager->SetTrajectory(new QweakSimTrajectory(aTrack));
118  }
119  else
120  {
121  //fpTrackingManager->SetStoreTrajectory(false);
122  //aTrack->SetTrackStatus(fStopAndKill);
123  fpTrackingManager->EventAborted();
124  return;
125  }
126  }
127  else if (TrackingFlag==1) //track primary electron and optical photon
128  {
129  if( (aTrack->GetParentID()==0) || (aTrack->GetDefinition() == G4OpticalPhoton::OpticalPhotonDefinition()) )
130  {
131  fpTrackingManager->SetStoreTrajectory(true);
132 
133  // store track
134  fpTrackingManager->SetTrajectory(new QweakSimTrajectory(aTrack));
135  }
136  else
137  {
138  //fpTrackingManager->SetStoreTrajectory(false);
139  //aTrack->SetTrackStatus(fStopAndKill);
140  fpTrackingManager->EventAborted();
141  return;
142  }
143 
144  }
145 
146  else if (TrackingFlag==2) //track primary electron and secondaries except optical photon
147  {
148  if( aTrack->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition() )
149  {
150  fpTrackingManager->SetStoreTrajectory(true);
151 
152  // store track
153  fpTrackingManager->SetTrajectory(new QweakSimTrajectory(aTrack));
154  }
155  else
156  {
157  //fpTrackingManager->SetStoreTrajectory(false);
158  //aTrack->SetTrackStatus(fStopAndKill);
159  fpTrackingManager->EventAborted();
160  return;
161  }
162 
163  }
164 
165  else // if (TrackingFlag==3) //track all particles
166  {
167  // store all trajectories (primary+secondaries)
168  fpTrackingManager->SetStoreTrajectory(true);
169  fpTrackingManager->SetTrajectory(new QweakSimTrajectory(aTrack));
170 
171  }
172  //
173  //-----------------------------------------------------------
174 
175 
176  // Check if the track already has track info
177  // in case of a primary event w/o user track info:
178  // create/store user track info
179  if( (aTrack->GetParentID()==0) && (aTrack->GetUserInformation()==0) )
180  {
181  // G4cout << G4endl << "###### Creating new Track User Information for primary" << G4endl << G4endl;
182 
183  // create user track info using current primary track pointer as an input
185 
186  // fill user track info with data stored in myUserInfo (class QweakUserInformation)
187 // anInfo->StorePrimaryQ2(myUserInfo->GetPrimaryQ2());
188 // anInfo->StoreCrossSection(myUserInfo->GetCrossSection());
189 // anInfo->StoreCrossSectionWeight(myUserInfo->GetCrossSectionWeight());
190 // anInfo->StorePrimaryEventNumber(myUserInfo->GetPrimaryEventNumber());
191 
192 /*
193  anInfo->StorePrimaryTrackID(myUserInfo->GetTrackID());
194 
195  int particleType = myUserInfo->GetPDGcode();
196  if (particleType==3)
197  anInfo->StorePrimaryParticle(G4Electron::ElectronDefinition());
198 
199  G4double theX = myUserInfo->GetOriginVertexPositionX();
200  G4double theY = myUserInfo->GetOriginVertexPositionY();
201  G4double theZ = myUserInfo->GetOriginVertexPositionZ();
202  anInfo->StorePrimaryPosition(G4ThreeVector(theX, theY, theZ));
203 
204  G4double mdX = myUserInfo->GetOriginVertexMomentumDirectionX();
205  G4double mdY = myUserInfo->GetOriginVertexMomentumDirectionY();
206  G4double mdZ = myUserInfo->GetOriginVertexMomentumDirectionZ();
207  G4double mdNormal = sqrt(mdX*mdX+mdY*mdY+mdZ*mdZ);
208  G4double Ek = myUserInfo->GetOriginVertexKineticEnergy();
209  G4double px = Ek*mdX/mdNormal;
210  G4double py = Ek*mdY/mdNormal;
211  G4double pz = Ek*mdZ/mdNormal;
212  anInfo->StorePrimaryMomentum(G4ThreeVector(px,py,pz));
213  anInfo->StorePrimaryEnergy(Ek);
214  anInfo->StorePrimaryKineticEnergy(Ek);
215  anInfo->StorePrimaryTime(myUserInfo->GetGlobalTime());
216 */
217 
218  // set the source track info (which is here identical to the primary track info)
219  anInfo->SetSourceTrackInformation(aTrack);
220  //anInfo->Print();
221 
222  // access the track pointer
223  G4Track* theTrack = (G4Track*) aTrack;
224 
225  // attach/expand track with user track info
226  theTrack->SetUserInformation(anInfo);
227 
228  //delete anInfo;
229  }
230 
231 
232 
233  // G4cout << G4endl << "###### Leaving QweakSimTrackingAction::PreUserTrackingAction() " << G4endl << G4endl;
234 
235 }
236 
237 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
238 void QweakSimTrackingAction::PostUserTrackingAction(const G4Track* /*aTrack*/)
239 {
240  //G4cout << G4endl << "###### Calling QweakSimTrackingAction::PostUserTrackingAction() " << G4endl << G4endl;
241 
242  //-------------------------------------------------------------------------
243  // attach user track info of the primary particle to the seconary particles
244  //-------------------------------------------------------------------------
245 
246  // get pointer to secondaries
247  G4TrackVector* secondaries = fpTrackingManager->GimmeSecondaries();
248 
249 
250  if(secondaries)
251  {
252 
253  // how many secondaries do we have ?
254  G4int nSeco = secondaries->size();
255 
256  //G4cout << "###### --> How many secondaries: " << nSeco << G4endl;
257 
258  if(nSeco>0)
259  {
260  //G4cout << "###### --> Within secondary loop " << G4endl;
261 
262  // get pointer to primary track info
263 // QweakSimTrackInformation* info = (QweakSimTrackInformation*)(aTrack->GetUserInformation());
264 
265 // // NEW: add the track information of the source track to the secondaries
266 // // (don't mix up source track info with unique primary track info)
267 // // nedded for figuring out if a primary electron or secondary particles
268 // // has created a cerenkov photon
269 
270 
271 
272 // //G4cout << "###### --> We got info = (QweakSimTrackInformation*)(aTrack->GetUserInformation() " << G4endl;
273 
274 
275 // for(G4int i=0;i<nSeco;i++)
276 // {
277 // // create secondary track info pointer for each secondary
278 // QweakSimTrackInformation* infoNew = new QweakSimTrackInformation(info);
279 
280 
281 // // copy primary track info to the current secondary
282 // (*secondaries)[i]->SetUserInformation(infoNew);
283 
284 // // check for cerenkov photon generation
285 // if ( (*secondaries)[i]->GetDefinition()==G4OpticalPhoton::OpticalPhotonDefinition() &&
286 // (*secondaries)[i]->GetCreatorProcess()->GetProcessName()=="Cerenkov" )
287 // {
288 // G4cout << " Cerenkov Photon generated here" << G4endl;
289 
290 // info->SetSourceTrackInformation(aTrack);
291 
292 // }
293 
294 // // tried it: G4 crashes if I delete infoNew
295 // //delete infoNew;
296 
297 // } // end for
298  } // end if (nSeco>0)
299  } // end if (secondaries)
300 
301  //G4cout << G4endl << "###### Leaving QweakSimTrackingAction::PostUserTrackingAction() " << G4endl << G4endl;
302 
303 }
304 
305 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
306 
Stores the information about the various tracks.
QweakSimTrackingAction(QweakSimUserInformation *)
void PreUserTrackingAction(const G4Track *aTrack)
void SetSourceTrackInformation(const G4Track *aTrack)
Scans the input file for /TrackingAction/xyz commands.
QweakSimTrackingActionMessenger * pTrackingActionMessenger
Class with additional track information like Q2.
void PostUserTrackingAction(const G4Track *aTrack)