ACTARSim
ACTAR TPC Simulation Reference Guide
ActarSimROOTAnalSciRing.cc
Go to the documentation of this file.
1 // - AUTHOR: Hector Alvarez-Pol 05/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 ActarSimROOTAnalSciRing
11 /// The ACTAR Scintillator detectorpart of the ROOT Analysis
12 /////////////////////////////////////////////////////////////////
13 
15 #include "ActarSimSciRingHit.hh"
17 
18 #include "G4ios.hh"
19 #include "G4RunManager.hh"
20 #include "G4SDManager.hh"
21 #include "G4Event.hh"
22 #include "G4Run.hh"
23 #include "G4Track.hh"
24 #include "G4ClassificationOfNewTrack.hh"
25 #include "G4TrackStatus.hh"
26 #include "G4Step.hh"
27 #include "G4Types.hh"
28 
29 //#include "G4PhysicalConstants.hh"
30 //#include "G4SystemOfUnits.hh"
31 
32 //ROOT INCLUDES
33 #include "TROOT.h"
34 //#include "TApplication.h"
35 //#include "TSystem.h"
36 //#include "TH1.h"
37 //#include "TH2.h"
38 //#include "TPad.h"
39 //#include "TCanvas.h"
40 #include "TTree.h"
41 #include "TFile.h"
42 #include "TClonesArray.h"
43 
44 //////////////////////////////////////////////////////////////////
45 /// Constructor
47  //The simulation file
49  simFile->cd();
50 
51  //The tree
53 
54  //TODO->Implement SimHits/Hits duality as R3B...
55 
56  //The clones array of SciRing Hits
57  sciRingHitCA = new TClonesArray("ActarSimSciRingHit",2);
58 
59  sciRingHitsBranch = eventTree->Branch("sciRingHits",&sciRingHitCA);
60 }
61 
62 //////////////////////////////////////////////////////////////////
63 /// Destructor. Makes nothing.
65 }
66 
67 //////////////////////////////////////////////////////////////////
68 /// Actions to perform in the SciRingntillator anal when generating the primaries
69 void ActarSimROOTAnalSciRing::GeneratePrimaries(const G4Event *anEvent){
70  if(anEvent){;} // to quiet the compiler
71 }
72 
73 //////////////////////////////////////////////////////////////////
74 /// Actions to perform in the SciRingntillator anal at the begining of the run
76  //Storing the runID
77  SetTheRunID(aRun->GetRunID());
78 
79  char newDirName[255];
80  sprintf(newDirName,"%s%i","Histos",aRun->GetRunID());
81  simFile->cd(newDirName);
82 
83  dirName = new char[255];
84 
85  sprintf(dirName,"%s","sciRing");
86  gDirectory->mkdir(dirName,dirName);
87  gDirectory->cd(dirName);
88 
89  simFile->cd();
90 }
91 
92 //////////////////////////////////////////////////////////////////
93 /// Actions to perform in the scintillator anal at the begining of the event
94 void ActarSimROOTAnalSciRing::BeginOfEventAction(const G4Event *anEvent) {
95  SetTheEventID(anEvent->GetEventID());
96 }
97 
98 //////////////////////////////////////////////////////////////////
99 /// Actions to perform in the scintillator anal at the beginning of the run
100 void ActarSimROOTAnalSciRing::EndOfEventAction(const G4Event *anEvent) {
101  FillingHits(anEvent);
102 }
103 
104 //////////////////////////////////////////////////////////////////
105 /// Actions to perform in the scintillator detector analysis after each step
107  if(aStep){;} // to quiet the compiler
108 }
109 
110 //////////////////////////////////////////////////////////////////
111 /// Defining the ActarSimSciHits from the ActarSimSciGeantHits
112 void ActarSimROOTAnalSciRing::FillingHits(const G4Event *anEvent) {
113  //Hit Container ID for ActarSimSciGeantHit
114  G4int hitsCollectionID =
115  G4SDManager::GetSDMpointer()->GetCollectionID("SciRingCollection");
116 
117  G4HCofThisEvent* HCofEvent = anEvent->GetHCofThisEvent();
118 
119  ActarSimSciRingGeantHitsCollection* hitsCollection =
120  (ActarSimSciRingGeantHitsCollection*) HCofEvent->GetHC(hitsCollectionID);
121 
122  //Number of R3BCalGeantHit (or steps) in the hitsCollection
123  G4int NbHits = hitsCollection->entries();
124  G4int NbHitsWithSomeEnergy = NbHits;
125  //G4cout << " NbHits: " << NbHits << G4endl;
126 
127  //We accept edep=0 GeantHits, we have to remove them
128  //from the total number of GeantHits accepted for creating a real CrystalHit
129  for (G4int i=0;i<NbHits;i++)
130  if((*hitsCollection)[i]->GetEdep()==0.)
131  NbHitsWithSomeEnergy--;
132  //G4cout << " NbHitsWithSomeEnergy: " << NbHitsWithSomeEnergy << G4endl;
133 
134  G4int NbCrystalsWithHit=NbHitsWithSomeEnergy;
135  G4String* nameWithSomeEnergy;
136  G4int* detIDWithSomeEnergy;
137  detIDWithSomeEnergy = new G4int[NbHitsWithSomeEnergy];
138  nameWithSomeEnergy = new G4String[NbHitsWithSomeEnergy];
139 
140  //keep the crystal identifiers of the GeantHits with some energy
141  G4int hitsWithEnergyCounter=0;
142  for (G4int i=0;i<NbHits;i++) {
143  if((*hitsCollection)[i]->GetEdep()>0.){
144  nameWithSomeEnergy[hitsWithEnergyCounter] = (*hitsCollection)[i]->GetDetName();
145  detIDWithSomeEnergy[hitsWithEnergyCounter] = (*hitsCollection)[i]->GetDetID();
146  //G4cout << "with some energy: name:" << nameWithSomeEnergy[hitsWithEnergyCounter]
147  // << " detID:"<<detIDWithSomeEnergy[hitsWithEnergyCounter] << G4endl;
148  hitsWithEnergyCounter++;
149  }
150  }
151 
152  //if(hitsWithEnergyCounter != NbHitsWithSomeEnergy) {
153  // G4cout << "ERROR in ActarSimROOTAnalSci::FillingHits(): hitsWithEnergyCounter!=NbHitsWithSomeEnergy" << G4endl;
154  // G4cout << " hitsWithEnergyCounter = " << hitsWithEnergyCounter << " while NbHitsWithSomeEnergy = " << NbHitsWithSomeEnergy << G4endl;
155  //}
156 
157  //We ask for the number of crystals with signal and
158  //this loop calculates them from the number of ActarSimSciGeantHits with energy
159  for (G4int i=0;i<NbHitsWithSomeEnergy;i++) {
160  for(G4int j=0;j<i;j++){
161  if(nameWithSomeEnergy[i] == nameWithSomeEnergy[j] && detIDWithSomeEnergy[i] == detIDWithSomeEnergy[j]){
162  NbCrystalsWithHit--;
163  break; //break stops the second for() sentence as soon as one repetition is found
164  }
165  }
166  }
167  //G4cout << " NbCrystalsWithHit: " << NbCrystalsWithHit << G4endl;
168 
169  G4String* nameWithHit;
170  G4int* detIDWithHit;
171  detIDWithHit = new G4int[NbCrystalsWithHit];
172  nameWithHit = new G4String[NbCrystalsWithHit];
173  hitsWithEnergyCounter=0;
174 
175  //keep the crystal identifiers of the final crystalHits
176  G4int hitsCounter=0;
177  for (G4int i=0;i<NbHits;i++) {
178  if ((*hitsCollection)[i]->GetEdep()>0.) {
179  if(hitsCounter==0) { //the first geantHit with energy always creates a newHit
180  nameWithHit[hitsCounter] = (*hitsCollection)[i]->GetDetName();
181  detIDWithHit[hitsCounter] = (*hitsCollection)[i]->GetDetID();
182  //G4cout << "with hit: name:" << nameWithHit[hitsCounter]
183  // << " detID:"<<detIDWithHit[hitsCounter] << G4endl;
184  hitsCounter++;
185  }
186  else {
187  for (G4int j=0;j<hitsCounter;j++) {
188  if ( (*hitsCollection)[i]->GetDetName() == nameWithHit[j] &&
189  (*hitsCollection)[i]->GetDetID() == detIDWithHit[j] ) {
190  break; //break stops the second for() sentence as soon as one repetition is found
191  }
192  if(j==hitsCounter-1){ //when the previous if was never true, there is no repetition
193  nameWithHit[hitsCounter] = (*hitsCollection)[i]->GetDetName();
194  detIDWithHit[hitsCounter] = (*hitsCollection)[i]->GetDetID();
195  //G4cout << "with hit: name:" << nameWithHit[hitsCounter]
196  // << " detID:"<<detIDWithHit[hitsCounter] << G4endl;
197  hitsCounter++;
198  }
199  }
200  }
201  }
202  }
203  //DELETE ME AFTER THE PRELIMINARY TESTS!
204  //if(hitsCounter != NbCrystalsWithHit) {
205  // G4cout << "ERROR in R3BROOTAnalCal::FillingHits(): hitsCounter!=NbCrystalsWithHit" << G4endl;
206  // G4cout << " hitsCounter = " << hitsCounter << " while NbCrystalsWithHit = " << NbCrystalsWithHit << G4endl;
207  //}
208 
209  //Let us create as many R3BCalCrystalHit as NbCrystalsWithHit
210  //TODO->Recover here the simhit/hit duality if needed!!
211  //if(((ActarSimROOTAnalysis*)gActarSimROOTAnalysis)->GetUseCrystalHitSim()==0){
212  theSciRingHit = new ActarSimSciRingHit*[NbCrystalsWithHit];
213  for (G4int i=0;i<NbCrystalsWithHit;i++)
215  //}
216  //else{ //In case that the simulated hits were used!
217  // theSciHit = new ActarSimSciHit*[NbCrystalsWithHit];
218  // for (G4int i=0;i<NbCrystalsWithHit;i++)
219  // theSciHit[i] = new ActarSimSciHitSim();
220  //}
221 
222  //Clear the ClonesArray before filling it
223  sciRingHitCA->Clear();
224 
225  G4bool counted = 0;
226  hitsCounter = 0; //now this variable is going to count the added GeantHits
227  G4String* name;
228  G4int* detID;
229  detID = new G4int[NbCrystalsWithHit];
230  name = new G4String[NbCrystalsWithHit];
231 
232  G4bool shouldThisBeStored = 0;
233  //Filling the ActarSimSciHit from the R3BCalGeantHit (or step)
234  for (G4int i=0;i<NbHits;i++) {
235  //do not accept GeantHits with edep=0
236  //if there is no other GeantHit with edep>0 in the same crystal!
237  shouldThisBeStored=0;
238  counted =0;
239  //G4cout << "ADDING THE HITS. GeantHit with name:" << (*hitsCollection)[i]->GetDetName()
240  // << " detID:"<< (*hitsCollection)[i]->GetDetID() << " edep:"<< (*hitsCollection)[i]->GetEdep() <<" under consideration"<< G4endl;
241  for (G4int j=0;j<NbCrystalsWithHit;j++) {
242  if( (*hitsCollection)[i]->GetDetName() == nameWithHit[j] &&
243  (*hitsCollection)[i]->GetDetID() == detIDWithHit[j] ) {
244  shouldThisBeStored=1;
245  break; //break stops the for() sentence as soon as one "energetic" partner is found
246  }
247  }
248  if(!shouldThisBeStored) continue; //so, continue with next, if there is no "energetic partner"
249  if(hitsCounter==0){ //only for the first geantHit accepted for storage
250  name[hitsCounter] = (*hitsCollection)[i]->GetDetName();
251  detID[hitsCounter] = (*hitsCollection)[i]->GetDetID();
252  AddCalCrystalHit(theSciRingHit[hitsCounter],(*hitsCollection)[i],0);
253  //G4cout << "ADD hit: name:" << name[hitsCounter]
254  // << " detID:"<<detID[hitsCounter] << " with code 0 (initial)" << G4endl;
255  hitsCounter++;
256  }
257  else { // from the second in advance, compare if a R3BCalCrystalHit
258  // in the same volume already exists
259  for (G4int j=0;j<hitsCounter;j++) {
260  if( (*hitsCollection)[i]->GetDetName() == name[j] &&
261  (*hitsCollection)[i]->GetDetID() == detID[j] ){ //identical Hit already present
262  AddCalCrystalHit(theSciRingHit[j],(*hitsCollection)[i],1);
263  //G4cout << "ADD hit: name:" << name[j]
264  // << " detID:"<< detID[j] << " with code 1" << G4endl;
265  counted = 1;
266  break; // stops the for() loop; the info is written in the first identical Hit
267  }
268  }
269  if(counted==0) { //No identical Hit present.
270  name[hitsCounter] = (*hitsCollection)[i]->GetDetName();
271  detID[hitsCounter] = (*hitsCollection)[i]->GetDetID();
272  AddCalCrystalHit(theSciRingHit[hitsCounter],(*hitsCollection)[i],0);
273  //G4cout << "ADD hit: name:" << name[hitsCounter]
274  // << " detID:"<<detID[hitsCounter] << " with code 0" << G4endl;
275  hitsCounter++;
276  }
277  }
278  }
279 
280  //BORRAME
281  //G4cout << G4endl <<"************************ END OF DEBUGGING *****************************************" << G4endl;
282  //G4cout <<"************************************************************************************" << G4endl;
283 
284  //at the end, fill the ClonesArray
285  //TODO-> Recover here the simhit/hit duality if needed!!
286  //if(((ActarSimROOTAnalysis*)gActarSimROOTAnalysis)->GetUseCrystalHitSim()==0){
287  for (G4int i=0;i<NbCrystalsWithHit;i++)
289  //}
290  //else{ //In case that the simulated hits were used!
291  // ActarSimSciHitSim* testSim;
292  // for (G4int i=0;i<NbCrystalsWithHit;i++) {
293  // testSim = (ActarSimSciHitSim*) (theSciHit[i]);
294  // new((*sciHitCA)[i])ActarSimSciHitSim(*testSim);
295  // }
296  //}
297 
298  // all branches at the same time!
299  //G4cout << " #@BITCOUNT "<< crystalHitsBranch->Fill() << G4endl;
300  //G4cout << " #@BITCOUNT "<< theR3BTree->Fill() << G4endl;
301 
302  delete [] detIDWithSomeEnergy;
303  delete [] nameWithSomeEnergy;
304  delete [] detIDWithHit;
305  delete [] nameWithHit;
306  delete [] detID;
307  delete [] name;
308  for (G4int i=0;i<NbCrystalsWithHit;i++) delete theSciRingHit[i];
309  delete [] theSciRingHit;
310 }
311 
312 //////////////////////////////////////////////////////////////////
313 /// Function to move the information from the ActarSimSciGeantHit (a step hit)
314 /// to ActarSimSciHit (an event hit) for the Darmstadt-Heidelberg Crystall Ball.
315 ///
316 /// Two modes are possible:
317 /// - mode == 0 : creation; the ActarSimSciHit is void and is
318 /// filled by the data from the ActarSimSciGeantHit
319 /// - mode == 1 : addition; the ActarSimSciHit was already created
320 /// by other ActarSimSciGeantHit and some data members are updated
323  G4int mode) {
324  if(mode == 0) { //creation
325  if( gHit->GetDetName() == "CsISector" ) cHit->SetType(1);
326  else G4cout << "ERROR in R3BROOTAnalCal::AddCalCrystalHit()." << G4endl
327  << "Unknown Detector Name: "<< gHit->GetDetName() << G4endl << G4endl;
328 
329  cHit->SetCopy(gHit->GetDetID());
330 
331  cHit->SetEnergy(gHit->GetEdep()/ CLHEP::MeV);
332  cHit->SetTime(gHit->GetToF() / CLHEP::ns);
333 
334  cHit->SetTrackID(gHit->GetTrackID());
335  cHit->SetEventID(GetTheEventID());
336  cHit->SetRunID(GetTheRunID());
337 
338  cHit->SetParticleID(gHit->GetParticleID());
339  cHit->SetParticleCharge(gHit->GetParticleCharge());
340  cHit->SetParticleMass(gHit->GetParticleMass());
341 
342  //TODO-> Recover here the simhit/hit duality if needed!!
343  /*
344  if(((ActarSimROOTAnalysis*)gActarSimROOTAnalysis)->GetUseCrystalHitSim()!=0){
345  if( fabs(gHit->GetLocalPos().z())> 120 )
346  ((ActarSimSciHitSim*)cHit)->SetEnergyPerZone(24,gHit->GetEdep()/ MeV);
347  else {
348  G4int bin = (G4int)((gHit->GetLocalPos().z()/10) + 12);
349  ((ActarSimSciHitSim*)cHit)->SetEnergyPerZone( bin,
350  gHit->GetEdep()/ MeV);
351  //G4cout << G4endl << "Initial posZ:" << gHit->GetLocalPos().z() << " bin " << bin
352  // << " E:" << ((ActarSimSciHitSim*)cHit)->GetEnergyPerZone( bin ) << G4endl << G4endl;
353  }
354  ((ActarSimSciHitSim*)cHit)->SetNbOfSteps(1);
355  ((ActarSimSciHitSim*)cHit)->SetTimeFirstStep(gHit->GetToF() / ns);
356  ((ActarSimSciHitSim*)cHit)->SetTimeLastStep(gHit->GetToF() / ns);
357  ((ActarSimSciHitSim*)cHit)->SetNbOfPrimaries(GetPrimNbOfParticles());
358  ((ActarSimSciHitSim*)cHit)->SetEnergyPrimary(GetPrimEnergy() / MeV);
359  ((ActarSimSciHitSim*)cHit)->SetThetaPrimary(GetPrimTheta() / rad);
360  ((ActarSimSciHitSim*)cHit)->SetPhiPrimary(GetPrimPhi() / rad);
361  if(gHit->GetProcessName()=="phot") ((ActarSimSciHitSim*)cHit)->SetFirstInteractionType(1);
362  else if(gHit->GetProcessName()=="compt") ((ActarSimSciHitSim*)cHit)->SetFirstInteractionType(2);
363  else if(gHit->GetProcessName()=="conv") ((ActarSimSciHitSim*)cHit)->SetFirstInteractionType(3);
364  else if(gHit->GetProcessName()=="msc") ((ActarSimSciHitSim*)cHit)->SetFirstInteractionType(4);
365  else if(gHit->GetProcessName()=="eBrem") ((ActarSimSciHitSim*)cHit)->SetFirstInteractionType(5);
366  else if(gHit->GetProcessName()=="Transportation") ((ActarSimSciHitSim*)cHit)->SetFirstInteractionType(6);
367  else ((ActarSimSciHitSim*)cHit)->SetFirstInteractionType(0);
368  ((ActarSimSciHitSim*)cHit)->SetFirstInteractionX(gHit->GetLocalPos().x());
369  ((ActarSimSciHitSim*)cHit)->SetFirstInteractionY(gHit->GetLocalPos().y());
370  ((ActarSimSciHitSim*)cHit)->SetFirstInteractionZ(gHit->GetLocalPos().z());
371  }
372  */
373  }
374  else if(mode==1){ //addition
375  cHit->SetEnergy(cHit->GetEnergy() + gHit->GetEdep()/ CLHEP::MeV);
376  if(gHit->GetToF()<cHit->GetTime()) cHit->SetTime(gHit->GetToF()/ CLHEP::ns);
377 
378  //TODO-> Recover here the simhit/hit duality if needed!!
379  /*
380  if(((ActarSimROOTAnalysis*)gActarSimROOTAnalysis)->GetUseCrystalHitSim()!=0){
381  if( fabs(gHit->GetLocalPos().z())> 120 )
382  ((R3BCalCrystalHitSim*)cHit)->SetEnergyPerZone(24,
383  ((R3BCalCrystalHitSim*)cHit)->GetEnergyPerZone(24)+gHit->GetEdep()/ MeV);
384  else {
385  G4int bin = (G4int)((gHit->GetLocalPos().z()/10) + 12);
386  ((R3BCalCrystalHitSim*)cHit)->SetEnergyPerZone( bin,
387  ((R3BCalCrystalHitSim*)cHit)->GetEnergyPerZone(bin)+gHit->GetEdep()/ MeV);
388  //G4cout << G4endl << "Adding posZ:" << gHit->GetLocalPos().z() << " bin " << bin << " E:"
389  // << ((R3BCalCrystalHitSim*)cHit)->GetEnergyPerZone( bin ) << G4endl << G4endl;
390  }
391  ((R3BCalCrystalHitSim*)cHit)->SetNbOfSteps(((R3BCalCrystalHitSim*)cHit)->GetNbOfSteps()+1);
392  if(gHit->GetToF()<((R3BCalCrystalHitSim*)cHit)->GetTime()){
393  ((R3BCalCrystalHitSim*)cHit)->SetTime(gHit->GetToF()/ ns);
394  ((R3BCalCrystalHitSim*)cHit)->SetTimeFirstStep(gHit->GetToF() / ns);
395  if(gHit->GetProcessName()=="phot") ((R3BCalCrystalHitSim*)cHit)->SetFirstInteractionType(1);
396  else if(gHit->GetProcessName()=="compt") ((R3BCalCrystalHitSim*)cHit)->SetFirstInteractionType(2);
397  else if(gHit->GetProcessName()=="conv") ((R3BCalCrystalHitSim*)cHit)->SetFirstInteractionType(3);
398  else if(gHit->GetProcessName()=="msc") ((R3BCalCrystalHitSim*)cHit)->SetFirstInteractionType(4);
399  else if(gHit->GetProcessName()=="eBrem") ((R3BCalCrystalHitSim*)cHit)->SetFirstInteractionType(5);
400  else if(gHit->GetProcessName()=="Transportation") ((R3BCalCrystalHitSim*)cHit)->SetFirstInteractionType(6);
401  else ((R3BCalCrystalHitSim*)cHit)->SetFirstInteractionType(0);
402  ((R3BCalCrystalHitSim*)cHit)->SetFirstInteractionX(gHit->GetLocalPos().x());
403  ((R3BCalCrystalHitSim*)cHit)->SetFirstInteractionY(gHit->GetLocalPos().y());
404  ((R3BCalCrystalHitSim*)cHit)->SetFirstInteractionZ(gHit->GetLocalPos().z());
405  }
406  else if(gHit->GetToF()>((R3BCalCrystalHitSim*)cHit)->GetTimeLastStep())
407  ((R3BCalCrystalHitSim*)cHit)->SetTimeLastStep(gHit->GetToF());
408  }
409 */
410  }
411 }
void BeginOfRunAction(const G4Run *)
Actions to perform in the SciRingntillator anal at the begining of the run.
void UserSteppingAction(const G4Step *)
Actions to perform in the scintillator detector analysis after each step.
void GeneratePrimaries(const G4Event *)
Actions to perform in the SciRingntillator anal when generating the primaries.
void FillingHits(const G4Event *anEvent)
Defining the ActarSimSciHits from the ActarSimSciGeantHits.
void SetParticleID(UInt_t pdgID)
G4THitsCollection< ActarSimSciRingGeantHit > ActarSimSciRingGeantHitsCollection
void SetParticleMass(Double_t pdgMass)
void SetEventID(UInt_t ev)
TFile * simFile
Local pointer to simFile.
TBranch * sciRingHitsBranch
Local branch for the scintillator ring hits.
TClonesArray * sciRingHitCA
ClonesArray of the hits in the scintillator ring.
void SetTrackID(UInt_t track)
void SetType(UInt_t ty)
ActarSimROOTAnalysis * gActarSimROOTAnalysis
Global pointer to this soliton.
void SetParticleCharge(Double_t pdgCharge)
void SetTime(Double_t ti)
void SetEnergy(Double_t ed)
void SetRunID(UInt_t run)
ActarSimSciRingHit ** theSciRingHit
Pointer to the hits in the scintillator ring.
void BeginOfEventAction(const G4Event *)
Actions to perform in the scintillator anal at the begining of the event.
void AddCalCrystalHit(ActarSimSciRingHit *, ActarSimSciRingGeantHit *, G4int)
TTree * eventTree
Local pointer to the event tree.
void SetCopy(UInt_t co)
void EndOfEventAction(const G4Event *)
Actions to perform in the scintillator anal at the beginning of the run.
~ActarSimROOTAnalSciRing()
Destructor. Makes nothing.