SToGS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros
SToGS_BaseROOTEventsActions.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 
29 
30 #include "SToGS_G4_TrackerSD.hh"
31 #include "SToGS_G4_CopClusterSD.hh"
33 
34 #include "TTree.h"
35 
36 #include "G4OpticalPhoton.hh"
37 #include "G4VProcess.hh"
38 #include "G4Event.hh"
39 #if G4MULTITHREADED
40 #include "G4Threading.hh"
41 #endif
42 #ifdef G4MULTITHREADED
43 #include "G4AutoLock.hh"
44 namespace { G4Mutex buildMutex = G4MUTEX_INITIALIZER; }
45 #endif
46 
48  SToGS::BaseROOTTreeRun(),
49  fTree(tree),
50  fPrimaryEvent(),
51  fEvent()
52 {
53  fTree->Branch("Pr.", &fPrimaryEvent);
54  fTree->Branch("Ev.", &fEvent);
55 }
56 
57 void SToGS::BaseROOTEventsRun::RecordEvent(const G4Event* evt)
58 {
60  G4int nb_hits_tot = 0, nb_hits_trac = 0, nb_hits_calo = 0;
61 
62  // check there is something to write in files
63  if( colltrackerID < 0 && collcaloID < 0 )
64  return;
65 
66  fPrimaryEvent.Clear("");
67  fEvent.Clear("");
68 
69  G4HCofThisEvent *HCE = evt->GetHCofThisEvent();
70  if(HCE)
71  {
72  THC = (SToGS::TrackerHitsCollection *)(HCE->GetHC(colltrackerID));
73  if ( colltrackerID > -1 ) {
74  nb_hits_trac += THC->entries();
75  }
76  CHC = (SToGS::CopClusterHitsCollection*)(HCE->GetHC(collcaloID));
77  if ( collcaloID > -1 ) {
78  nb_hits_calo += CHC->entries();
79  }
80  nb_hits_tot = nb_hits_trac + nb_hits_calo;
81  }
82 
83  // write primary part
84  G4int istart = 0, iend = evt->GetNumberOfPrimaryVertex();
85  G4int K = 0; G4double H = 0.0;
86  for( G4int i = istart; i < iend ; i++ ){
87 
88  G4int jstart = 0, jend = evt->GetPrimaryVertex(i)->GetNumberOfParticle(); // to get the next vertex
89  for( G4int j = jstart; j < jend; j++ ) {
90 
91  SBRPHit *hit = fPrimaryEvent.AddHit();
92 
93  K++;
94  G4PrimaryParticle *prim = evt->GetPrimaryVertex(i)->GetPrimary(j); // get the next particle for vertex #i
95 
96  hit->fPDG = prim->GetPDGcode();
97  //
98  //hit->fE = (std::sqrt( prim->GetMomentum().mag2() + prim->GetMass() * prim->GetMass()) - prim->GetMass())/CLHEP::keV;
99  hit->fE = prim->GetKineticEnergy()/CLHEP::MeV;
100  //
101  hit->fX = evt->GetPrimaryVertex(i)->GetX0()/CLHEP::cm;
102  hit->fY = evt->GetPrimaryVertex(i)->GetY0()/CLHEP::cm;
103  hit->fZ = evt->GetPrimaryVertex(i)->GetZ0()/CLHEP::cm;
104  //
105  hit->fT = evt->GetPrimaryVertex(i)->GetT0();
106  //
107 
108  hit->fPX = prim->GetMomentumDirection().x();
109  hit->fPY = prim->GetMomentumDirection().y();
110  hit->fPZ = prim->GetMomentumDirection().z();
111 
112  hit->fFlag = prim->GetTrackID()-1;
113 
114  H += (std::sqrt( prim->GetMomentum().mag2() + prim->GetMass() * prim->GetMass()) - prim->GetMass())/CLHEP::MeV;
115  //H += (std::sqrt( prim->GetMomentum().mag2() + prim->GetMass() * prim->GetMass()) - prim->GetMass())/CLHEP::keV;
116  }
117  }
118  fPrimaryEvent.SetEMult(H,K);
119 
120  if ( nb_hits_trac ) {
121  int n_hit = nb_hits_trac;
122 
123  K = n_hit; H = 0.0;
124  for (int i = 0 ; i < n_hit; i++) {
125 
126  SBRHit *hit = fEvent.AddHit();
127 
128  if ( hit == 0x0 ) {
129  G4cout << " Cannot add more that " << i << " hits in the collection " << G4endl;
130  break;
131  }
132 
133  hit->fPDG = (*THC)[i]->GetPDGcode();
134 
135  hit->fE = (*THC)[i]->GetEdep()/CLHEP::MeV;
136  //hit->fE = (*THC)[i]->GetEdep()/CLHEP::keV;
137  hit->fX = (*THC)[i]->GetPos().x()/CLHEP::cm;
138  hit->fY = (*THC)[i]->GetPos().y()/CLHEP::cm;
139  hit->fZ = (*THC)[i]->GetPos().z()/CLHEP::cm;
140  hit->fT = (*THC)[i]->GetToF();
141 
142  hit->fFlag = (*THC)[i]->GetPrimaryID()-1 ;
143  hit->fUID = (*THC)[i]->GetDetID();
144 
145  H += (*THC)[i]->GetEdep()/CLHEP::MeV;
146  //H += (*THC)[i]->GetEdep()/CLHEP::keV;
147  }
148  fEvent.SetEMult(H,K);
149  }
150  if ( nb_hits_calo )
151  {
152  int n_hit = nb_hits_calo;
153 
154  K = n_hit; H = 0.0;
155  for (int i = 0 ; i < n_hit; i++) {
156 
157  SBRHit *hit = fEvent.AddHit();
158  if ( hit == 0x0 ) {
159  G4cout << " Cannot add more that " << i << " hits in the collection " << G4endl;
160  break;
161  }
162  hit->fPDG = -1 ;
163 
164  hit->fE = (*CHC)[i]->GetEdep()/CLHEP::MeV;
165  //hit->fE = (*CHC)[i]->GetEdep()/CLHEP::keV;
166  hit->fX = (*CHC)[i]->GetPos().x()/CLHEP::cm;
167  hit->fY = (*CHC)[i]->GetPos().y()/CLHEP::cm;
168  hit->fZ = (*CHC)[i]->GetPos().z()/CLHEP::cm;
169  hit->fT = (*CHC)[i]->GetToF();
170 
171  hit->fFlag = (*CHC)[i]->GetNbHits() ;
172  hit->fUID = (*CHC)[i]->GetDetID();
173 
174  H += (*CHC)[i]->GetEdep()/CLHEP::MeV;
175  //H += (*CHC)[i]->GetEdep()/CLHEP::keV;
176  }
177  fEvent.SetEMult(H,K);
178  }
179  if ( fRecordOption == 1 ) {
180  if ( fEvent.GetMultTot() > 0 ) fTree->Fill();
181  }
182  else fTree->Fill();
183 }
184 
186  SToGS::BaseROOTTreeAction(conffile),
187  fisOptical(0),
188  fOpticalEventBeg(0x0),
189  fOpticalEventEnd(0x0)
190 {
191  // open the ascii file
192  std::ifstream file(conffile.data());
193  if ( file.is_open() == false ) {
194  G4cout << " ** SToGS::BaseROOTEventsUserAction WARNING ** Cannot open file " << conffile << ", default parameters to be used "<< G4endl;
195  }
196  else {
197  std::string aline; getline(file,aline);
198  while ( file.good() && !file.eof() ) {
199 
200  if ( aline[0] == '#' ){
201  getline(file,aline);
202  continue;
203  } // this line is a comment
204 
205  std::string key, what, option; std::istringstream decode(aline); decode.clear();
206  decode >> key;
207 
208  if ( key == "scintillation:" ) { // to switch on following scintillation photons
209  decode >> fisOptical ;
210  }
211  getline(file,aline);
212  }
213  file.close();
214  }
215 }
216 
218 {
219 #ifdef G4MULTITHREADED
220  G4AutoLock lock(&buildMutex);
221 #endif
222 
223  G4Run* therun = 0x0; G4cout << " In SToGS::BaseROOTEventsUserAction, Generate a new Run " << G4endl;
224 
225  if ( fTree == 0x0 ) {
226  fTree = new TTree(fTreeName.data(),fTreeTitle.data());
227  fTree->SetDirectory(0x0);
228  // fTree->SetMaxTreeSize(fMaxEvents);
229 
230  if (fisOptical > 0) {
231  fOpticalEventBeg = new SBROpticalEvent(); fOpticalEventEnd = new SBROpticalEvent();
232  //
233  fTree->Branch("OpticalEvBegin.", fOpticalEventBeg);
234  fTree->Branch("OpticalEvEnd.", fOpticalEventEnd);
235  }
236  }
237 
238  // creates a new run. As the file is open bu AsciiRun no need to open something for master ... maybe one day to keep some globals ?
239 #if G4MULTITHREADED
240  if ( G4Threading::IsWorkerThread() ) {
241  // if ( IsMaster() ) {
242  BaseROOTEventsRun *loc_therun = new BaseROOTEventsRun(fTree);
243  therun = loc_therun;
244  }
245  else therun = G4UserRunAction::GenerateRun();
246 #else
247  therun = new BaseROOTEventsRun(fTree);
248 #endif
249  return therun;
250 }
251 
253 {
255  if (fisOptical > 0) {
256  delete fOpticalEventBeg; delete fOpticalEventEnd;
257  }
258 }
259 
261 {
263 
264  // clear events
265  if ( fisOptical > 0 ) {
266  fOpticalEventBeg->Clear("");
267  fOpticalEventEnd->Clear("");
268  }
269 }
270 
272 {
273  if ( fisOptical < 1 ) {
274  return;
275  }
276  //Count what process generated the optical photons
277  if(atrack->GetDefinition()==G4OpticalPhoton::OpticalPhotonDefinition()){
278  // particle is optical photon
279  if(atrack->GetParentID()>0){
280  // particle is secondary
281  if(atrack->GetCreatorProcess()->GetProcessName()=="Scintillation") {
282  // G4cout << "Scintillation Pre " << atrack->GetTrackID() << " " << atrack->GetTrackLength() << G4endl;
283 
284  SToGS::PrimaryTrackInformation *pinfo = (SToGS::PrimaryTrackInformation *)atrack->GetUserInformation();
285  G4int pid = atrack->GetTrackID();
286  if ( pinfo ) {
287  pid = pinfo->GetPrimaryID();
288  }
289 
290  SBROpticalHit *hit = fOpticalEventBeg->AddHit();
291 
292  hit->fX = atrack->GetPosition().x();
293  hit->fY = atrack->GetPosition().y();
294  hit->fZ = atrack->GetPosition().z();
295 
296  hit->fTA = atrack->GetGlobalTime();
297  hit->fTL = atrack->GetLocalTime();
298 
299  hit->fLength = atrack->GetTrackLength();
300  hit->fNbSteps = atrack->GetCurrentStepNumber();
301 
302  hit->fPrimaryID = pid-1;
303  hit->fSecondaryID = atrack->GetVolume()->GetCopyNo(); // hit->fUID
304 
305  }
306  }
307  }
308 }
310 {
311  if ( fisOptical < 1 ) {
312  return;
313  }
314  //Count what process generated the optical photons
315  if(atrack->GetDefinition()==G4OpticalPhoton::OpticalPhotonDefinition()){
316  // particle is optical photon
317  if(atrack->GetParentID()>0){
318  // particle is secondary
319  if(atrack->GetCreatorProcess()->GetProcessName()=="Scintillation") {
320  // G4cout << "Scintillation Post " << atrack->GetTrackID() << " " << atrack->GetTrackLength() << G4endl;
321 
322  SToGS::PrimaryTrackInformation *pinfo = (SToGS::PrimaryTrackInformation *)atrack->GetUserInformation();
323  G4int pid = atrack->GetTrackID();
324  if ( pinfo ) {
325  pid = pinfo->GetPrimaryID();
326  }
327 
328  SBROpticalHit *hit = fOpticalEventEnd->AddHit();
329 
330  hit->fX = atrack->GetPosition().x();
331  hit->fY = atrack->GetPosition().y();
332  hit->fZ = atrack->GetPosition().z();
333 
334  hit->fTA = atrack->GetGlobalTime();
335  hit->fTL = atrack->GetLocalTime();
336 
337  hit->fLength = atrack->GetTrackLength();
338  hit->fNbSteps = atrack->GetCurrentStepNumber();
339 
340  hit->fPrimaryID = pid-1; // PrimaryID[(*THC)[i]->GetTrackID()]-1 ;
341  hit->fSecondaryID = atrack->GetVolume()->GetCopyNo(); // hit->fUID
342  }
343  }
344  }
345 }
346 
347 
348 
G4THitsCollection< TrackerHit > TrackerHitsCollection
SToGS Base Root Polarized Hit.
Double32_t fPY
virtual void PreUserTrackingAction(const G4Track *)
G4THitsCollection< CopClusterHit > CopClusterHitsCollection
virtual void EndOfRunAction(const G4Run *)
Double32_t fPZ
virtual void RecordEvent(const G4Event *evt)
Double32_t fX
BaseROOTEventsUserAction(G4String conffile="setup/SToGS_Tree_actions.conf")
Double32_t fT
G4int GetPrimaryID() const
Get the primary ID.
virtual void EndOfRunAction(const G4Run *)
Double32_t fZ
virtual void BeginOfEventAction(const G4Event *)
G4int fisOptical
true if one has also photons from scintillations
virtual void BeginOfEventAction(const G4Event *)
Double32_t fY
virtual void PostUserTrackingAction(const G4Track *)
User Track Information that kept the origine of any track.
Double32_t fPX
Double32_t fE