CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Classes | Public Member Functions | Private Attributes
HFShowerParam Class Reference

#include <HFShowerParam.h>

Classes

struct  Hit
 

Public Member Functions

std::vector< HitgetHits (const G4Step *aStep, double weight, bool &isKilled)
 
 HFShowerParam (const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
 
virtual ~HFShowerParam ()
 

Private Attributes

double aperture_
 
bool applyFidCut_
 
double attLMeanInv_
 
double edMin_
 
TH2F * em_2d_1_
 
TH2F * em_2d_2_
 
TH1F * em_lateral_1_
 
TH1F * em_lateral_2_
 
TH1F * em_long_1_
 
TH1F * em_long_1_tuned_
 
TH1F * em_long_2_
 
TH1F * em_long_gflash_
 
TH1F * em_long_sl_
 
G4int emPDG_
 
G4int epPDG_
 
bool equalizeTimeShift_
 
std::unique_ptr< HFFibrefibre_
 
bool fillHisto_
 
G4int gammaPDG_
 
std::unique_ptr< HFGflashgflash_
 
std::vector< double > gpar_
 
const HcalDDDSimConstantshcalConstants_
 
TH1F * hzvem_
 
TH1F * hzvhad_
 
bool onlyLong_
 
bool parametrizeLast_
 
double pePerGeV_
 
double ref_index_
 
std::unique_ptr< HFShowerLibraryshowerLibrary_
 
bool trackEM_
 

Detailed Description

Definition at line 25 of file HFShowerParam.h.

Constructor & Destructor Documentation

HFShowerParam::HFShowerParam ( const std::string &  name,
const HcalDDDSimConstants hcons,
const HcalSimulationParameters hps,
edm::ParameterSet const &  p 
)

Definition at line 29 of file HFShowerParam.cc.

References aperture_, applyFidCut_, attLMeanInv_, funct::cos(), edMin_, em_2d_1_, em_2d_2_, em_lateral_1_, em_lateral_2_, em_long_1_, em_long_1_tuned_, em_long_2_, em_long_gflash_, em_long_sl_, equalizeTimeShift_, fibre_, fillHisto_, HcalDDDSimConstants::getGparHF(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), gflash_, gpar_, hcalConstants_, hzvem_, hzvhad_, edm::Service< T >::isAvailable(), TFileDirectory::make(), TFileService::mkdir(), mergeVDriftHistosByStation::name, onlyLong_, AlCaHLTBitMon_ParallelJobs::p, parametrizeLast_, pePerGeV_, ref_index_, showerLibrary_, compare::tfile, and trackEM_.

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  equalizeTimeShift_ = m_HF2.getParameter<bool>("EqualizeTimeShift");
42  onlyLong_ = m_HF2.getParameter<bool>("OnlyLong");
43  ref_index_ = m_HF.getParameter<double>("RefIndex");
44  double lambdaMean = m_HF.getParameter<double>("LambdaMean");
45  aperture_ = cos(asin(m_HF.getParameter<double>("Aperture")));
46  applyFidCut_ = m_HF.getParameter<bool>("ApplyFiducialCut");
47  parametrizeLast_ = m_HF.getUntrackedParameter<bool>("ParametrizeLast", false);
49 
50  edm::LogVerbatim("HFShower") << "HFShowerParam::Use of shower library is set to " << useShowerLibrary
51  << " Use of Gflash is set to " << useGflash << " P.E. per GeV " << pePerGeV_
52  << ", ref. index of fibre " << ref_index_ << ", Track EM Flag " << trackEM_ << ", edMin "
53  << edMin_ << " GeV, use of Short fibre info in"
54  << " shower library set to " << !(onlyLong_)
55  << " equalize flag for time shift in fibre is set to " << equalizeTimeShift_
56  << ", use of parametrization for last part set to " << parametrizeLast_
57  << ", Mean lambda " << lambdaMean << ", aperture (cutoff) " << aperture_
58  << ", Application of Fiducial Cut " << applyFidCut_;
59 
60 #ifdef plotDebug
62  if (tfile.isAvailable()) {
63  fillHisto_ = true;
64  edm::LogVerbatim("HFShower") << "HFShowerParam::Save histos in directory "
65  << "ProfileFromParam";
66  TFileDirectory showerDir = tfile->mkdir("ProfileFromParam");
67  hzvem_ = showerDir.make<TH1F>("hzvem", "Longitudinal Profile (EM Part);Number of PE", 330, 0.0, 1650.0);
68  hzvhad_ = showerDir.make<TH1F>("hzvhad", "Longitudinal Profile (Had Part);Number of PE", 330, 0.0, 1650.0);
69  em_2d_1_ = showerDir.make<TH2F>(
70  "em_2d_1", "Lateral Profile vs. Shower Depth;cm;Events", 800, 800.0, 1600.0, 100, 50.0, 150.0);
71  em_long_1_ =
72  showerDir.make<TH1F>("em_long_1", "Longitudinal Profile;Radiation Length;Number of Spots", 800, 800.0, 1600.0);
73  em_long_1_tuned_ = showerDir.make<TH1F>(
74  "em_long_1_tuned", "Longitudinal Profile;Radiation Length;Number of Spots", 800, 800.0, 1600.0);
75  em_lateral_1_ = showerDir.make<TH1F>("em_lateral_1", "Lateral Profile;cm;Events", 100, 50.0, 150.0);
76  em_2d_2_ = showerDir.make<TH2F>(
77  "em_2d_2", "Lateral Profile vs. Shower Depth;cm;Events", 800, 800.0, 1600.0, 100, 50.0, 150.0);
78  em_long_2_ =
79  showerDir.make<TH1F>("em_long_2", "Longitudinal Profile;Radiation Length;Number of Spots", 800, 800.0, 1600.0);
80  em_lateral_2_ = showerDir.make<TH1F>("em_lateral_2", "Lateral Profile;cm;Events", 100, 50.0, 150.0);
81  em_long_gflash_ = showerDir.make<TH1F>(
82  "em_long_gflash", "Longitudinal Profile From GFlash;cm;Number of Spots", 800, 800.0, 1600.0);
83  em_long_sl_ = showerDir.make<TH1F>(
84  "em_long_sl", "Longitudinal Profile From Shower Library;cm;Number of Spots", 800, 800.0, 1600.0);
85  } else {
86  fillHisto_ = false;
87  edm::LogVerbatim("HFShower") << "HFShowerParam::No file is available for saving histos so the "
88  << "flag is set to false";
89  }
90 #endif
91 
92  if (useShowerLibrary)
93  showerLibrary_ = std::make_unique<HFShowerLibrary>(name, hcalConstants_, hps, p);
94  else
95  showerLibrary_.reset(nullptr);
96  if (useGflash)
97  gflash_ = std::make_unique<HFGflash>(p);
98  else
99  gflash_.reset(nullptr);
100  fibre_ = std::make_unique<HFFibre>(name, hcalConstants_, hps, p);
101  attLMeanInv_ = fibre_->attLength(lambdaMean);
102  edm::LogVerbatim("HFShower") << "att. length used for (lambda=" << lambdaMean
103  << ") = " << 1 / (attLMeanInv_ * CLHEP::cm) << " cm";
104 }
std::unique_ptr< HFGflash > gflash_
Definition: HFShowerParam.h:47
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
double aperture_
Definition: HFShowerParam.h:49
TH1F * em_long_sl_
Definition: HFShowerParam.h:55
double attLMeanInv_
Definition: HFShowerParam.h:49
double pePerGeV_
Definition: HFShowerParam.h:49
TH1F * em_long_2_
Definition: HFShowerParam.h:53
TH1F * em_long_1_tuned_
Definition: HFShowerParam.h:54
TH1F * em_long_1_
Definition: HFShowerParam.h:53
TH1F * em_long_gflash_
Definition: HFShowerParam.h:54
TH1F * em_lateral_2_
Definition: HFShowerParam.h:53
std::vector< double > gpar_
Definition: HFShowerParam.h:52
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
bool isAvailable() const
Definition: Service.h:40
std::unique_ptr< HFFibre > fibre_
Definition: HFShowerParam.h:46
T * make(const Args &...args) const
make new ROOT object
TH1F * em_lateral_1_
Definition: HFShowerParam.h:53
const std::vector< double > & getGparHF() const
bool equalizeTimeShift_
Definition: HFShowerParam.h:50
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
bool parametrizeLast_
Definition: HFShowerParam.h:50
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::unique_ptr< HFShowerLibrary > showerLibrary_
Definition: HFShowerParam.h:45
double ref_index_
Definition: HFShowerParam.h:49
const HcalDDDSimConstants * hcalConstants_
Definition: HFShowerParam.h:44
tuple tfile
Definition: compare.py:324
HFShowerParam::~HFShowerParam ( )
virtual

