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 
13 
14 #include "G4VPhysicalVolume.hh"
15 #include "G4Step.hh"
16 #include "G4Track.hh"
17 #include "G4ParticleTable.hh"
18 #include "Randomize.hh"
19 #include "CLHEP/Units/PhysicalConstants.h"
20 #include "CLHEP/Units/SystemOfUnits.h"
21 
22 //ROOT
23 #include "TROOT.h"
24 #include "TFile.h"
25 #include "TTree.h"
26 #include "TBranchObject.h"
27 
28 //#define DebugLog
29 
31  : hf(nullptr), evtInfo(nullptr), emBranch(nullptr), hadBranch(nullptr),
32  nMomBin(0), totEvents(0), evtPerBin(0),
33  nBinsE(0),nBinsEta(0),nBinsPhi(0),
34  nEvtPerBinE(0),nEvtPerBinEta(0),nEvtPerBinPhi(0),
35  SLenergies(),SLetas(),SLphis() {
36 
37  initFile(p);
38 }
39 
40 //=============================================================================================
41 
43  if (hf) hf->Close();
44 }
45 
46 
47 //=============================================================================================
48 
51 //
52 // Init TFile and associated TBranch's of CASTOR Root file
53 // holding library events
54 //
56 
57  //
58  // Read PSet for Castor shower library
59  //
60 
61  edm::ParameterSet m_CS = p.getParameter<edm::ParameterSet>("CastorShowerLibrary");
62  edm::FileInPath fp = m_CS.getParameter<edm::FileInPath>("FileName");
63  std::string pTreeName = fp.fullPath();
64  std::string branchEvInfo = m_CS.getUntrackedParameter<std::string>("BranchEvt");
65  std::string branchEM = m_CS.getUntrackedParameter<std::string>("BranchEM");
66  std::string branchHAD = m_CS.getUntrackedParameter<std::string>("BranchHAD");
67  verbose = m_CS.getUntrackedParameter<bool>("Verbosity",false);
68 
69  // Open TFile
70  if (pTreeName.find(".") == 0) pTreeName.erase(0,2);
71  const char* nTree = pTreeName.c_str();
72  hf = TFile::Open(nTree);
73 
74  // Check that TFile has been successfully opened
75  if (!hf->IsOpen()) {
76  edm::LogError("CastorShower") << "CastorShowerLibrary: opening " << nTree << " failed";
77  throw cms::Exception("Unknown", "CastorShowerLibrary")
78  << "Opening of " << pTreeName << " fails\n";
79  } else {
80  edm::LogInfo("CastorShower") << "CastorShowerLibrary: opening " << nTree << " successfully";
81  }
82 
83  // Check for the TBranch holding EventInfo in "Events" TTree
84  TTree* event = (TTree *) hf ->Get("CastorCherenkovPhotons");
85  if (event) {
86  evtInfo = (TBranchObject *) event->GetBranch(branchEvInfo.c_str());
87  if (evtInfo) {
89  } else {
90  edm::LogError("CastorShower") << "CastorShowerLibrary: " << branchEvInfo.c_str()
91  << " Branch does not exit in Event";
92  throw cms::Exception("Unknown", "CastorShowerLibrary") << "Event information absent\n";
93  }
94  } else {
95  edm::LogError("CastorShower") << "CastorShowerLibrary: Events Tree does not exist";
96  throw cms::Exception("Unknown", "CastorShowerLibrary") << "Events tree absent\n";
97  }
98 
99  // Get EM and HAD Branchs
100  emBranch = (TBranchObject *) event->GetBranch(branchEM.c_str());
101  if (verbose) emBranch->Print();
102  hadBranch = (TBranchObject *) event->GetBranch(branchHAD.c_str());
103  if (verbose) hadBranch->Print();
104  edm::LogInfo("CastorShower") << "CastorShowerLibrary: Branch " << branchEM
105  << " has " << emBranch->GetEntries()
106  << " entries and Branch " << branchHAD
107  << " has " << hadBranch->GetEntries()
108  << " entries";
109 
110 }
111 
112 //=============================================================================================
113 
116 //
117 // Get EventInfo from the "TBranch* branch" of Root file
118 // holding library events
119 //
120 // Based on HFShowerLibrary::loadEventInfo
121 //
123 
125  branch->SetAddress(&eventInfo);
126  branch->GetEntry(0);
127  // Initialize shower library general parameters
128 
135  SLetas = eventInfo->Eta.getBin();
138  SLphis = eventInfo->Phi.getBin();
139 
140  // Convert from GeV to MeV
141  for (unsigned int i=0; i<SLenergies.size(); i++) { SLenergies[i] *= CLHEP::GeV; }
142 
143  edm::LogInfo("CastorShower") << " CastorShowerLibrary::loadEventInfo : "
144  << "\n \n Total number of events : " << totEvents
145  << "\n Number of bins (E) : " << nBinsE
146  << "\n Number of events/bin (E) : " << nEvtPerBinE
147  << "\n Number of bins (Eta) : " << nBinsEta
148  << "\n Number of events/bin (Eta) : " << nEvtPerBinEta
149  << "\n Number of bins (Phi) : " << nBinsPhi
150  << "\n Number of events/bin (Phi) : " << nEvtPerBinPhi << "\n";
151 
152  for (unsigned int i=0; i<nBinsE; i++)
153  edm::LogInfo("CastorShower") << "CastorShowerLibrary: SLenergies[" << i << "] = "
154  << SLenergies[i]/CLHEP::GeV << " GeV";
155  for (unsigned int i=0; i<nBinsEta; i++)
156  edm::LogInfo("CastorShower") << "CastorShowerLibrary: SLetas[" << i << "] = "
157  << SLetas[i];
158  for (unsigned int i=0; i<nBinsPhi; i++)
159  edm::LogInfo("CastorShower") << "CastorShowerLibrary: SLphis[" << i << "] = "
160  << SLphis[i] << " rad";
161 }
162 
163 //=============================================================================================
164 
166 
167  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
168  G4Track * track = aStep->GetTrack();
169 
170  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
171 
173  hit.Clear();
174 
175  ok = false;
177  if (!isEM && !G4TrackToParticleID::isStableHadronIon(track)) {
178  return hit;
179  }
180  ok = true;
181 
182  double pin = preStepPoint->GetTotalEnergy();
183  double zint = hitPoint.z();
184  double R=sqrt(hitPoint.x()*hitPoint.x() + hitPoint.y()*hitPoint.y());
185  double theta = atan2(R,std::abs(zint));
186  double phiin = atan2(hitPoint.y(),hitPoint.x());
187  double etain = -std::log(std::tan((CLHEP::pi-theta)*0.5));
188 
189  // Replace "interpolation/extrapolation" by new method "select" that just randomly
190  // selects a record from the appropriate energy bin and fills its content to
191  // "showerEvent" instance of class CastorShowerEvent
192 
193  if (isEM) {
194  select(0, pin, etain, phiin);
195  } else {
196  select(1, pin, etain, phiin);
197  }
198 
199  hit = (*showerEvent);
200  return hit;
201 }
202 
203 //=============================================================================================
204 
207 //
208 // Retrieve event # "record" from the library and stores it
209 // into vector<CastorShowerHit> showerHit;
210 //
211 // Based on HFShowerLibrary::getRecord
212 //
213 // Modified 02/02/09 by W. Carvalho
214 //
216 
217 #ifdef DebugLog
218  LogDebug("CastorShower") << "CastorShowerLibrary::getRecord: ";
219 #endif
220  int nrc = record;
222  if (type > 0) {
223  hadBranch->SetAddress(&showerEvent);
224  hadBranch->GetEntry(nrc);
225  } else {
226  emBranch->SetAddress(&showerEvent);
227  emBranch->GetEntry(nrc);
228  }
229 
230 #ifdef DebugLog
231  int nHit = showerEvent->getNhit();
232  LogDebug("CastorShower") << "CastorShowerLibrary::getRecord: Record " << record
233  << " of type " << type << " with " << nHit
234  << " CastorShowerHits";
235 
236 #endif
237 
238 }
239 
240 //=======================================================================================
241 
242 void CastorShowerLibrary::select(int type, double pin, double etain, double phiin) {
244 //
245 // Selects an event from the library based on
246 //
247 // type: 0 --> em
248 // >0 --> had
249 // pin : momentum
250 // etain: eta (if not given, disregard the eta binning
251 // phiin: phi (if not given, disregard the phi binning
252 //
253 // Created 30/01/09 by W. Carvalho
254 //
256 
257  int irec; // to hold record number
258  double r = G4UniformRand(); // to randomly select within an energy bin (r=[0,1])
259 
260  // Randomly select a record from the right energy bin in the library,
261  // based on track momentum (pin)
262 
263  int ienergy = FindEnergyBin(pin);
264  int ieta = FindEtaBin(etain);
265 #ifdef DebugLog
266  if (verbose) edm::LogInfo("CastorShower") << " ienergy = " << ienergy ;
267  if (verbose) edm::LogInfo("CastorShower") << " ieta = " << ieta;
268 #endif
269 
270  int iphi;
271  const double phiMin = 0.;
272  const double phiMax = M_PI/4.; // 45 * (pi/180) rad
273  if(phiin < phiMin) phiin = phiin + M_PI ;
274  if(phiin >= phiMin && phiin < phiMax) {
275  iphi = FindPhiBin(phiin) ;
276  } else {
277  double remainder = fmod(phiin, phiMax) ;
278  phiin = phiMin + remainder ;
279  iphi = FindPhiBin(phiin) ;
280  }
281 #ifdef DebugLog
282  if (verbose) edm::LogInfo("CastorShower") << " iphi = " << iphi;
283 #endif
284  if (ienergy==-1) ienergy=0; // if pin < first bin, choose an event in the first one
285  if (ieta==-1) ieta=0; // idem for eta
286  if (iphi!=-1) irec = int(nEvtPerBinE*ienergy+nEvtPerBinEta*ieta+nEvtPerBinPhi*(iphi+r));
287  else irec = int(nEvtPerBinE*(ienergy+r));
288 
289 #ifdef DebugLog
290  edm::LogInfo("CastorShower") << "CastorShowerLibrary:: Select record " << irec
291  << " of type " << type ;
292 #endif
293 
294  // Retrieve record number "irec" from the library
295  getRecord (type, irec);
296 
297 }
298 
299 //=======================================================================================
300 
302  //
303  // returns the integer index of the energy bin, taken from SLenergies vector
304  // returns -1 if ouside valid range
305  //
306  if (energy >= SLenergies.back()) return SLenergies.size()-1;
307 
308  unsigned int i = 0;
309  for(;i<SLenergies.size()-1;i++)
310  if (energy >= SLenergies.at(i) && energy < SLenergies.at(i+1)) return (int)i;
311 
312  // now i points to the last but 1 bin
313  if (energy>=SLenergies.at(i)) return (int)i;
314  // energy outside bin range
315  return -1;
316 }
317 
318 //=======================================================================================
319 
321  //
322  // returns the integer index of the eta bin, taken from SLetas vector
323  // returns -1 if ouside valid range
324  //
325  if (eta>=SLetas.back()) return SLetas.size()-1;
326  unsigned int i = 0;
327  for(;i<SLetas.size()-1;i++)
328  if (eta >= SLetas.at(i) && eta < SLetas.at(i+1)) return (int)i;
329  // now i points to the last but 1 bin
330  if (eta>=SLetas.at(i)) return (int)i;
331  // eta outside bin range
332  return -1;
333 }
334 
335 //=======================================================================================
336 
338  //
339  // returns the integer index of the phi bin, taken from SLphis vector
340  // returns -1 if ouside valid range
341  //
342  // needs protection in case phi is outside range -pi,pi
343  //
344  if (phi>=SLphis.back()) return SLphis.size()-1;
345  unsigned int i = 0;
346  for(;i<SLphis.size()-1;i++)
347  if (phi >= SLphis.at(i) && phi < SLphis.at(i+1)) return (int)i;
348  // now i points to the last but 1 bin
349  if (phi>=SLphis.at(i)) return (int)i;
350  // phi outside bin range
351  return -1;
352 }
353 
354 //=======================================================================================
#define LogDebug(id)
void Clear(Option_t *option="") override
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const double GeV
Definition: MathUtil.h:16
void select(int, double, double=0, double=9)
CastorShowerLibrary(const std::string &name, edm::ParameterSet const &p)
JetCorrectorParameters::Record record
Definition: classes.h:7
Geom::Theta< T > theta() const
#define nullptr
std::vector< double > SLetas
const Double_t pi
CastorShowerEvent getShowerHits(const G4Step *, bool &)
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
unsigned int getNEvts()
#define M_PI
unsigned int getNEvtPerBin()
static bool isStableHadronIon(const G4Track *)
void initFile(edm::ParameterSet const &)
std::vector< double > SLenergies
CastorShowerLibraryInfo * eventInfo
void loadEventInfo(TBranchObject *)
CastorShowerEvent * showerEvent
TBranchObject * evtInfo
static bool isGammaElectronPositron(int pdgCode)
double getBin(int i)
unsigned int getNhit()
Definition: event.py:1