CMS 3D CMS Logo

CastorShowerLibrary.cc
Go to the documentation of this file.
1 // File: CastorShowerLibrary.cc
3 // Description: Shower library for CASTOR calorimeter
4 // Adapted from HFShowerLibrary
5 //
6 // Wagner Carvalho
7 //
9 
12 
13 #include "G4VPhysicalVolume.hh"
14 #include "G4Step.hh"
15 #include "G4Track.hh"
16 #include "G4ParticleTable.hh"
17 #include "Randomize.hh"
18 #include "G4PhysicalConstants.hh"
19 #include "CLHEP/Units/GlobalSystemOfUnits.h"
20 
21 //#define DebugLog
22 
24  : hf(nullptr), evtInfo(nullptr), emBranch(nullptr), hadBranch(nullptr),
25  nMomBin(0), totEvents(0), evtPerBin(0),
26  nBinsE(0),nBinsEta(0),nBinsPhi(0),
27  nEvtPerBinE(0),nEvtPerBinEta(0),nEvtPerBinPhi(0),
28  SLenergies(),SLetas(),SLphis() {
29 
30  initFile(p);
31 
32 }
33 
34 //=============================================================================================
35 
36 
38  if (hf) hf->Close();
39 }
40 
41 
42 //=============================================================================================
43 
46 //
47 // Init TFile and associated TBranch's of CASTOR Root file
48 // holding library events
49 //
51 
52  //
53  // Read PSet for Castor shower library
54  //
55 
56  edm::ParameterSet m_CS = p.getParameter<edm::ParameterSet>("CastorShowerLibrary");
57  edm::FileInPath fp = m_CS.getParameter<edm::FileInPath>("FileName");
58  std::string pTreeName = fp.fullPath();
59  std::string branchEvInfo = m_CS.getUntrackedParameter<std::string>("BranchEvt");
60  std::string branchEM = m_CS.getUntrackedParameter<std::string>("BranchEM");
61  std::string branchHAD = m_CS.getUntrackedParameter<std::string>("BranchHAD");
62  verbose = m_CS.getUntrackedParameter<bool>("Verbosity",false);
63 
64  // Open TFile
65  if (pTreeName.find(".") == 0) pTreeName.erase(0,2);
66  const char* nTree = pTreeName.c_str();
67  hf = TFile::Open(nTree);
68 
69  // Check that TFile has been successfully opened
70  if (!hf->IsOpen()) {
71  edm::LogError("CastorShower") << "CastorShowerLibrary: opening " << nTree << " failed";
72  throw cms::Exception("Unknown", "CastorShowerLibrary")
73  << "Opening of " << pTreeName << " fails\n";
74  } else {
75  edm::LogInfo("CastorShower") << "CastorShowerLibrary: opening " << nTree << " successfully";
76  }
77 
78  // Check for the TBranch holding EventInfo in "Events" TTree
79  TTree* event = (TTree *) hf ->Get("CastorCherenkovPhotons");
80  if (event) {
81  evtInfo = (TBranchObject *) event->GetBranch(branchEvInfo.c_str());
82  if (evtInfo) {
84  } else {
85  edm::LogError("CastorShower") << "CastorShowerLibrary: " << branchEvInfo.c_str()
86  << " Branch does not exit in Event";
87  throw cms::Exception("Unknown", "CastorShowerLibrary") << "Event information absent\n";
88  }
89  } else {
90  edm::LogError("CastorShower") << "CastorShowerLibrary: Events Tree does not exist";
91  throw cms::Exception("Unknown", "CastorShowerLibrary") << "Events tree absent\n";
92  }
93 
94  // Get EM and HAD Branchs
95  emBranch = (TBranchObject *) event->GetBranch(branchEM.c_str());
96  if (verbose) emBranch->Print();
97  hadBranch = (TBranchObject *) event->GetBranch(branchHAD.c_str());
98  if (verbose) hadBranch->Print();
99  edm::LogInfo("CastorShower") << "CastorShowerLibrary: Branch " << branchEM
100  << " has " << emBranch->GetEntries()
101  << " entries and Branch " << branchHAD
102  << " has " << hadBranch->GetEntries()
103  << " entries";
104 
105 }
106 
107 //=============================================================================================
108 
111 //
112 // Get EventInfo from the "TBranch* branch" of Root file
113 // holding library events
114 //
115 // Based on HFShowerLibrary::loadEventInfo
116 //
118 
120  branch->SetAddress(&eventInfo);
121  branch->GetEntry(0);
122  // Initialize shower library general parameters
123 
125 // nMomBin = eventInfo->Energy.getNBins();
126 // evtPerBin = eventInfo->Energy.getNEvtPerBin();
127 // pmom = eventInfo->Energy.getBin();
133  SLetas = eventInfo->Eta.getBin();
136  SLphis = eventInfo->Phi.getBin();
137 
138  // Convert from GeV to MeV
139  for (unsigned int i=0; i<SLenergies.size(); i++) SLenergies[i] *= GeV;
140 
141  edm::LogInfo("CastorShower") << " CastorShowerLibrary::loadEventInfo : "
142  << "\n \n Total number of events : " << totEvents
143  << "\n Number of bins (E) : " << nBinsE
144  << "\n Number of events/bin (E) : " << nEvtPerBinE
145  << "\n Number of bins (Eta) : " << nBinsEta
146  << "\n Number of events/bin (Eta) : " << nEvtPerBinEta
147  << "\n Number of bins (Phi) : " << nBinsPhi
148  << "\n Number of events/bin (Phi) : " << nEvtPerBinPhi << "\n";
149 
150  for (unsigned int i=0; i<nBinsE; i++)
151  edm::LogInfo("CastorShower") << "CastorShowerLibrary: SLenergies[" << i << "] = "
152  << SLenergies[i]/GeV << " GeV";
153  for (unsigned int i=0; i<nBinsEta; i++)
154  edm::LogInfo("CastorShower") << "CastorShowerLibrary: SLetas[" << i << "] = "
155  << SLetas[i];
156  for (unsigned int i=0; i<nBinsPhi; i++)
157  edm::LogInfo("CastorShower") << "CastorShowerLibrary: SLphis[" << i << "] = "
158  << SLphis[i] << " rad";
159 }
160 
161 //=============================================================================================
162 
163 void CastorShowerLibrary::initParticleTable(G4ParticleTable * theParticleTable) {
165 //
166 // Set particle codes according to PDG encoding
167 //
168 // Based on HFShowerLibrary::initRun
169 //
171 
172  G4String particleName;
173 
174  edm::LogInfo("CastorShower") << "CastorShowerLibrary::initParticleTable"
175  << " *** Accessing PDGEncoding ***" ;
176 
177  emPDG = theParticleTable->FindParticle(particleName="e-")->GetPDGEncoding();
178  epPDG = theParticleTable->FindParticle(particleName="e+")->GetPDGEncoding();
179  gammaPDG = theParticleTable->FindParticle(particleName="gamma")->GetPDGEncoding();
180  pi0PDG = theParticleTable->FindParticle(particleName="pi0")->GetPDGEncoding();
181  etaPDG = theParticleTable->FindParticle(particleName="eta")->GetPDGEncoding();
182  nuePDG = theParticleTable->FindParticle(particleName="nu_e")->GetPDGEncoding();
183  numuPDG = theParticleTable->FindParticle(particleName="nu_mu")->GetPDGEncoding();
184  nutauPDG = theParticleTable->FindParticle(particleName="nu_tau")->GetPDGEncoding();
185  anuePDG = theParticleTable->FindParticle(particleName="anti_nu_e")->GetPDGEncoding();
186  anumuPDG = theParticleTable->FindParticle(particleName="anti_nu_mu")->GetPDGEncoding();
187  anutauPDG = theParticleTable->FindParticle(particleName="anti_nu_tau")->GetPDGEncoding();
188  geantinoPDG = theParticleTable->FindParticle(particleName="geantino")->GetPDGEncoding();
189  mumPDG = theParticleTable->FindParticle(particleName="mu-")->GetPDGEncoding();
190  mupPDG = theParticleTable->FindParticle(particleName="mu+")->GetPDGEncoding();
191 
192 //#ifdef DebugLog
193  LogDebug("CastorShower") << "CastorShowerLibrary: Particle codes for e- = " << emPDG
194  << ", e+ = " << epPDG << ", gamma = " << gammaPDG
195  << ", pi0 = " << pi0PDG << ", eta = " << etaPDG
196  << ", geantino = " << geantinoPDG << "\n nu_e = "
197  << nuePDG << ", nu_mu = " << numuPDG << ", nu_tau = "
198  << nutauPDG << ", anti_nu_e = " << anuePDG
199  << ", anti_nu_mu = " << anumuPDG << ", anti_nu_tau = "
200  << anutauPDG << ", mu- = " << mumPDG << ", mu+ = "
201  << mupPDG;
202 //#endif
203 
204  edm::LogInfo("CastorShower") << " ***** Successfully called: "
205  << "CastorShowerLibrary::initParticleTable() ***** " ;
206 
207 }
208 
209 
210 //=============================================================================================
211 
213 
214  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
215  G4Track * track = aStep->GetTrack();
216 
217  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
218  G4String partType = track->GetDefinition()->GetParticleName();
219  int parCode = track->GetDefinition()->GetPDGEncoding();
220 
222  hit.Clear();
223 
224  ok = false;
225  if (parCode == pi0PDG || parCode == etaPDG || parCode == nuePDG ||
226  parCode == numuPDG || parCode == nutauPDG || parCode == anuePDG ||
227  parCode == anumuPDG || parCode == anutauPDG || parCode == geantinoPDG)
228  return hit;
229  ok = true;
230 
231  double pin = preStepPoint->GetTotalEnergy();
232  double zint = hitPoint.z();
233  double R=sqrt(hitPoint.x()*hitPoint.x() + hitPoint.y()*hitPoint.y());
234  double theta = atan2(R,std::abs(zint));
235  double phiin = atan2(hitPoint.y(),hitPoint.x());
236  double etain = -std::log(std::tan((pi-theta)*0.5));
237 
238  // Replace "interpolation/extrapolation" by new method "select" that just randomly
239  // selects a record from the appropriate energy bin and fills its content to
240  // "showerEvent" instance of class CastorShowerEvent
241 
242  if (parCode == emPDG || parCode == epPDG || parCode == gammaPDG ) {
243  select(0, pin, etain, phiin);
244  // if (pin<SLenergies[nBinsE-1]) {
245  // interpolate(0, pin);
246  // } else {
247  // extrapolate(0, pin);
248  // }
249  } else {
250  select(1, pin, etain, phiin);
251  // if (pin<SLenergies[nBinsE-1]) {
252  // interpolate(1, pin);
253  // } else {
254  // extrapolate(1, pin);
255  // }
256  }
257 
258  hit = (*showerEvent);
259  return hit;
260 
261 }
262 
263 //=============================================================================================
264 
267 //
268 // Retrieve event # "record" from the library and stores it
269 // into vector<CastorShowerHit> showerHit;
270 //
271 // Based on HFShowerLibrary::getRecord
272 //
273 // Modified 02/02/09 by W. Carvalho
274 //
276 
277 #ifdef DebugLog
278  LogDebug("CastorShower") << "CastorShowerLibrary::getRecord: ";
279 #endif
280  int nrc = record;
282  if (type > 0) {
283  hadBranch->SetAddress(&showerEvent);
284  hadBranch->GetEntry(nrc);
285  } else {
286  emBranch->SetAddress(&showerEvent);
287  emBranch->GetEntry(nrc);
288  }
289 
290 #ifdef DebugLog
291  int nHit = showerEvent->getNhit();
292  LogDebug("CastorShower") << "CastorShowerLibrary::getRecord: Record " << record
293  << " of type " << type << " with " << nHit
294  << " CastorShowerHits";
295 
296 #endif
297 
298 }
299 
300 //=======================================================================================
301 
302 void CastorShowerLibrary::select(int type, double pin, double etain, double phiin) {
304 //
305 // Selects an event from the library based on
306 //
307 // type: 0 --> em
308 // >0 --> had
309 // pin : momentum
310 // etain: eta (if not given, disregard the eta binning
311 // phiin: phi (if not given, disregard the phi binning
312 //
313 // Created 30/01/09 by W. Carvalho
314 //
316 
317  int irec; // to hold record number
318  double r = G4UniformRand(); // to randomly select within an energy bin (r=[0,1])
319 
320  // Randomly select a record from the right energy bin in the library,
321  // based on track momentum (pin)
322 
323  // pin < SLenergies[MIN]
324 /*
325  if (pin<SLenergies[0]) {
326  edm::LogWarning("CastorShower") << "CastorShowerLibrary: Warning, pin = " << pin
327  << " less than minimum SLenergies " << SLenergies[0] << " in library."
328  << " For the moment, selecting hit from the first bin" ;
329  irec = int(nEvtPerBinE*r) + 1;
330  // pin > SLenergies[MAX]
331  } else if (pin>SLenergies[nBinsE]) {
332 
333 // This part needs rethinking because the last element of SLenergies is no longer the upper limit of the bins
334  edm::LogWarning("CastorShower") << "CastorShowerLibrary: Warning, pin = " << pin
335  << " greater than maximum SLenergies " << SLenergies[nBinsE] << " in library."
336  << " For the moment, selecting hit from the last bin";
337  irec = (nBinsE-1)*nEvtPerBinE + int(evtPerBin*r) + 1;
338  // SLenergies[MIN] < pin < SLenergies[MAX]
339  } else {
340  for (unsigned int j=0; j<nBinsE; j++) {
341  if (pin >= SLenergies[j] && pin < SLenergies[j+1]) {
342  irec = j*evtPerBin + int(evtPerBin*r) + 1 ;
343  if (irec<0) {
344  edm::LogWarning("CastorShower") << "CastorShowerLibrary:: Illegal irec = "
345  << irec << " now set to 0";
346  irec = 0;
347  } else if (irec > totEvents) {
348  edm::LogWarning("CastorShower") << "CastorShowerLibrary:: Illegal irec = "
349  << irec << " now set to "<< totEvents;
350  irec = totEvents;
351  }
352  }
353  }
354  }
355 */
356  int ienergy = FindEnergyBin(pin);
357  int ieta = FindEtaBin(etain);
358 #ifdef DebugLog
359  if (verbose) edm::LogInfo("CastorShower") << " ienergy = " << ienergy ;
360  if (verbose) edm::LogInfo("CastorShower") << " ieta = " << ieta;
361 #endif
362 
363  int iphi;
364  double phiMin = 0. ;
365  double phiMax = M_PI/4 ; // 45 * (pi/180) rad
366  if(phiin < phiMin) phiin = phiin + M_PI ;
367  if(phiin >= phiMin && phiin < phiMax) {
368  iphi = FindPhiBin(phiin) ;
369  } else {
370  double remainder = fmod(phiin , M_PI/4) ;
371  phiin = phiMin + remainder ;
372  iphi = FindPhiBin(phiin) ;
373  }
374 #ifdef DebugLog
375  if (verbose) edm::LogInfo("CastorShower") << " iphi = " << iphi;
376 #endif
377  if (ienergy==-1) ienergy=0; // if pin < first bin, choose an event in the first one
378  if (ieta==-1) ieta=0; // idem for eta
379  if (iphi!=-1) irec = int(nEvtPerBinE*ienergy+nEvtPerBinEta*ieta+nEvtPerBinPhi*(iphi+r));
380  else irec = int(nEvtPerBinE*(ienergy+r));
381 
382 #ifdef DebugLog
383  edm::LogInfo("CastorShower") << "CastorShowerLibrary:: Select record " << irec
384  << " of type " << type ;
385 #endif
386 
387  // Retrieve record number "irec" from the library
388  getRecord (type, irec);
389 
390 }
392  //
393  // returns the integer index of the energy bin, taken from SLenergies vector
394  // returns -1 if ouside valid range
395  //
396  if (energy >= SLenergies.back()) return SLenergies.size()-1;
397 
398  unsigned int i = 0;
399  for(;i<SLenergies.size()-1;i++)
400  if (energy >= SLenergies.at(i) && energy < SLenergies.at(i+1)) return (int)i;
401 
402  // now i points to the last but 1 bin
403  if (energy>=SLenergies.at(i)) return (int)i;
404  // energy outside bin range
405  return -1;
406 }
408  //
409  // returns the integer index of the eta bin, taken from SLetas vector
410  // returns -1 if ouside valid range
411  //
412  if (eta>=SLetas.back()) return SLetas.size()-1;
413  unsigned int i = 0;
414  for(;i<SLetas.size()-1;i++)
415  if (eta >= SLetas.at(i) && eta < SLetas.at(i+1)) return (int)i;
416  // now i points to the last but 1 bin
417  if (eta>=SLetas.at(i)) return (int)i;
418  // eta outside bin range
419  return -1;
420 }
422  //
423  // returns the integer index of the phi bin, taken from SLphis vector
424  // returns -1 if ouside valid range
425  //
426  // needs protection in case phi is outside range -pi,pi
427  //
428  if (phi>=SLphis.back()) return SLphis.size()-1;
429  unsigned int i = 0;
430  for(;i<SLphis.size()-1;i++)
431  if (phi >= SLphis.at(i) && phi < SLphis.at(i+1)) return (int)i;
432  // now i points to the last but 1 bin
433  if (phi>=SLphis.at(i)) return (int)i;
434  // phi outside bin range
435  return -1;
436 }
#define LogDebug(id)
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
void initParticleTable(G4ParticleTable *)
T getUntrackedParameter(std::string const &, T const &) const
const double GeV
Definition: MathUtil.h:16
void select(int, double, double=0, double=9)
JetCorrectorParameters::Record record
Definition: classes.h:7
Geom::Theta< T > theta() const
#define nullptr
std::vector< double > SLetas
const Double_t pi
unsigned int getNBins()
TBranchObject * emBranch
std::vector< double > SLphis
TBranchObject * hadBranch
T sqrt(T t)
Definition: SSEVec.h:18
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CastorShowerLibrary(std::string &name, edm::ParameterSet const &p)
unsigned int getNEvts()
#define M_PI
unsigned int getNEvtPerBin()
void initFile(edm::ParameterSet const &)
std::vector< double > SLenergies
CastorShowerLibraryInfo * eventInfo
void loadEventInfo(TBranchObject *)
CastorShowerEvent getShowerHits(G4Step *, bool &)
CastorShowerEvent * showerEvent
TBranchObject * evtInfo
double getBin(int i)
unsigned int getNhit()
Definition: event.py:1