Definition at line 106 of file HFShowerParam.cc.

106 {}

Member Function Documentation

std::vector< HFShowerParam::Hit > HFShowerParam::getHits ( const G4Step *  aStep,
double  weight,
bool &  isKilled 
)

Definition at line 108 of file HFShowerParam.cc.

References funct::abs(), aperture_, applyFidCut_, attLMeanInv_, HFShowerParam::Hit::depth, HLT_FULL_cff::depth, HFShowerParam::Hit::edep, edMin_, em_2d_1_, em_2d_2_, em_lateral_1_, em_lateral_2_, em_long_1_, em_long_2_, em_long_gflash_, em_long_sl_, equalizeTimeShift_, funct::exp(), fibre_, fillHisto_, GeV, gflash_, gpar_, hzvem_, hzvhad_, mps_fire::i, cuy::ii, G4TrackToParticleID::isGammaElectronPositron(), convertSQLiteXML::ok, onlyLong_, parametrizeLast_, fed_dqm_sourceclient-live_cfg::path, pePerGeV_, HFFibreFiducial::PMTNumber(), HFShowerParam::Hit::position, position, funct::pow(), diffTwoXMLs::r1, diffTwoXMLs::r2, ref_index_, alignCSCRings::s, showerLibrary_, mathSSE::sqrt(), AlCaHLTBitMon_QueryRunRegistry::string, HFShowerParam::Hit::time, HLT_FULL_cff::track, trackEM_, histoStyle::weight, z, and gpuVertexFinder::zv.

