CMS 3D CMS Logo

HFShowerParam.cc
Go to the documentation of this file.
1 // File: HFShowerParam.cc
3 // Description: Parametrized version of HF hits
5 
9 
13 
14 #include "G4VPhysicalVolume.hh"
15 #include "G4Step.hh"
16 #include "G4Track.hh"
17 #include "G4NavigationHistory.hh"
18 #include "Randomize.hh"
19 
20 #include "CLHEP/Units/GlobalPhysicalConstants.h"
21 #include "CLHEP/Units/GlobalSystemOfUnits.h"
22 
23 #include <iostream>
24 
25 //#define EDM_ML_DEBUG
26 //#define plotDebug
27 //#define mkdebug
28 
30  const HcalDDDSimConstants* hcons,
31  const HcalSimulationParameters* hps,
32  edm::ParameterSet const& p)
33  : hcalConstants_(hcons), fillHisto_(false) {
34  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
35  edm::ParameterSet m_HF2 = m_HF.getParameter<edm::ParameterSet>("HFShowerBlock");
36  pePerGeV_ = m_HF.getParameter<double>("PEPerGeV");
37  trackEM_ = m_HF.getParameter<bool>("TrackEM");
38  bool useShowerLibrary = m_HF.getParameter<bool>("UseShowerLibrary");
39  bool useGflash = m_HF.getParameter<bool>("UseHFGflash");
40  edMin_ = m_HF.getParameter<double>("EminLibrary");
41  onlyLong_ = m_HF2.getParameter<bool>("OnlyLong");
42  ref_index_ = m_HF.getParameter<double>("RefIndex");
43  double lambdaMean = m_HF.getParameter<double>("LambdaMean");
44  aperture_ = cos(asin(m_HF.getParameter<double>("Aperture")));
45  applyFidCut_ = m_HF.getParameter<bool>("ApplyFiducialCut");
46  parametrizeLast_ = m_HF.getUntrackedParameter<bool>("ParametrizeLast", false);
48 
49  edm::LogVerbatim("HFShower") << "HFShowerParam::Use of shower library is set to " << useShowerLibrary
50  << " Use of Gflash is set to " << useGflash << " P.E. per GeV " << pePerGeV_
51  << ", ref. index of fibre " << ref_index_ << ", Track EM Flag " << trackEM_ << ", edMin "
52  << edMin_ << " GeV, use of Short fibre info in"
53  << " shower library set to " << !(onlyLong_)
54  << ", use of parametrization for last part set to " << parametrizeLast_
55  << ", Mean lambda " << lambdaMean << ", aperture (cutoff) " << aperture_
56  << ", Application of Fiducial Cut " << applyFidCut_;
57 
58 #ifdef plotDebug
60  if (tfile.isAvailable()) {
61  fillHisto_ = true;
62  edm::LogVerbatim("HFShower") << "HFShowerParam::Save histos in directory "
63  << "ProfileFromParam";
64  TFileDirectory showerDir = tfile->mkdir("ProfileFromParam");
65  hzvem_ = showerDir.make<TH1F>("hzvem", "Longitudinal Profile (EM Part);Number of PE", 330, 0.0, 1650.0);
66  hzvhad_ = showerDir.make<TH1F>("hzvhad", "Longitudinal Profile (Had Part);Number of PE", 330, 0.0, 1650.0);
67  em_2d_1_ = showerDir.make<TH2F>(
68  "em_2d_1", "Lateral Profile vs. Shower Depth;cm;Events", 800, 800.0, 1600.0, 100, 50.0, 150.0);
69  em_long_1_ =
70  showerDir.make<TH1F>("em_long_1", "Longitudinal Profile;Radiation Length;Number of Spots", 800, 800.0, 1600.0);
71  em_long_1_tuned_ = showerDir.make<TH1F>(
72  "em_long_1_tuned", "Longitudinal Profile;Radiation Length;Number of Spots", 800, 800.0, 1600.0);
73  em_lateral_1_ = showerDir.make<TH1F>("em_lateral_1", "Lateral Profile;cm;Events", 100, 50.0, 150.0);
74  em_2d_2_ = showerDir.make<TH2F>(
75  "em_2d_2", "Lateral Profile vs. Shower Depth;cm;Events", 800, 800.0, 1600.0, 100, 50.0, 150.0);
76  em_long_2_ =
77  showerDir.make<TH1F>("em_long_2", "Longitudinal Profile;Radiation Length;Number of Spots", 800, 800.0, 1600.0);
78  em_lateral_2_ = showerDir.make<TH1F>("em_lateral_2", "Lateral Profile;cm;Events", 100, 50.0, 150.0);
79  em_long_gflash_ = showerDir.make<TH1F>(
80  "em_long_gflash", "Longitudinal Profile From GFlash;cm;Number of Spots", 800, 800.0, 1600.0);
81  em_long_sl_ = showerDir.make<TH1F>(
82  "em_long_sl", "Longitudinal Profile From Shower Library;cm;Number of Spots", 800, 800.0, 1600.0);
83  } else {
84  fillHisto_ = false;
85  edm::LogVerbatim("HFShower") << "HFShowerParam::No file is available for saving histos so the "
86  << "flag is set to false";
87  }
88 #endif
89 
90  if (useShowerLibrary)
91  showerLibrary_ = std::make_unique<HFShowerLibrary>(name, hcalConstants_, hps, p);
92  else
93  showerLibrary_.reset(nullptr);
94  if (useGflash)
95  gflash_ = std::make_unique<HFGflash>(p);
96  else
97  gflash_.reset(nullptr);
98  fibre_ = std::make_unique<HFFibre>(name, hcalConstants_, hps, p);
99  attLMeanInv_ = fibre_->attLength(lambdaMean);
100  edm::LogVerbatim("HFShower") << "att. length used for (lambda=" << lambdaMean
101  << ") = " << 1 / (attLMeanInv_ * CLHEP::cm) << " cm";
102 }
103 
105 
106 std::vector<HFShowerParam::Hit> HFShowerParam::getHits(const G4Step* aStep, double weight, bool& isKilled) {
107  auto const preStepPoint = aStep->GetPreStepPoint();
108  auto const track = aStep->GetTrack();
110  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
111  double zv = std::abs(hitPoint.z()) - gpar_[4] - 0.5 * gpar_[1];
112  G4ThreeVector localPoint = G4ThreeVector(hitPoint.x(), hitPoint.y(), zv);
113 
114  double pin = (preStepPoint->GetTotalEnergy()) / CLHEP::GeV;
115  double zint = hitPoint.z();
116  double zz = std::abs(zint) - gpar_[4];
117 
118 #ifdef EDM_ML_DEBUG
119  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits " << track->GetDefinition()->GetParticleName()
120  << " of energy " << pin << " GeV Pos x,y,z = " << hitPoint.x() << "," << hitPoint.y()
121  << "," << zint << " (" << zz << "," << localPoint.z() << ", "
122  << (localPoint.z() + 0.5 * gpar_[1]) << ") Local " << localPoint;
123 #endif
124  std::vector<HFShowerParam::Hit> hits;
126  hit.position = hitPoint;
127 
128  // look for other charged particles
129  bool other = false;
130  double pBeta = track->GetDynamicParticle()->GetTotalMomentum() / track->GetDynamicParticle()->GetTotalEnergy();
131  double dirz = (track->GetDynamicParticle()->GetMomentumDirection()).z();
132  if (hitPoint.z() < 0)
133  dirz *= -1.;
134 #ifdef EDM_ML_DEBUG
135  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits Momentum "
136  << track->GetDynamicParticle()->GetMomentumDirection() << " HitPoint " << hitPoint
137  << " dirz " << dirz;
138 #endif
139  if (!isEM && track->GetDefinition()->GetPDGCharge() != 0 && pBeta > (1 / ref_index_) &&
140  aStep->GetTotalEnergyDeposit() > 0.) {
141  other = true;
142  }
143 
144  // take only e+-/gamma/or special particles
145  if (isEM || other) {
146  // Leave out the last part
147  double edep = 0.;
148  if ((!trackEM_) && ((zz < (gpar_[1] - gpar_[2])) || parametrizeLast_) && (!other)) {
149  edep = pin;
150  isKilled = true;
151  } else if ((track->GetDefinition()->GetPDGCharge() != 0) && (pBeta > (1 / ref_index_)) && (dirz > aperture_)) {
152  edep = (aStep->GetTotalEnergyDeposit()) / GeV;
153  }
154  std::string path = "ShowerLibrary";
155 #ifdef EDM_ML_DEBUG
156  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits edep = " << edep << " weight " << weight << " final "
157  << edep * weight << ", Kill = " << isKilled << ", pin = " << pin
158  << ", edMin = " << edMin_ << " Other " << other;
159 #endif
160  edep *= weight;
161  if (edep > 0) {
162  if ((showerLibrary_.get() || gflash_.get()) && isKilled && pin > edMin_ && (!other)) {
163  if (showerLibrary_.get()) {
164  std::vector<HFShowerLibrary::Hit> hitSL = showerLibrary_->getHits(aStep, isKilled, weight, onlyLong_);
165  for (unsigned int i = 0; i < hitSL.size(); i++) {
166  bool ok = true;
167 #ifdef EDM_ML_DEBUG
168  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits applyFidCut = " << applyFidCut_;
169 #endif
170  if (applyFidCut_) { // @@ For showerlibrary no z-cut for Short (no z)
171  int npmt = HFFibreFiducial::PMTNumber(hitSL[i].position);
172  if (npmt <= 0)
173  ok = false;
174  }
175  if (ok) {
176  hit.position = hitSL[i].position;
177  hit.depth = hitSL[i].depth;
178  hit.time = hitSL[i].time;
179  hit.edep = 1;
180  hits.push_back(hit);
181 #ifdef plotDebug
182  if (fillHisto_) {
183  double zv = std::abs(hit.position.z()) - gpar_[4];
184  hzvem_->Fill(zv);
185  em_long_sl_->Fill(hit.position.z() / CLHEP::cm);
186  double sq = sqrt(pow(hit.position.x() / CLHEP::cm, 2) + pow(hit.position.y() / CLHEP::cm, 2));
187  double zp = hit.position.z() / CLHEP::cm;
188  if (hit.depth == 1) {
189  em_2d_1_->Fill(zp, sq);
190  em_lateral_1_->Fill(sq);
191  em_long_1_->Fill(zp);
192  } else if (hit.depth == 2) {
193  em_2d_2_->Fill(zp, sq);
194  em_lateral_2_->Fill(sq);
195  em_long_2_->Fill(zp);
196  }
197  }
198 #endif
199 #ifdef EDM_ML_DEBUG
200  edm::LogVerbatim("HFShower")
201  << "HFShowerParam: Hit at depth " << hit.depth << " with edep " << hit.edep << " Time " << hit.time;
202 #endif
203  }
204  }
205  } else { // GFlash clusters with known z
206  std::vector<HFGflash::Hit> hitSL = gflash_->gfParameterization(aStep, onlyLong_);
207  for (unsigned int i = 0; i < hitSL.size(); ++i) {
208  bool ok = true;
209  G4ThreeVector pe_effect(hitSL[i].position.x(), hitSL[i].position.y(), hitSL[i].position.z());
210  double zv = std::abs(pe_effect.z()) - gpar_[4];
211  //depth
212  int depth = 1;
213  int npmt = 0;
214  if (zv < 0. || zv > gpar_[1]) {
215 #ifdef mkdebug
216  edm::LogVerbatim("HFShower") << "-#Zcut-HFShowerParam::getHits:z=" << zv << ",m=" << gpar_[1];
217 #endif
218  ok = false;
219  }
220  if (ok && applyFidCut_) {
221  npmt = HFFibreFiducial::PMTNumber(pe_effect);
222 #ifdef EDM_ML_DEBUG
223  edm::LogVerbatim("HFShower") << "HFShowerParam::getHits:#PMT= " << npmt << ",z = " << zv;
224 #endif
225  if (npmt <= 0) {
226 #ifdef EDM_ML_DEBUG
227  edm::LogVerbatim("HFShower") << "-#PMT=0 cut-HFShowerParam::getHits: npmt = " << npmt;
228 #endif
229  ok = false;
230  } else if (npmt > 24) { // a short fibre
231  if (zv > gpar_[0]) {
232  depth = 2;
233  } else {
234 #ifdef EDM_ML_DEBUG
235  edm::LogVerbatim("HFShower") << "-SHORT cut-HFShowerParam::getHits:zMin=" << gpar_[0];
236 #endif
237  ok = false;
238  }
239  }
240 #ifdef EDM_ML_DEBUG
241  edm::LogVerbatim("HFShower")
242  << "HFShowerParam: npmt " << npmt << " zv " << std::abs(pe_effect.z()) << ":" << gpar_[4] << ":" << zv
243  << ":" << gpar_[0] << " ok " << ok << " depth " << depth;
244 #endif
245  } else {
246  if (G4UniformRand() > 0.5)
247  depth = 2;
248  if (depth == 2 && zv < gpar_[0])
249  ok = false;
250  }
251  //attenuation
252  double dist = fibre_->zShift(localPoint, depth, 0); // distance to PMT
253  double r1 = G4UniformRand();
254 #ifdef EDM_ML_DEBUG
255  edm::LogVerbatim("HFShower") << "HFShowerParam:Distance to PMT (" << npmt << ") " << dist
256  << ", exclusion flag " << (r1 > exp(-attLMeanInv_ * zv));
257 #endif
258  if (r1 > exp(-attLMeanInv_ * dist))
259  ok = false;
260  if (ok) {
261  double r2 = G4UniformRand();
262 #ifdef EDM_ML_DEBUG
263  edm::LogVerbatim("HFShower")
264  << "HFShowerParam:Extra exclusion " << r2 << ">" << weight << " " << (r2 > weight);
265 #endif
266  if (r2 < weight) {
267  double time = fibre_->tShift(localPoint, depth, 0);
268 
269  hit.position = hitSL[i].position;
270  hit.depth = depth;
271  hit.time = time + hitSL[i].time;
272  hit.edep = 1;
273  hits.push_back(hit);
274 #ifdef plotDebug
275  if (fillHisto_) {
276  em_long_gflash_->Fill(pe_effect.z() / CLHEP::cm, hitSL[i].edep);
277  hzvem_->Fill(zv);
278  double sq = sqrt(pow(hit.position.x() / CLHEP::cm, 2) + pow(hit.position.y() / CLHEP::cm, 2));
279  double zp = hit.position.z() / CLHEP::cm;
280  if (hit.depth == 1) {
281  em_2d_1_->Fill(zp, sq);
282  em_lateral_1_->Fill(s);
283  em_long_1_->Fill(zp);
284  } else if (hit.depth == 2) {
285  em_2d_2_->Fill(zp, sq);
286  em_lateral_2_->Fill(sq);
287  em_long_2_->Fill(zp);
288  }
289  }
290 #endif
291 #ifdef EDM_ML_DEBUG
292  edm::LogVerbatim("HFShower")
293  << "HFShowerParam: Hit at depth " << hit.depth << " with edep " << hit.edep << " Time " << hit.time;
294 #endif
295  }
296  }
297  }
298  }
299  } else {
300  path = "Rest";
301  edep *= pePerGeV_;
302  double tSlice = (aStep->GetPostStepPoint()->GetGlobalTime());
303  double time = fibre_->tShift(localPoint, 1, 0); // remaining part
304  bool ok = true;
305  if (applyFidCut_) { // @@ For showerlibrary no z-cut for Short (no z)
306  int npmt = HFFibreFiducial::PMTNumber(hitPoint);
307  if (npmt <= 0)
308  ok = false;
309  }
310 #ifdef EDM_ML_DEBUG
311  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits hitPoint " << hitPoint << " flag " << ok;
312 #endif
313  if (ok) {
314  hit.depth = 1;
315  hit.time = tSlice + time;
316  hit.edep = edep;
317  hits.push_back(hit);
318 #ifdef EDM_ML_DEBUG
319  edm::LogVerbatim("HFShower") << "HFShowerParam: Hit at depth 1 with edep " << edep << " Time " << tSlice
320  << ":" << time << ":" << hit.time;
321 #endif
322 #ifdef plotDebug
323  double zv = std::abs(hitPoint.z()) - gpar_[4];
324  if (fillHisto_) {
325  hzvhad_->Fill(zv);
326  }
327 #endif
328  if (zz >= gpar_[0]) {
329  time = fibre_->tShift(localPoint, 2, 0);
330  hit.depth = 2;
331  hit.time = tSlice + time;
332  hits.push_back(hit);
333 #ifdef EDM_ML_DEBUG
334  edm::LogVerbatim("HFShower") << "HFShowerParam: Hit at depth 2 with edep " << edep << " Time " << tSlice
335  << ":" << time << hit.time;
336 #endif
337 #ifdef plotDebug
338  if (fillHisto_) {
339  hzvhad_->Fill(zv);
340  }
341 #endif
342  }
343  }
344  }
345 #ifdef EDM_ML_DEBUG
346  for (unsigned int ii = 0; ii < hits.size(); ++ii) {
347  double zv = std::abs(hits[ii].position.z());
348  if (zv > 12790)
349  edm::LogVerbatim("HFShower") << "HFShowerParam: Abnormal hit along " << path << " in "
350  << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName() << " at "
351  << hits[ii].position << " zz " << zv << " Edep " << edep << " due to "
352  << track->GetDefinition()->GetParticleName() << " time " << hit.time;
353  }
354  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits kill (" << isKilled << ") track " << track->GetTrackID()
355  << " at " << hitPoint << " and deposit " << edep << " " << hits.size() << " times"
356  << " ZZ " << zz << " " << gpar_[0];
357 #endif
358  }
359  }
360  return hits;
361 }
HFShowerParam::onlyLong_
bool onlyLong_
Definition: HFShowerParam.h:50
mps_fire.i
i
Definition: mps_fire.py:428
geometryCSVtoXML.zz
zz
Definition: geometryCSVtoXML.py:19
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11713
HFShowerParam::getHits
std::vector< Hit > getHits(const G4Step *aStep, double weight, bool &isKilled)
Definition: HFShowerParam.cc:106
funct::false
false
Definition: Factorize.h:29
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
Calorimetry_cff.useShowerLibrary
useShowerLibrary
Definition: Calorimetry_cff.py:269
hit::y
double y
Definition: SiStripHitEffFromCalibTree.cc:90
TFileDirectory::make
T * make(const Args &... args) const
make new ROOT object
Definition: TFileDirectory.h:53
HFFibreFiducial.h
mps_merge.weight
weight
Definition: mps_merge.py:88
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
HFShowerParam.h
protons_cff.time
time
Definition: protons_cff.py:39
HcalDDDSimConstants::getGparHF
const std::vector< double > & getGparHF() const
Definition: HcalDDDSimConstants.h:45
HFShowerParam::em_long_gflash_
TH1F * em_long_gflash_
Definition: HFShowerParam.h:54
gpuVertexFinder::zv
float *__restrict__ zv
Definition: gpuFitVertices.h:26
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
HFShowerParam::~HFShowerParam
virtual ~HFShowerParam()
Definition: HFShowerParam.cc:104
TFileDirectory
Definition: TFileDirectory.h:24
HFShowerParam::fibre_
std::unique_ptr< HFFibre > fibre_
Definition: HFShowerParam.h:46
HFShowerParam::hzvhad_
TH1F * hzvhad_
Definition: HFShowerParam.h:54
convertSQLiteXML.ok
bool ok
Definition: convertSQLiteXML.py:98
HFFibreFiducial::PMTNumber
int PMTNumber(const G4ThreeVector &pe_effect)
Definition: HFFibreFiducial.cc:9
G4TrackToParticleID::isGammaElectronPositron
static bool isGammaElectronPositron(int pdgCode)
Definition: G4TrackToParticleID.cc:17
hit::x
double x
Definition: SiStripHitEffFromCalibTree.cc:89
HcalDDDSimConstants
Definition: HcalDDDSimConstants.h:24
HFShowerParam::hzvem_
TH1F * hzvem_
Definition: HFShowerParam.h:54
HFShowerParam::HFShowerParam
HFShowerParam(const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
Definition: HFShowerParam.cc:29
alignCSCRings.s
s
Definition: alignCSCRings.py:92
HFShowerParam::gpar_
std::vector< double > gpar_
Definition: HFShowerParam.h:52
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
HFShowerParam::applyFidCut_
bool applyFidCut_
Definition: HFShowerParam.h:50
HFShowerParam::em_long_1_
TH1F * em_long_1_
Definition: HFShowerParam.h:53
Service.h
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
HFShowerParam::showerLibrary_
std::unique_ptr< HFShowerLibrary > showerLibrary_
Definition: HFShowerParam.h:45
tfile
Definition: tfile.py:1
DDAxes::z
trackingPlots.other
other
Definition: trackingPlots.py:1460
LEDCalibrationChannels.depth
depth
Definition: LEDCalibrationChannels.py:65
hit::z
double z
Definition: SiStripHitEffFromCalibTree.cc:91
HFShowerParam::em_lateral_2_
TH1F * em_lateral_2_
Definition: HFShowerParam.h:53
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HFShowerParam::em_lateral_1_
TH1F * em_lateral_1_
Definition: HFShowerParam.h:53
TFileService.h
edm::ParameterSet
Definition: ParameterSet.h:47
HFShowerParam::em_long_1_tuned_
TH1F * em_long_1_tuned_
Definition: HFShowerParam.h:54
HFShowerParam::trackEM_
bool trackEM_
Definition: HFShowerParam.h:50
HFShowerParam::em_long_sl_
TH1F * em_long_sl_
Definition: HFShowerParam.h:55
GeV
const double GeV
Definition: MathUtil.h:16
diffTwoXMLs.r2
r2
Definition: diffTwoXMLs.py:73
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
HFShowerParam::hcalConstants_
const HcalDDDSimConstants * hcalConstants_
Definition: HFShowerParam.h:44
edm::Service< TFileService >
HFShowerParam::parametrizeLast_
bool parametrizeLast_
Definition: HFShowerParam.h:50
HFShowerParam::edMin_
double edMin_
Definition: HFShowerParam.h:49
HFShowerParam::gflash_
std::unique_ptr< HFGflash > gflash_
Definition: HFShowerParam.h:47
compare.tfile
tfile
Definition: compare.py:325
HFShowerParam::ref_index_
double ref_index_
Definition: HFShowerParam.h:49
HFShowerParam::em_long_2_
TH1F * em_long_2_
Definition: HFShowerParam.h:53
G4TrackToParticleID.h
diffTwoXMLs.r1
r1
Definition: diffTwoXMLs.py:53
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
HFShowerParam::fillHisto_
bool fillHisto_
Definition: HFShowerParam.h:48
HFShowerParam::pePerGeV_
double pePerGeV_
Definition: HFShowerParam.h:49
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HFShowerParam::attLMeanInv_
double attLMeanInv_
Definition: HFShowerParam.h:49
Exception.h
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
castor_dqm_sourceclient_file_cfg.path
path
Definition: castor_dqm_sourceclient_file_cfg.py:37
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
HFShowerParam::Hit
Definition: HFShowerParam.h:34
HcalSimulationParameters
Definition: HcalSimulationParameters.h:6
cuy.ii
ii
Definition: cuy.py:590
HFShowerParam::aperture_
double aperture_
Definition: HFShowerParam.h:49
weight
Definition: weight.py:1
hit
Definition: SiStripHitEffFromCalibTree.cc:88
HFShowerParam::em_2d_1_
TH2F * em_2d_1_
Definition: HFShowerParam.h:56
HFShowerParam::em_2d_2_
TH2F * em_2d_2_
Definition: HFShowerParam.h:56