CMS 3D CMS Logo

HcalTB02Analysis.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HcalTestBeam
4 // Class : HcalTB02Analysis
5 //
6 // Implementation:
7 // Main analysis class for Hcal Test Beam 2002 Analysis
8 //
9 // Original Author:
10 // Created: Sun May 21 10:14:34 CEST 2006
11 //
12 
13 // system include files
14 #include <cmath>
15 #include <iostream>
16 #include <iomanip>
17 #include <memory>
18 #include <map>
19 #include <vector>
20 #include <string>
21 
22 // user include files
23 // to retreive hits
33 
38 
42 
43 #include "G4HCofThisEvent.hh"
44 #include "G4SDManager.hh"
45 #include "G4Step.hh"
46 #include "G4Track.hh"
47 #include "G4ThreeVector.hh"
48 #include "G4VProcess.hh"
49 
50 #include <CLHEP/Random/RandGaussQ.h>
51 #include <CLHEP/Random/Randomize.h>
52 #include <CLHEP/Units/GlobalSystemOfUnits.h>
53 #include <CLHEP/Units/GlobalPhysicalConstants.h>
54 
55 //#define EDM_ML_DEBUG
56 
57 namespace CLHEP {
58  class HepRandomEngine;
59 }
60 
61 class HcalTB02Analysis : public SimProducer, public Observer<const BeginOfEvent*>, public Observer<const EndOfEvent*> {
62 public:
64  HcalTB02Analysis(const HcalTB02Analysis&) = delete; // stop default
65  const HcalTB02Analysis& operator=(const HcalTB02Analysis&) = delete;
66  ~HcalTB02Analysis() override;
67 
68  void produce(edm::Event&, const edm::EventSetup&) override;
69 
70 private:
71  // observer methods
72  void update(const BeginOfEvent* evt) override;
73  void update(const EndOfEvent* evt) override;
74 
76  void clear();
77  void finish();
78 
79 private:
80  // Private Tuples
81  std::unique_ptr<HcalTB02Histo> histo;
82 
83  // to read from parameter set
85  const bool hcalOnly;
86  const std::vector<std::string> names;
87 
88  //To be saved
89  std::map<int, float> energyInScints, energyInCrystals;
90  std::map<int, float> primaries;
95  int maxTime;
100 };
101 
102 //
103 // constructors and destructor
104 //
105 
107  : m_Anal(p.getParameter<edm::ParameterSet>("HcalTB02Analysis")),
108  hcalOnly(m_Anal.getUntrackedParameter<bool>("HcalClusterOnly", true)),
109  names(m_Anal.getParameter<std::vector<std::string> >("Names")) {
110  produces<HcalTB02HistoClass>();
111 
112  edm::LogVerbatim("HcalTBSim") << "HcalTB02Analysis:: Initialised as observer of "
113  << "BeginOfJob/BeginOfEvent/EndOfEvent with "
114  << "Parameter values:\n \thcalOnly = " << hcalOnly;
115 
116  histo = std::make_unique<HcalTB02Histo>(m_Anal);
117 }
118 
120 
121 //
122 // member functions
123 //
124 
126  std::unique_ptr<HcalTB02HistoClass> product(new HcalTB02HistoClass);
127  fillEvent(*product);
128  e.put(std::move(product));
129 }
130 
132 #ifdef EDM_ML_DEBUG
133  edm::LogVerbatim("HcalTBSim") << "HcalTB02Analysis: =====> Begin of event = " << (*evt)()->GetEventID();
134 #endif
135  clear();
136 }
137 
139  CLHEP::HepRandomEngine* engine = G4Random::getTheEngine();
140  CLHEP::RandGaussQ randGauss(*engine);
141 
142  // Look for the Hit Collection
143 #ifdef EDM_ML_DEBUG
144  edm::LogVerbatim("HcalTBSim") << "HcalTB02Analysis::Fill event " << (*evt)()->GetEventID();
145 #endif
146  // access to the G4 hit collections
147  G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
148  int ihit = 0;
149 
150  // HCAL
151  std::string sd = names[0];
152  int HCHCid = G4SDManager::GetSDMpointer()->GetCollectionID(sd);
153  CaloG4HitCollection* theHCHC = (CaloG4HitCollection*)allHC->GetHC(HCHCid);
154 #ifdef EDM_ML_DEBUG
155  edm::LogVerbatim("HcalTBSim") << "HcalTB02Analysis :: Hit Collection for " << sd << " of ID " << HCHCid
156  << " is obtained at " << theHCHC;
157 #endif
158  int nentries = 0;
159  nentries = theHCHC->entries();
160  if (nentries == 0)
161  return;
162 
163  int xentries = 0;
164  int XTALSid = 0;
165  CaloG4HitCollection* theXTHC = nullptr;
166 
167  if (!hcalOnly) {
168  // XTALS
169  sd = names[1];
170  XTALSid = G4SDManager::GetSDMpointer()->GetCollectionID(sd);
171  theXTHC = (CaloG4HitCollection*)allHC->GetHC(XTALSid);
172 #ifdef EDM_ML_DEBUG
173  edm::LogVerbatim("HcalTBSim") << "HcalTB02Analysis :: Hit Collection for " << sd << " of ID " << XTALSid
174  << " is obtained at " << theXTHC;
175 #endif
176  xentries = theXTHC->entries();
177  if (xentries == 0)
178  return;
179  }
180 
181 #ifdef EDM_ML_DEBUG
182  edm::LogVerbatim("HcalTBSim") << "HcalTB02Analysis :: There are " << nentries << " HCal hits, and" << xentries
183  << " xtal hits";
184 #endif
185  float ETot = 0., xETot = 0.;
186  int scintID = 0, xtalID = 0;
187 
188  // HCAL
189  std::unique_ptr<HcalTB02HcalNumberingScheme> org(new HcalTB02HcalNumberingScheme());
190 
191  if (HCHCid >= 0 && theHCHC != nullptr) {
192  for (ihit = 0; ihit < nentries; ihit++) {
193  CaloG4Hit* aHit = (*theHCHC)[ihit];
194  scintID = aHit->getUnitID();
195  int layer = org->getlayerID(scintID);
196  float enEm = aHit->getEM();
197  float enhad = aHit->getHadr();
198 
199  if (layer == 0) {
200  enEm = enEm / 2.;
201  enhad = enhad / 2.;
202  }
203 
204  energyInScints[scintID] += enEm + enhad;
205  primaries[aHit->getTrackID()] += enEm + enhad;
206  float time = aHit->getTimeSliceID();
207  if (time > maxTime)
208  maxTime = static_cast<int>(time);
209  histo->fillAllTime(time);
210  }
211 
212  //
213  // Profile
214  //
215 
216  float TowerEne[8][18], TowerEneCF[8][18], LayerEne[19], EnRing[100];
217  for (int iphi = 0; iphi < 8; iphi++) {
218  for (int jeta = 0; jeta < 18; jeta++) {
219  TowerEne[iphi][jeta] = 0.;
220  TowerEneCF[iphi][jeta] = 0.;
221  }
222  }
223 
224  for (int ilayer = 0; ilayer < 19; ilayer++)
225  LayerEne[ilayer] = 0.;
226  for (int iring = 0; iring < 100; iring++)
227  EnRing[iring] = 0.;
228 
229  for (std::map<int, float>::iterator is = energyInScints.begin(); is != energyInScints.end(); is++) {
230  ETot = (*is).second;
231 
232  int layer = org->getlayerID((*is).first);
233 
234  if ((layer != 17) && (layer != 18)) {
235  float eta = org->getetaID((*is).first);
236  float phi = org->getphiID((*is).first);
237 
238  SEnergy += ETot;
239  int iphi = static_cast<int>(phi);
240  int ieta = static_cast<int>(eta);
241  TowerEne[iphi][ieta] += ETot;
242 
243  TowerEneCF[iphi][ieta] += ETot * (1. + 0.1 * randGauss.fire());
244  double dR = 0.08727 * std::sqrt((eta - 8.) * (eta - 8.) + (phi - 3.) * (phi - 3.));
245  EnRing[static_cast<int>(dR / 0.01)] += ETot;
246  }
247 
248  LayerEne[layer] += ETot;
249  }
250 
251  for (int ilayer = 0; ilayer < 19; ilayer++) {
252  histo->fillProfile(ilayer, LayerEne[ilayer] / GeV);
253  }
254 
255  for (int iring = 0; iring < 100; iring++)
256  histo->fillTransProf(0.01 * iring + 0.005, EnRing[iring] / GeV);
257 
258  for (int iphi = 0; iphi < 8; iphi++) {
259  for (int jeta = 0; jeta < 18; jeta++) {
260  SEnergyN += TowerEneCF[iphi][jeta] + 3. * randGauss.fire(); // QGSP
261 
262  double Rand = 3. * randGauss.fire(); // QGSP
263 
264  if ((iphi >= 0) && (iphi < 7)) {
265  if ((jeta >= 5) && (jeta < 12)) {
266  E7x7Matrix += TowerEne[iphi][jeta];
267  E7x7MatrixN += TowerEneCF[iphi][jeta] + Rand;
268 
269  if ((iphi >= 1) && (iphi < 6)) {
270  if ((jeta >= 6) && (jeta < 11)) {
271  E5x5Matrix += TowerEne[iphi][jeta];
272  E5x5MatrixN += TowerEneCF[iphi][jeta] + Rand;
273  }
274  }
275  }
276  }
277  }
278  }
279 
280  //
281  // Find Primary info:
282  //
283  int trackID = 0;
284  G4PrimaryParticle* thePrim = nullptr;
285  G4int nvertex = (*evt)()->GetNumberOfPrimaryVertex();
286 #ifdef EDM_ML_DEBUG
287  edm::LogVerbatim("HcalTBSim") << "HcalTB02Analysis :: Event has " << nvertex << " vertex";
288 #endif
289  if (nvertex == 0)
290  edm::LogWarning("HcalTBSim") << "HcalTB02Analysis:: End Of Event "
291  << "ERROR: no vertex";
292 
293  for (int i = 0; i < nvertex; i++) {
294  G4PrimaryVertex* avertex = (*evt)()->GetPrimaryVertex(i);
295  if (avertex == nullptr) {
296  edm::LogWarning("HcalTBSim") << "HcalTB02Analysis:: End Of Event "
297  << "ERROR: pointer to vertex = 0";
298  } else {
299 #ifdef EDM_ML_DEBUG
300  int npart = avertex->GetNumberOfParticle();
301  edm::LogVerbatim("HcalTBSim") << "HcalTB02Analysis::Vertex number :" << i << " with " << npart << " particles";
302 #endif
303  if (thePrim == nullptr)
304  thePrim = avertex->GetPrimary(trackID);
305  }
306  }
307 
308  double px = 0., py = 0., pz = 0.;
309 
310  if (thePrim != nullptr) {
311  px = thePrim->GetPx();
312  py = thePrim->GetPy();
313  pz = thePrim->GetPz();
314  pInit = std::sqrt(pow(px, 2.) + pow(py, 2.) + pow(pz, 2.));
315  if (pInit == 0) {
316  edm::LogWarning("HcalTBSim") << "HcalTB02Analysis:: End Of Event "
317  << " ERROR: primary has p=0 ";
318  } else {
319  float costheta = pz / pInit;
320  float theta = acos(std::min(std::max(costheta, float(-1.)), float(1.)));
321  eta = -std::log(std::tan(theta / 2));
322  if (px != 0)
323  phi = std::atan(py / px);
324  }
325  particleType = thePrim->GetPDGcode();
326  } else {
327 #ifdef EDM_ML_DEBUG
328  edm::LogVerbatim("HcalTBSim") << "HcalTB02Analysis:: End Of Event ERROR: could not find primary ";
329 #endif
330  }
331 
332  CaloG4Hit* firstHit = (*theHCHC)[0];
333  incidentEnergy = (firstHit->getIncidentEnergy() / GeV);
334 
335  } // number of Hits > 0
336 
337  if (!hcalOnly) {
338  // XTALS
339 
340  if (XTALSid >= 0 && theXTHC != nullptr) {
341  for (int xihit = 0; xihit < xentries; xihit++) {
342  CaloG4Hit* xaHit = (*theXTHC)[xihit];
343 
344  float xenEm = xaHit->getEM();
345  float xenhad = xaHit->getHadr();
346  xtalID = xaHit->getUnitID();
347 
348  energyInCrystals[xtalID] += xenEm + xenhad;
349  }
350 
351  float xCrysEne[7][7];
352  for (int irow = 0; irow < 7; irow++) {
353  for (int jcol = 0; jcol < 7; jcol++) {
354  xCrysEne[irow][jcol] = 0.;
355  }
356  }
357 
358  for (std::map<int, float>::iterator is = energyInCrystals.begin(); is != energyInCrystals.end(); is++) {
359  int xtalID = (*is).first;
360  xETot = (*is).second;
361 
362  int irow = static_cast<int>(xtalID / 100.);
363  int jcol = static_cast<int>(xtalID - 100. * irow);
364 
365  xSEnergy += xETot;
366  xCrysEne[irow][jcol] = xETot;
367 
368  float dR = std::sqrt(0.01619 * 0.01619 * (jcol - 3) * (jcol - 3) + 0.01606 * 0.01606 * (irow - 3) * (irow - 3));
369  histo->fillTransProf(dR, xETot * 1.05);
370 
371  if ((irow > 0) && (irow < 6)) {
372  if ((jcol > 0) && (jcol < 6)) {
373  xE5x5Matrix += xCrysEne[irow][jcol];
374  xE5x5MatrixN += xCrysEne[irow][jcol] + 108.5 * randGauss.fire();
375 
376  if ((irow > 1) && (irow < 5)) {
377  if ((jcol > 1) && (jcol < 5)) {
378  xE3x3Matrix += xCrysEne[irow][jcol];
379  xE3x3MatrixN += xCrysEne[irow][jcol] + 108.5 * randGauss.fire();
380  }
381  }
382  }
383  }
384  }
385 
386  if (!hcalOnly) {
387  // assert(theXTHC);
388  if (theXTHC != nullptr) {
389  CaloG4Hit* xfirstHit = (*theXTHC)[0];
390  xIncidentEnergy = xfirstHit->getIncidentEnergy() / GeV;
391  }
392  }
393 
394  } // number of Hits > 0
395  }
396 
397  int iEvt = (*evt)()->GetEventID();
398  if (iEvt < 10)
399  edm::LogVerbatim("HcalTBSim") << " Event " << iEvt;
400  else if ((iEvt < 100) && (iEvt % 10 == 0))
401  edm::LogVerbatim("HcalTBSim") << " Event " << iEvt;
402  else if ((iEvt < 1000) && (iEvt % 100 == 0))
403  edm::LogVerbatim("HcalTBSim") << " Event " << iEvt;
404  else if ((iEvt < 10000) && (iEvt % 1000 == 0))
405  edm::LogVerbatim("HcalTBSim") << " Event " << iEvt;
406 }
407 
409  //Beam information
410  product.set_Nprim(float(primaries.size()));
411  product.set_partType(particleType);
412  product.set_Einit(pInit / GeV);
413  product.set_eta(eta);
414  product.set_phi(phi);
415  product.set_Eentry(incidentEnergy);
416 
417  //Calorimeter energy
418  product.set_ETot(SEnergy / GeV);
419  product.set_E7x7(E7x7Matrix / GeV);
420  product.set_E5x5(E5x5Matrix / GeV);
421  product.set_ETotN(SEnergyN / GeV);
422  product.set_E7x7N(E7x7MatrixN / GeV);
423  product.set_E5x5N(E5x5MatrixN / GeV);
424  product.set_NUnit(float(energyInScints.size()));
425  product.set_Ntimesli(float(maxTime));
426 
427  //crystal information
428  product.set_xEentry(xIncidentEnergy);
429  product.set_xNUnit(float(energyInCrystals.size()));
430  product.set_xETot(xSEnergy / GeV);
431  product.set_xETotN(xSEnergyN / GeV);
432  product.set_xE5x5(xE5x5Matrix / GeV);
433  product.set_xE3x3(xE3x3Matrix / GeV);
434  product.set_xE5x5N(xE5x5MatrixN / GeV);
435  product.set_xE3x3N(xE3x3MatrixN / GeV);
436 }
437 
439  primaries.clear();
440  particleType = 0;
441  pInit = eta = phi = incidentEnergy = 0;
442 
444  E7x7MatrixN = E5x5MatrixN = 0;
445  energyInScints.clear();
446  maxTime = 0;
447 
448  xIncidentEnergy = 0;
449  energyInCrystals.clear();
452 }
453 
455  /*
456  //Profile
457  std::ofstream oFile;
458  oFile.open("profile.dat");
459  float st[19] = {0.8,
460  0.4, 0.4, 0.4, 0.4, 0.4,
461  0.4, 0.4, 0.4, 0.4, 0.4,
462  0.4, 0.4, 0.4, 0.4, 0.4,
463  0.8, 1.0, 1.0};
464 
465  //cm of material (brass) in front of scintillator layer i:
466 
467  float w[19] = {7.45, //iron !
468  6.00, 6.00, 6.00, 6.00, 6.00, //brass
469  6.00, 6.00, 6.00, 6.00, 6.60, //brass
470  6.60, 6.60, 6.60, 6.60, 6.60, //brass
471  8.90, 20.65, 19.5}; //brass,iron !
472 
473  for (int ilayer = 0; ilayer<19; ilayer++) {
474 
475  // Histogram mean and sigma calculated from the ROOT histos
476  edm::LogVerbatim("HcalTBSim") << "Layer number: " << ilayer << " Mean = "
477  << histo->getMean(ilayer) << " sigma = "
478  << histo->getRMS(ilayer) << " LThick= "
479  << w[ilayer] << " SThick= " << st[ilayer];
480 
481  oFile << ilayer << " " << histo->getMean(ilayer) << " "
482  << histo->getRMS(ilayer) << " " << w[ilayer] << " " << st[ilayer]
483  << std::endl;
484 
485  }
486  oFile.close();
487  */
488 }
489 
void set_xE3x3N(float v)
Log< level::Info, true > LogVerbatim
void set_E5x5N(float v)
void set_Ntimesli(float v)
int getTrackID() const
Definition: CaloG4Hit.h:64
void set_NUnit(float v)
#define DEFINE_SIMWATCHER(type)
void set_Eentry(float v)
std::map< int, float > energyInScints
void set_xEentry(float v)
void set_xETotN(float v)
void set_Nprim(float v)
void set_E7x7(float v)
~HcalTB02Analysis() override
double getEM() const
Definition: CaloG4Hit.h:55
void set_ETotN(float v)
double npart
Definition: HydjetWrapper.h:48
void set_Einit(float v)
const HcalTB02Analysis & operator=(const HcalTB02Analysis &)=delete
void produce(edm::Event &, const edm::EventSetup &) override
void set_xETot(float v)
void set_E7x7N(float v)
const std::string names[nVars_]
double getHadr() const
Definition: CaloG4Hit.h:58
void fillEvent(HcalTB02HistoClass &)
const edm::ParameterSet m_Anal
void set_xE3x3(float v)
HcalTB02Analysis(const edm::ParameterSet &p)
void update(const BeginOfEvent *evt) override
This routine will be called when the appropriate signal arrives.
void set_xE5x5N(float v)
std::unique_ptr< HcalTB02Histo > histo
T sqrt(T t)
Definition: SSEVec.h:19
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
std::map< int, float > primaries
const std::vector< std::string > names
int getTimeSliceID() const
Definition: CaloG4Hit.h:68
void set_E5x5(float v)
void set_partType(float v)
uint32_t getUnitID() const
Definition: CaloG4Hit.h:66
HLT enums.
G4THitsCollection< CaloG4Hit > CaloG4HitCollection
void set_ETot(float v)
void set_xE5x5(float v)
Log< level::Warning, false > LogWarning
double getIncidentEnergy() const
Definition: CaloG4Hit.h:61
Geom::Theta< T > theta() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511
std::map< int, float > energyInCrystals
void set_xNUnit(float v)