108  {
109  auto const preStepPoint = aStep->GetPreStepPoint();
110  auto const track = aStep->GetTrack();
112  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
113  double zv = std::abs(hitPoint.z()) - gpar_[4] - 0.5 * gpar_[1];
114  G4ThreeVector localPoint = G4ThreeVector(hitPoint.x(), hitPoint.y(), zv);
115 
116  double pin = (preStepPoint->GetTotalEnergy()) / CLHEP::GeV;
117  double zint = hitPoint.z();
118  double zz = std::abs(zint) - gpar_[4];
119 
120 #ifdef EDM_ML_DEBUG
121  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits " << track->GetDefinition()->GetParticleName()
122  << " of energy " << pin << " GeV Pos x,y,z = " << hitPoint.x() << "," << hitPoint.y()
123  << "," << zint << " (" << zz << "," << localPoint.z() << ", "
124  << (localPoint.z() + 0.5 * gpar_[1]) << ") Local " << localPoint;
125 #endif
126  std::vector<HFShowerParam::Hit> hits;
128  hit.position = hitPoint;
129 
130  // look for other charged particles
131  bool other = false;
132  double pBeta = track->GetDynamicParticle()->GetTotalMomentum() / track->GetDynamicParticle()->GetTotalEnergy();
133  double dirz = (track->GetDynamicParticle()->GetMomentumDirection()).z();
134  if (hitPoint.z() < 0)
135  dirz *= -1.;
136 #ifdef EDM_ML_DEBUG
137  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits Momentum "
138  << track->GetDynamicParticle()->GetMomentumDirection() << " HitPoint " << hitPoint
139  << " dirz " << dirz;
140 #endif
141  if (!isEM && track->GetDefinition()->GetPDGCharge() != 0 && pBeta > (1 / ref_index_) &&
142  aStep->GetTotalEnergyDeposit() > 0.) {
143  other = true;
144  }
145 
146  // take only e+-/gamma/or special particles
147  if (isEM || other) {
148  // Leave out the last part
149  double edep = 0.;
150  if ((!trackEM_) && ((zz < (gpar_[1] - gpar_[2])) || parametrizeLast_) && (!other)) {
151  edep = pin;
152  isKilled = true;
153  } else if ((track->GetDefinition()->GetPDGCharge() != 0) && (pBeta > (1 / ref_index_)) && (dirz > aperture_)) {
154  edep = (aStep->GetTotalEnergyDeposit()) / GeV;
155  }
156  std::string path = "ShowerLibrary";
157 #ifdef EDM_ML_DEBUG
158  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits edep = " << edep << " weight " << weight << " final "
159  << edep * weight << ", Kill = " << isKilled << ", pin = " << pin
160  << ", edMin = " << edMin_ << " Other " << other;
161 #endif
162  edep *= weight;
163  if (edep > 0) {
164  if ((showerLibrary_.get() || gflash_.get()) && isKilled && pin > edMin_ && (!other)) {
165  if (showerLibrary_.get()) {
166  std::vector<HFShowerLibrary::Hit> hitSL = showerLibrary_->getHits(aStep, isKilled, weight, onlyLong_);
167  for (unsigned int i = 0; i < hitSL.size(); i++) {
168  bool ok = true;
169 #ifdef EDM_ML_DEBUG
170  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits applyFidCut = " << applyFidCut_;
171 #endif
172  if (applyFidCut_) { // @@ For showerlibrary no z-cut for Short (no z)
173  int npmt = HFFibreFiducial::PMTNumber(hitSL[i].position);
174  if (npmt <= 0)
175  ok = false;
176  }
177  if (ok) {
178  hit.position = hitSL[i].position;
179  hit.depth = hitSL[i].depth;
180  hit.time = hitSL[i].time;
181  hit.edep = 1;
182  hits.push_back(hit);
183 #ifdef plotDebug
184  if (fillHisto_) {
185  double zv = std::abs(hit.position.z()) - gpar_[4];
186  hzvem_->Fill(zv);
187  em_long_sl_->Fill(hit.position.z() / CLHEP::cm);
188  double sq = sqrt(pow(hit.position.x() / CLHEP::cm, 2) + pow(hit.position.y() / CLHEP::cm, 2));
189  double zp = hit.position.z() / CLHEP::cm;
190  if (hit.depth == 1) {
191  em_2d_1_->Fill(zp, sq);
192  em_lateral_1_->Fill(sq);
193  em_long_1_->Fill(zp);
194  } else if (hit.depth == 2) {
195  em_2d_2_->Fill(zp, sq);
196  em_lateral_2_->Fill(sq);
197  em_long_2_->Fill(zp);
198  }
199  }
200 #endif
201 #ifdef EDM_ML_DEBUG
202  edm::LogVerbatim("HFShower")
203  << "HFShowerParam: Hit at depth " << hit.depth << " with edep " << hit.edep << " Time " << hit.time;
204 #endif
205  }
206  }
207  } else { // GFlash clusters with known z
208  std::vector<HFGflash::Hit> hitSL = gflash_->gfParameterization(aStep, onlyLong_);
209  for (unsigned int i = 0; i < hitSL.size(); ++i) {
210  bool ok = true;
211  G4ThreeVector pe_effect(hitSL[i].position.x(), hitSL[i].position.y(), hitSL[i].position.z());
212  double zv = std::abs(pe_effect.z()) - gpar_[4];
213  //depth
214  int depth = 1;
215  int npmt = 0;
216  if (zv < 0. || zv > gpar_[1]) {
217 #ifdef mkdebug
218  edm::LogVerbatim("HFShower") << "-#Zcut-HFShowerParam::getHits:z=" << zv << ",m=" << gpar_[1];
219 #endif
220  ok = false;
221  }
222  if (ok && applyFidCut_) {
223  npmt = HFFibreFiducial::PMTNumber(pe_effect);
224 #ifdef EDM_ML_DEBUG
225  edm::LogVerbatim("HFShower") << "HFShowerParam::getHits:#PMT= " << npmt << ",z = " << zv;
226 #endif
227  if (npmt <= 0) {
228 #ifdef EDM_ML_DEBUG
229  edm::LogVerbatim("HFShower") << "-#PMT=0 cut-HFShowerParam::getHits: npmt = " << npmt;
230 #endif
231  ok = false;
232  } else if (npmt > 24) { // a short fibre
233  if (zv > gpar_[0]) {
234  depth = 2;
235  } else {
236 #ifdef EDM_ML_DEBUG
237  edm::LogVerbatim("HFShower") << "-SHORT cut-HFShowerParam::getHits:zMin=" << gpar_[0];
238 #endif
239  ok = false;
240  }
241  }
242 #ifdef EDM_ML_DEBUG
243  edm::LogVerbatim("HFShower")
244  << "HFShowerParam: npmt " << npmt << " zv " << std::abs(pe_effect.z()) << ":" << gpar_[4] << ":" << zv
245  << ":" << gpar_[0] << " ok " << ok << " depth " << depth;
246 #endif
247  } else {
248  if (G4UniformRand() > 0.5)
249  depth = 2;
250  if (depth == 2 && zv < gpar_[0])
251  ok = false;
252  }
253  //attenuation
254  double dist = fibre_->zShift(localPoint, depth, 0); // distance to PMT
255  double r1 = G4UniformRand();
256 #ifdef EDM_ML_DEBUG
257  edm::LogVerbatim("HFShower") << "HFShowerParam:Distance to PMT (" << npmt << ") " << dist
258  << ", exclusion flag " << (r1 > exp(-attLMeanInv_ * zv));
259 #endif
260  if (r1 > exp(-attLMeanInv_ * dist))
261  ok = false;
262  if (ok) {
263  double r2 = G4UniformRand();
264 #ifdef EDM_ML_DEBUG
265  edm::LogVerbatim("HFShower")
266  << "HFShowerParam:Extra exclusion " << r2 << ">" << weight << " " << (r2 > weight);
267 #endif
268  if (r2 < weight) {
269  double time = (equalizeTimeShift_) ? (fibre_->tShift(localPoint, depth, -1))
270  : (fibre_->tShift(localPoint, depth, 0));
271 
272  hit.position = hitSL[i].position;
273  hit.depth = depth;
274  hit.time = time + hitSL[i].time;
275  hit.edep = 1;
276  hits.push_back(hit);
277 #ifdef plotDebug
278  if (fillHisto_) {
279  em_long_gflash_->Fill(pe_effect.z() / CLHEP::cm, hitSL[i].edep);
280  hzvem_->Fill(zv);
281  double sq = sqrt(pow(hit.position.x() / CLHEP::cm, 2) + pow(hit.position.y() / CLHEP::cm, 2));
282  double zp = hit.position.z() / CLHEP::cm;
283  if (hit.depth == 1) {
284  em_2d_1_->Fill(zp, sq);
285  em_lateral_1_->Fill(s);
286  em_long_1_->Fill(zp);
287  } else if (hit.depth == 2) {
288  em_2d_2_->Fill(zp, sq);
289  em_lateral_2_->Fill(sq);
290  em_long_2_->Fill(zp);
291  }
292  }
293 #endif
294 #ifdef EDM_ML_DEBUG
295  edm::LogVerbatim("HFShower")
296  << "HFShowerParam: Hit at depth " << hit.depth << " with edep " << hit.edep << " Time " << hit.time;
297 #endif
298  }
299  }
300  }
301  }
302  } else {
303  path = "Rest";
304  edep *= pePerGeV_;
305  double tSlice = (aStep->GetPostStepPoint()->GetGlobalTime());
306  double time = (equalizeTimeShift_) ? (fibre_->tShift(localPoint, 1, -1))
307  : (fibre_->tShift(localPoint, 1, 0)); // remaining part
308  bool ok = true;
309  if (applyFidCut_) { // @@ For showerlibrary no z-cut for Short (no z)
310  int npmt = HFFibreFiducial::PMTNumber(hitPoint);
311  if (npmt <= 0)
312  ok = false;
313  }
314 #ifdef EDM_ML_DEBUG
315  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits hitPoint " << hitPoint << " flag " << ok;
316 #endif
317  if (ok) {
318  hit.depth = 1;
319  hit.time = tSlice + time;
320  hit.edep = edep;
321  hits.push_back(hit);
322 #ifdef EDM_ML_DEBUG
323  edm::LogVerbatim("HFShower") << "HFShowerParam: Hit at depth 1 with edep " << edep << " Time " << tSlice
324  << ":" << time << ":" << hit.time;
325 #endif
326 #ifdef plotDebug
327  double zv = std::abs(hitPoint.z()) - gpar_[4];
328  if (fillHisto_) {
329  hzvhad_->Fill(zv);
330  }
331 #endif
332  if (zz >= gpar_[0]) {
333  time = (equalizeTimeShift_) ? (fibre_->tShift(localPoint, 2, -1)) : (fibre_->tShift(localPoint, 2, 0));
334  hit.depth = 2;
335  hit.time = tSlice + time;
336  hits.push_back(hit);
337 #ifdef EDM_ML_DEBUG
338  edm::LogVerbatim("HFShower") << "HFShowerParam: Hit at depth 2 with edep " << edep << " Time " << tSlice
339  << ":" << time << hit.time;
340 #endif
341 #ifdef plotDebug
342  if (fillHisto_) {
343  hzvhad_->Fill(zv);
344  }
345 #endif
346  }
347  }
348  }
349 #ifdef EDM_ML_DEBUG
350  for (unsigned int ii = 0; ii < hits.size(); ++ii) {
351  double zv = std::abs(hits[ii].position.z());
352  if (zv > 12790)
353  edm::LogVerbatim("HFShower") << "HFShowerParam: Abnormal hit along " << path << " in "
354  << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName() << " at "
355  << hits[ii].position << " zz " << zv << " Edep " << edep << " due to "
356  << track->GetDefinition()->GetParticleName() << " time " << hit.time;
357  }
358  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits kill (" << isKilled << ") track " << track->GetTrackID()
359  << " at " << hitPoint << " and deposit " << edep << " " << hits.size() << " times"
360  << " ZZ " << zz << " " << gpar_[0];
361 #endif
362  }
363  }
364  return hits;
365 }
std::unique_ptr< HFGflash > gflash_
Definition: HFShowerParam.h:47
Log< level::Info, true > LogVerbatim
const double GeV
Definition: MathUtil.h:16
double aperture_
Definition: HFShowerParam.h:49
TH1F * em_long_sl_
Definition: HFShowerParam.h:55
double attLMeanInv_
Definition: HFShowerParam.h:49
double pePerGeV_
Definition: HFShowerParam.h:49
float *__restrict__ zv
TH1F * em_long_2_
Definition: HFShowerParam.h:53
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
int ii
Definition: cuy.py:589
TH1F * em_long_1_
Definition: HFShowerParam.h:53
TH1F * em_long_gflash_
Definition: HFShowerParam.h:54
TH1F * em_lateral_2_
Definition: HFShowerParam.h:53
std::vector< double > gpar_
Definition: HFShowerParam.h:52
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::unique_ptr< HFFibre > fibre_
Definition: HFShowerParam.h:46
TH1F * em_lateral_1_
Definition: HFShowerParam.h:53
G4ThreeVector position
Definition: HFShowerParam.h:36
int PMTNumber(const G4ThreeVector &pe_effect)
bool equalizeTimeShift_
Definition: HFShowerParam.h:50
bool parametrizeLast_
Definition: HFShowerParam.h:50
std::unique_ptr< HFShowerLibrary > showerLibrary_
Definition: HFShowerParam.h:45
static int position[264][3]
Definition: ReadPGInfo.cc:289
static bool isGammaElectronPositron(int pdgCode)
double ref_index_
Definition: HFShowerParam.h:49
int weight
Definition: histoStyle.py:51
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

