ACTARSim
ACTAR TPC Simulation Reference Guide
ActarSimROOTAnalSilRing.cc
Go to the documentation of this file.
1 // - AUTHOR: Hector Alvarez-Pol 05/2005
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 ActarSimROOTAnalSilRing
11 /// The ACTAR SilRingicon detectorpart of the ROOT Analysis
12 /////////////////////////////////////////////////////////////////
13 
15 
16 #include "G4ios.hh"
17 #include "G4RunManager.hh"
18 #include "G4SDManager.hh"
19 #include "G4Event.hh"
20 #include "G4Run.hh"
21 #include "G4Track.hh"
22 #include "G4ClassificationOfNewTrack.hh"
23 #include "G4TrackStatus.hh"
24 #include "G4Step.hh"
25 #include "G4Types.hh"
26 
27 //#include "G4PhysicalConstants.hh"
28 //#include "G4SystemOfUnits.hh"
29 
30 //ROOT INCLUDES
31 #include "TROOT.h"
32 //#include "TApplication.h"
33 //#include "TSystem.h"
34 //#include "TH1.h"
35 //#include "TH2.h"
36 //#include "TPad.h"
37 //#include "TCanvas.h"
38 #include "TTree.h"
39 #include "TFile.h"
40 #include "TClonesArray.h"
41 
42 //////////////////////////////////////////////////////////////////
43 /// Constructor
45  //The simulation file
47  simFile->cd();
48 
49  //The tree
51 
52  //The clones array of Sci Hits
53  //TODO-> Repeat here the trick SimpleHit/HitSim/Hit that is existing in ActarSimAnalSilRing.cc
54  silRingHitCA = new TClonesArray("ActarSimSilRingHit",2);
55 
56  silRingHitsBranch = eventTree->Branch("silRingHits",&silRingHitCA);
57  silRingHitsBranch->SetAutoDelete(kTRUE);
58 }
59 
60 //////////////////////////////////////////////////////////////////
61 /// Destructor. Makes nothing.
63 }
64 
65 //////////////////////////////////////////////////////////////////
66 /// Actions to perform in the silicon anal when generating the primaries
67 void ActarSimROOTAnalSilRing::GeneratePrimaries(const G4Event *anEvent){
68  if(anEvent){;} /* keep the compiler "quiet" */
69 }
70 
71 //////////////////////////////////////////////////////////////////
72 /// Actions to perform in the silicon anal at the begining of the run
74  if (aRun){;} /* keep the compiler "quiet" */
75 
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","silRing");
86  gDirectory->mkdir(dirName,dirName);
87  gDirectory->cd(dirName);
88 
89  simFile->cd();
90 }
91 
92 //////////////////////////////////////////////////////////////////
93 /// Actions to perform in the silicon anal at the begining of the event
94 void ActarSimROOTAnalSilRing::BeginOfEventAction(const G4Event *anEvent) {
95  SetTheEventID(anEvent->GetEventID());
96 }
97 
98 //////////////////////////////////////////////////////////////////
99 /// Actions to perform in the silicon anal at the beginning of the run
100 /// Defining the ActarSimSilHit from the ActarSimSilGeantHits
101 void ActarSimROOTAnalSilRing::EndOfEventAction(const G4Event *anEvent) {
102  FillingHits(anEvent);
103 }
104 
105 //////////////////////////////////////////////////////////////////
106 /// Defining the ActarSimSilHits from the ActarSimSilGeantHits
107 ///
108 /// A simple algorithm checking the number of primaries
109 /// reaching the Sil and calculating their mean parameters
110 /// taking into account their energy deposition on the silicons
111 /// NOTE that only primaries can produce Hits following this scheme
112 void ActarSimROOTAnalSilRing::FillingHits(const G4Event *anEvent) {
113  //Hit Container ID for ActarSimSilGeantHit
114  G4int hitsCollectionID =
115  G4SDManager::GetSDMpointer()->GetCollectionID("SilRingCollection");
116 
117  G4HCofThisEvent* HCofEvent = anEvent->GetHCofThisEvent();
118 
119  ActarSimSilRingGeantHitsCollection* hitsCollection =
120  (ActarSimSilRingGeantHitsCollection*) HCofEvent->GetHC(hitsCollectionID);
121 
122  //Number of ActarSimSilRingGeantHit (or steps) in the hitsCollection
123  G4int NbHits = hitsCollection->entries();
124 
125  G4int indepHits = 0; //number of Hits
126  G4int* trackIDtable; //stores the trackID of primary particles for each (valid) GeantHit
127  trackIDtable = new G4int[NbHits];
128  G4int* detIDtable; //stores the detIDs for each (valid) GeantHit
129  detIDtable = new G4int[NbHits];
130  G4int* IDtable; //stores the order in previous array for each (valid) GeantHit
131  IDtable = new G4int[NbHits];
132 
133  for (G4int i=0;i<NbHits;i++) {
134  if((*hitsCollection)[i]->GetParentID()==0) { //step from primary
135  if(indepHits==0) { //only for the first Hit
136  trackIDtable[indepHits] = (*hitsCollection)[i]->GetTrackID();
137  detIDtable[indepHits] = (*hitsCollection)[i]->GetDetID();
138  IDtable[i] = indepHits;
139  indepHits++;
140  }
141  else { // this part is never reached. Maybe because there is always only one indepHits that has parentID equals to 0.
142  for(G4int j=0; j<indepHits;j++) {
143  if( (*hitsCollection)[i]->GetTrackID() == trackIDtable[j] &&
144  (*hitsCollection)[i]->GetDetID() == detIDtable[j]) { //checking trackID and detID
145  IDtable[i] = j;
146  break; //not a new Hit
147  }
148  if(j==indepHits-1){ //we got the last hit and there was no match!
149  trackIDtable[indepHits] = (*hitsCollection)[i]->GetTrackID();
150  detIDtable[indepHits] = (*hitsCollection)[i]->GetDetID();
151  IDtable[i] = indepHits;
152  indepHits++;
153  }
154  }
155  }
156  }
157  }
158 
159  //Let us create as many ActarSimSilRingHit as independent primary particles
160  theSilRingHit = new ActarSimSilRingHit*[indepHits];
161  for (G4int i=0;i<indepHits;i++)
163 
164  //Clear the ClonesArray before filling it
165  silRingHitCA->Clear();
166 
167  //a variable to check if the Hit was already created
168  G4int* existing;
169  existing = new G4int[indepHits];
170  for(G4int i=0;i<indepHits;i++) existing[i] = 0;
171 
172  for(G4int i=0;i<NbHits;i++) {
173  if( (*hitsCollection)[i]->GetParentID()==0 ) { //step from primary
174  //the IDtable[i] contains the order in the indepHits list
175  if( existing[IDtable[i]]==0) { //if the indepHits does not exist
176  AddSilRingHit(theSilRingHit[IDtable[i]],(*hitsCollection)[i],0);
177  existing[IDtable[i]] = 1;
178  }
179  else
180  AddSilRingHit(theSilRingHit[IDtable[i]],(*hitsCollection)[i],1);
181  }
182  }
183 
184  //at the end, fill the ClonesArray
185  for (G4int i=0;i<indepHits;i++)
187 
188  delete [] trackIDtable;
189  delete [] IDtable;
190  delete [] existing;
191  delete [] detIDtable;
192  for (G4int i=0;i<indepHits;i++) delete theSilRingHit[i];
193  delete [] theSilRingHit;
194 }
195 
196 //////////////////////////////////////////////////////////////////
197 /// Function to move the information from the ActarSimSilGeantHit (a step hit)
198 /// to ActarSimSilHit (a simple primary representation in the silicon)
199 /// Two modes are possible:
200 /// mode == 0 : creation; the ActarSimSilHit is void and is
201 /// filled by the data from the ActarSimSilGeantHit
202 /// mode == 1 : addition; the ActarSimSilHit was already created
203 /// by other ActarSimSilGeantHit and some data members are updated
206  G4int mode) {
207  if(mode == 0) { //creation
208  cHit->SetDetectorID(gHit->GetDetID());
209 
210  cHit->SetXPos(gHit->GetLocalPrePos().x()/CLHEP::mm);
211  cHit->SetYPos(gHit->GetLocalPrePos().y()/CLHEP::mm);
212  cHit->SetZPos(gHit->GetLocalPrePos().z()/CLHEP::mm);
213 
214  cHit->SetTime(gHit->GetToF()/CLHEP::ns);
215  cHit->SetEnergy(gHit->GetEdep()/CLHEP::MeV);
216  cHit->SetEBeforeSil(gHit->GetEBeforeSil()/CLHEP::MeV);
217  cHit->SetEAfterSil(gHit->GetEAfterSil()/CLHEP::MeV);
218 
219  cHit->SetTrackID(gHit->GetTrackID());
220  cHit->SetEventID(GetTheEventID());
221  cHit->SetRunID(GetTheRunID());
222 
223  cHit->SetParticleID(gHit->GetParticleID());
224  cHit->SetParticleCharge(gHit->GetParticleCharge());
225  cHit->SetParticleMass(gHit->GetParticleMass());
226 
227  cHit->SetStepsContributing(1);
228  }
229 
230  else if(mode==1){ //addition
231  cHit->SetEnergy(cHit->GetEnergy() + gHit->GetEdep());
232  //taking the smaller outgoing energy of the geantHits
233  if(cHit->GetEAfterSil()>gHit->GetEAfterSil()) cHit->SetEAfterSil(gHit->GetEAfterSil()/CLHEP::MeV);
234  //taking the larger incoming energy of the geantHits
235  if(cHit->GetEBeforeSil()<gHit->GetEBeforeSil()) cHit->SetEBeforeSil(gHit->GetEBeforeSil()/CLHEP::MeV);
236 
238  // The mean value of a distribution {x_i} can also be computed iteratively
239  // if the values x_i are drawn one-by-one. After a new value x, the new mean is:
240  // mean(t) = mean(t-1) + (1/t)(x-mean(t-1))
241  cHit->SetXPos(cHit->GetXPos() +
242  (gHit->GetLocalPrePos().x()-cHit->GetXPos())/((G4double)cHit->GetStepsContributing()));
243  cHit->SetYPos(cHit->GetYPos() +
244  (gHit->GetLocalPrePos().y()-cHit->GetYPos())/((G4double)cHit->GetStepsContributing()));
245  cHit->SetZPos(cHit->GetZPos() +
246  (gHit->GetLocalPrePos().z()-cHit->GetZPos())/((G4double)cHit->GetStepsContributing()));
247 
248  //taking the shorter time of the geantHits
249  if(cHit->GetTime()>gHit->GetToF()) cHit->SetTime(gHit->GetToF()/CLHEP::ns);
250  }
251 }
252 
253 //////////////////////////////////////////////////////////////////
254 /// Actions to perform in the ACTAR gas detector analysis after each step
256  if(aStep){;} /* keep the compiler "quiet" */
257 }
void SetXPos(Double_t x)
void SetParticleCharge(Double_t pdgCharge)
void SetDetectorID(Int_t det)
void SetTime(Double_t ti)
G4THitsCollection< ActarSimSilRingGeantHit > ActarSimSilRingGeantHitsCollection
void BeginOfRunAction(const G4Run *)
Actions to perform in the silicon anal at the begining of the run.
void SetStepsContributing(UInt_t step)
void SetParticleMass(Double_t pdgMass)
~ActarSimROOTAnalSilRing()
Destructor. Makes nothing.
void SetEAfterSil(Double_t ea)
void AddSilRingHit(ActarSimSilRingHit *cHit, ActarSimSilRingGeantHit *gHit, G4int mode)
void BeginOfEventAction(const G4Event *)
Actions to perform in the silicon anal at the begining of the event.
void SetYPos(Double_t y)
void SetEBeforeSil(Double_t eb)
void SetEnergy(Double_t ed)
TFile * simFile
Local pointer to simFile.
void SetRunID(UInt_t run)
void SetEventID(UInt_t ev)
void SetZPos(Double_t z)
ActarSimROOTAnalysis * gActarSimROOTAnalysis
Global pointer to this soliton.
void FillingHits(const G4Event *)
TBranch * silRingHitsBranch
Local branch for the silicon hits.
ActarSimSilRingHit ** theSilRingHit
Pointer to the hits in the silicon ring.
void EndOfEventAction(const G4Event *)
void UserSteppingAction(const G4Step *)
Actions to perform in the ACTAR gas detector analysis after each step.
void GeneratePrimaries(const G4Event *)
Actions to perform in the silicon anal when generating the primaries.
TTree * eventTree
Local pointer to the event tree.
void SetTrackID(UInt_t tr)
TClonesArray * silRingHitCA
ClonesArray of the hits in the silicon ring.
void SetParticleID(UInt_t pdgID)