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