Member Data Documentation

double HFShowerParam::aperture_
private

Definition at line 49 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

bool HFShowerParam::applyFidCut_
private

Definition at line 50 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

double HFShowerParam::attLMeanInv_
private

Definition at line 49 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

double HFShowerParam::edMin_
private

Definition at line 49 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH2F* HFShowerParam::em_2d_1_
private

Definition at line 56 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH2F * HFShowerParam::em_2d_2_
private

Definition at line 56 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F * HFShowerParam::em_lateral_1_
private

Definition at line 53 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F * HFShowerParam::em_lateral_2_
private

Definition at line 53 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F* HFShowerParam::em_long_1_
private

Definition at line 53 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F * HFShowerParam::em_long_1_tuned_
private

Definition at line 54 of file HFShowerParam.h.

Referenced by HFShowerParam().

TH1F * HFShowerParam::em_long_2_
private

Definition at line 53 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F * HFShowerParam::em_long_gflash_
private

Definition at line 54 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F* HFShowerParam::em_long_sl_
private

Definition at line 55 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

G4int HFShowerParam::emPDG_
private

Definition at line 51 of file HFShowerParam.h.

G4int HFShowerParam::epPDG_
private

Definition at line 51 of file HFShowerParam.h.

bool HFShowerParam::equalizeTimeShift_
private

Definition at line 50 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

std::unique_ptr<HFFibre> HFShowerParam::fibre_
private

Definition at line 46 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

bool HFShowerParam::fillHisto_
private

Definition at line 48 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

G4int HFShowerParam::gammaPDG_
private

Definition at line 51 of file HFShowerParam.h.

std::unique_ptr<HFGflash> HFShowerParam::gflash_
private

Definition at line 47 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

std::vector<double> HFShowerParam::gpar_
private

Definition at line 52 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

const HcalDDDSimConstants* HFShowerParam::hcalConstants_
private

Definition at line 44 of file HFShowerParam.h.

Referenced by HFShowerParam().

TH1F* HFShowerParam::hzvem_
private

Definition at line 54 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F * HFShowerParam::hzvhad_
private

Definition at line 54 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

bool HFShowerParam::onlyLong_
private

Definition at line 50 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

bool HFShowerParam::parametrizeLast_
private

Definition at line 50 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

double HFShowerParam::pePerGeV_
private

Definition at line 49 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

double HFShowerParam::ref_index_
private

Definition at line 49 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

std::unique_ptr<HFShowerLibrary> HFShowerParam::showerLibrary_
private

Definition at line 45 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

bool HFShowerParam::trackEM_
private

Definition at line 50 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().