CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
HFShowerParam Class Reference

#include <HFShowerParam.h>

Classes

struct  Hit
 

Public Member Functions

std::vector< HitgetHits (G4Step *aStep, double weight)
 
 HFShowerParam (std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
 
void initRun (G4ParticleTable *, HcalDDDSimConstants *)
 
virtual ~HFShowerParam ()
 

Private Member Functions

std::vector< double > getDDDArray (const std::string &, const DDsvalues_type &)
 

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
 
HFFibrefibre
 
bool fillHisto
 
G4int gammaPDG
 
HFGflashgflash
 
std::vector< double > gpar
 
TH1F * hzvem
 
TH1F * hzvhad
 
bool onlyLong
 
bool parametrizeLast
 
double pePerGeV
 
double ref_index
 
HFShowerLibraryshowerLibrary
 
bool trackEM
 

Detailed Description

Definition at line 27 of file HFShowerParam.h.

Constructor & Destructor Documentation

HFShowerParam::HFShowerParam ( std::string &  name,
const DDCompactView cpv,
edm::ParameterSet const &  p 
)

Definition at line 30 of file HFShowerParam.cc.

References aperture, applyFidCut, HFFibre::attLength(), 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, fibre, fillHisto, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), gflash, g4SimHits_cfi::HFGflash, g4SimHits_cfi::HFShowerLibrary, hzvem, hzvhad, edm::Service< T >::isAvailable(), LogDebug, TFileDirectory::make(), TFileService::mkdir(), onlyLong, parametrizeLast, pePerGeV, ref_index, showerLibrary, compare::tfile, trackEM, and g4SimHits_cfi::useShowerLibrary.

31  : showerLibrary(nullptr),
32  fibre(nullptr), gflash(nullptr),
33  fillHisto(false) {
34  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
35  pePerGeV = m_HF.getParameter<double>("PEPerGeV");
36  trackEM = m_HF.getParameter<bool>("TrackEM");
37  bool useShowerLibrary = m_HF.getParameter<bool>("UseShowerLibrary");
38  bool useGflash = m_HF.getParameter<bool>("UseHFGflash");
39  edMin = m_HF.getParameter<double>("EminLibrary");
40  onlyLong = m_HF.getParameter<bool>("OnlyLong");
41  ref_index = m_HF.getParameter<double>("RefIndex");
42  double lambdaMean = m_HF.getParameter<double>("LambdaMean");
43  aperture = cos(asin(m_HF.getParameter<double>("Aperture")));
44  applyFidCut = m_HF.getParameter<bool>("ApplyFiducialCut");
45  parametrizeLast = m_HF.getUntrackedParameter<bool>("ParametrizeLast",false);
46  edm::LogInfo("HFShower") << "HFShowerParam::Use of shower library is set to "
47  << useShowerLibrary << " Use of Gflash is set to "
48  << useGflash << " P.E. per GeV " << pePerGeV
49  << ", ref. index of fibre " << ref_index
50  << ", Track EM Flag " << trackEM << ", edMin "
51  << edMin << " GeV, use of Short fibre info in"
52  << " shower library set to " << !(onlyLong)
53  << ", use of parametrization for last part set to "
54  << parametrizeLast << ", Mean lambda " << lambdaMean
55  << ", aperture (cutoff) " << aperture
56  << ", Application of Fiducial Cut " << applyFidCut;
57 
58 #ifdef plotDebug
60  if (tfile.isAvailable()) {
61  fillHisto = true;
62  LogDebug("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>("em_2d_1","Lateral Profile vs. Shower Depth;cm;Events",800,800.0,1600.0,100,50.0,150.0);
68  em_long_1 = showerDir.make<TH1F>("em_long_1","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
69  em_long_1_tuned = showerDir.make<TH1F>("em_long_1_tuned","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
70  em_lateral_1 = showerDir.make<TH1F>("em_lateral_1","Lateral Profile;cm;Events",100,50.0,150.0);
71  em_2d_2 = showerDir.make<TH2F>("em_2d_2","Lateral Profile vs. Shower Depth;cm;Events",800,800.0,1600.0,100,50.0,150.0);
72  em_long_2 = showerDir.make<TH1F>("em_long_2","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
73  em_lateral_2 = showerDir.make<TH1F>("em_lateral_2","Lateral Profile;cm;Events",100,50.0,150.0);
74  em_long_gflash = showerDir.make<TH1F>("em_long_gflash","Longitudinal Profile From GFlash;cm;Number of Spots",800,800.0,1600.0);
75  em_long_sl = showerDir.make<TH1F>("em_long_sl","Longitudinal Profile From Shower Library;cm;Number of Spots",800,800.0,1600.0);
76  } else {
77  fillHisto = false;
78  edm::LogInfo("HFShower") << "HFShowerParam::No file is available for "
79  << "saving histos so the flag is set to false";
80  }
81 #endif
82 
84  if (useGflash) gflash = new HFGflash(p);
85  fibre = new HFFibre(name, cpv, p);
86  attLMeanInv = fibre->attLength(lambdaMean);
87  edm::LogInfo("HFShower") << "att. length used for (lambda=" << lambdaMean
88  << ") = " << 1/(attLMeanInv*cm) << " cm";
89 }
#define LogDebug(id)
T getParameter(std::string const &) const
double attLMeanInv
Definition: HFShowerParam.h:55
T getUntrackedParameter(std::string const &, T const &) const
bool parametrizeLast
Definition: HFShowerParam.h:56
double attLength(double lambda)
Definition: HFFibre.cc:94
HFShowerLibrary * showerLibrary
Definition: HFShowerParam.h:52
TH1F * em_lateral_2
Definition: HFShowerParam.h:60
HFFibre * fibre
Definition: HFShowerParam.h:53
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
bool isAvailable() const
Definition: Service.h:46
T * make(const Args &...args) const
make new ROOT object
HFGflash * gflash
Definition: HFShowerParam.h:54
double ref_index
Definition: HFShowerParam.h:55
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
TH1F * em_long_gflash
Definition: HFShowerParam.h:61
TH1F * em_long_sl
Definition: HFShowerParam.h:62
TH1F * em_long_1
Definition: HFShowerParam.h:60
TH1F * em_long_1_tuned
Definition: HFShowerParam.h:61
TH1F * em_lateral_1
Definition: HFShowerParam.h:60
TH1F * em_long_2
Definition: HFShowerParam.h:60
HFShowerParam::~HFShowerParam ( )
virtual

Definition at line 91 of file HFShowerParam.cc.

References fibre, gflash, and showerLibrary.

91  {
92  if (fibre) delete fibre;
93  if (gflash) delete gflash;
94  if (showerLibrary) delete showerLibrary;
95 }
HFShowerLibrary * showerLibrary
Definition: HFShowerParam.h:52
HFFibre * fibre
Definition: HFShowerParam.h:53
HFGflash * gflash
Definition: HFShowerParam.h:54

Member Function Documentation

std::vector< double > HFShowerParam::getDDDArray ( const std::string &  str,
const DDsvalues_type sv 
)
private

Definition at line 403 of file HFShowerParam.cc.

References DDfetch(), DDValue::doubles(), Exception, LogDebug, harvestTrackValidationPlots::str, and relativeConstraints::value.

405 {
406 #ifdef DebugLog
407  LogDebug("HFShower") << "HFShowerParam:getDDDArray called for " << str;
408 #endif
409  DDValue value(str);
410  if (DDfetch(&sv,value))
411  {
412 #ifdef DebugLog
413  LogDebug("HFShower") << value;
414 #endif
415  const std::vector<double> & fvec = value.doubles();
416  int nval = fvec.size();
417  if (nval < 2) {
418  edm::LogError("HFShower") << "HFShowerParam : # of " << str
419  << " bins " << nval << " < 2 ==> illegal";
420  throw cms::Exception("Unknown", "HFShowerParam") << "nval < 2 for array "
421  << str << "\n";
422  }
423  return fvec;
424  } else {
425  edm::LogError("HFShower") << "HFShowerParam : cannot get array " << str;
426  throw cms::Exception("Unknown", "HFShowerParam") << "cannot get array "
427  << str << "\n";
428  }
429 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
Definition: value.py:1
std::vector< HFShowerParam::Hit > HFShowerParam::getHits ( G4Step *  aStep,
double  weight 
)

Definition at line 117 of file HFShowerParam.cc.

References funct::abs(), aperture, applyFidCut, attLMeanInv, gather_cfg::cout, particleFlowClusterECALTimeSelected_cfi::depth, HFShowerParam::Hit::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, emPDG, epPDG, JetChargeProducer_cfi::exp, fibre, fillHisto, gammaPDG, HFShowerLibrary::getHits(), GetVolume(), GeV, gflash, HFGflash::gfParameterization(), gpar, hfClusterShapes_cfi::hits, hzvem, hzvhad, mps_fire::i, cuy::ii, GetRecoTauVFromDQM_MC_cff::kk, convertSQLiteXML::ok, onlyLong, trackingPlots::other, parametrizeLast, callgraph::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, ntuplemaker::time, HiIsolationCommonParameters_cff::track, trackEM, HFFibre::tShift(), mps_merge::weight, z, HFFibre::zShift(), and geometryCSVtoXML::zz.

Referenced by HCalSD::getFromParam().

118  {
119  G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
120  G4Track * track = aStep->GetTrack();
121  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
122  G4int particleCode = track->GetDefinition()->GetPDGEncoding();
123  double zv = std::abs(hitPoint.z()) - gpar[4] - 0.5*gpar[1];
124  G4ThreeVector localPoint = G4ThreeVector(hitPoint.x(),hitPoint.y(),zv);
125 
126  double pin = (preStepPoint->GetTotalEnergy())/GeV;
127  double zint = hitPoint.z();
128  double zz = std::abs(zint) - gpar[4];
129 
130 #ifdef DebugLog
131  edm::LogInfo("HFShower") << "HFShowerParam: getHits "
132  << track->GetDefinition()->GetParticleName()
133  << " of energy " << pin << " GeV"
134  << " Pos x,y,z = " << hitPoint.x() << ","
135  << hitPoint.y() << "," << zint << " (" << zz << ","
136  << localPoint.z() << ", "
137  << (localPoint.z()+0.5*gpar[1]) << ") Local "
138  << localPoint;
139 #endif
140  std::vector<HFShowerParam::Hit> hits;
142  hit.position = hitPoint;
143 
144  // look for other charged particles
145  bool other = false;
146  double pBeta = track->GetDynamicParticle()->GetTotalMomentum() / track->GetDynamicParticle()->GetTotalEnergy();
147  double dirz = (track->GetDynamicParticle()->GetMomentumDirection()).z();
148  if (hitPoint.z() < 0) dirz *= -1.;
149 #ifdef DebugLog
150  edm::LogInfo("HFShower") << "HFShowerParam: getHits Momentum "
151  <<track->GetDynamicParticle()->GetMomentumDirection()
152  << " HitPoint " << hitPoint << " dirz " << dirz;
153 #endif
154  if (particleCode != emPDG && particleCode != epPDG && particleCode != gammaPDG ) {
155  if (track->GetDefinition()->GetPDGCharge() != 0 && pBeta > (1/ref_index) &&
156  aStep->GetTotalEnergyDeposit() > 0) other = true;
157  }
158 
159  // take only e+-/gamma/or special particles
160  if (particleCode == emPDG || particleCode == epPDG ||
161  particleCode == gammaPDG || other) {
162  // Leave out the last part
163  double edep = 0.;
164  bool kill = false;
165  if ((!trackEM) && ((zz<(gpar[1]-gpar[2])) || parametrizeLast) && (!other)){
166  edep = pin;
167  kill = true;
168  } else if ((track->GetDefinition()->GetPDGCharge() != 0) &&
169  (pBeta > (1/ref_index)) && (dirz > aperture)) {
170  edep = (aStep->GetTotalEnergyDeposit())/GeV;
171  }
172  std::string path = "ShowerLibrary";
173 #ifdef DebugLog
174  edm::LogInfo("HFShower") << "HFShowerParam: getHits edep = " << edep
175  << " weight " << weight << " final " <<edep*weight
176  << ", Kill = " << kill << ", pin = " << pin
177  << ", edMin = " << edMin << " Other " << other;
178 #endif
179  edep *= weight;
180  if (edep > 0) {
181  if ((showerLibrary || gflash) && kill && pin > edMin && (!other)) {
182  if (showerLibrary) {
183  std::vector<HFShowerLibrary::Hit> hitSL = showerLibrary->getHits(aStep,kill,weight,onlyLong);
184  for (unsigned int i=0; i<hitSL.size(); i++) {
185  bool ok = true;
186 #ifdef DebugLog
187  edm::LogInfo("HFShower") << "HFShowerParam: getHits applyFidCut = " << applyFidCut;
188 #endif
189  if (applyFidCut) { // @@ For showerlibrary no z-cut for Short (no z)
190  int npmt = HFFibreFiducial:: PMTNumber(hitSL[i].position);
191  if (npmt <= 0) ok = false;
192  }
193  if (ok) {
194  hit.position = hitSL[i].position;
195  hit.depth = hitSL[i].depth;
196  hit.time = hitSL[i].time;
197  hit.edep = 1;
198  hits.push_back(hit);
199 #ifdef plotDebug
200  if (fillHisto) {
201  double zv = std::abs(hit.position.z()) - gpar[4];
202  hzvem->Fill(zv);
203  em_long_sl->Fill(hit.position.z()/cm);
204  double sq = sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2));
205  double zp = hit.position.z()/cm;
206  if (hit.depth == 1) {
207  em_2d_1->Fill(zp, sq);
208  em_lateral_1->Fill(sq);
209  em_long_1->Fill(zp);
210  } else if (hit.depth == 2) {
211  em_2d_2->Fill(zp, sq);
212  em_lateral_2->Fill(sq);
213  em_long_2->Fill(zp);
214  }
215  }
216 #endif
217 #ifdef DebugLog
218  edm::LogInfo("HFShower") << "HFShowerParam: Hit at depth "
219  << hit.depth << " with edep " << hit.edep
220  << " Time " << hit.time;
221 #endif
222  }
223  }
224  } else { // GFlash clusters with known z
225  std::vector<HFGflash::Hit>hitSL=gflash->gfParameterization(aStep,kill, onlyLong);
226  for (unsigned int i=0; i<hitSL.size(); ++i) {
227  bool ok = true;
228  G4ThreeVector pe_effect(hitSL[i].position.x(), hitSL[i].position.y(),
229  hitSL[i].position.z());
230  double zv = std::abs(pe_effect.z()) - gpar[4];
231  //depth
232  int depth = 1;
233  int npmt = 0;
234  if (zv < 0. || zv > gpar[1]) {
235 #ifdef mkdebug
236  std::cout<<"-#Zcut-HFShowerParam::getHits:z="<<zv<<",m="<<gpar[1]<<std::endl;
237 #endif
238  ok = false;
239  }
240  if (ok && applyFidCut) {
241  npmt = HFFibreFiducial:: PMTNumber(pe_effect);
242 #ifdef DebugLog
243  edm::LogInfo("HFShower") << "HFShowerParam::getHits:#PMT= "
244  << npmt << ",z = " << zv;
245 #endif
246  if (npmt <= 0) {
247 #ifdef DebugLog
248  edm::LogInfo("HFShower") << "-#PMT=0 cut-HFShowerParam::"
249  << "getHits: npmt = " << npmt;
250 #endif
251  ok = false;
252  } else if (npmt > 24) { // a short fibre
253  if (zv > gpar[0]) {
254  depth = 2;
255  } else {
256 #ifdef DebugLog
257  edm::LogInfo("HFShower") << "-SHORT cut-HFShowerParam::"
258  << "getHits:zMin=" << gpar[0];
259 #endif
260  ok = false;
261  }
262  }
263 #ifdef DebugLog
264  edm::LogInfo("HFShower") << "HFShowerParam: npmt " << npmt
265  << " zv " << std::abs(pe_effect.z())
266  << ":" << gpar[4] << ":" << zv << ":"
267  << gpar[0] << " ok " << ok << " depth "
268  << depth;
269 #endif
270  } else {
271  if (G4UniformRand() > 0.5) depth = 2;
272  if (depth == 2 && zv < gpar[0]) ok = false;
273  }
274  //attenuation
275  double dist = fibre->zShift(localPoint,depth,0); // distance to PMT
276  double r1 = G4UniformRand();
277 #ifdef DebugLog
278  edm::LogInfo("HFShower") << "HFShowerParam:Distance to PMT (" <<npmt
279  << ") " << dist << ", exclusion flag "
280  << (r1 > exp(-attLMeanInv*zv));
281 #endif
282  if (r1 > exp(-attLMeanInv*dist)) ok = false;
283  if (ok) {
284  double r2 = G4UniformRand();
285 #ifdef DebugLog
286  edm::LogInfo("HFShower") << "HFShowerParam:Extra exclusion "
287  << r2 << ">" << weight << " "
288  << (r2 > weight);
289 #endif
290  if (r2 < weight) {
291  double time = fibre->tShift(localPoint,depth,0);
292 
293  hit.position = hitSL[i].position;
294  hit.depth = depth;
295  hit.time = time + hitSL[i].time;
296  hit.edep = 1;
297  hits.push_back(hit);
298 #ifdef plotDebug
299  if (fillHisto) {
300  em_long_gflash->Fill(pe_effect.z()/cm, hitSL[i].edep);
301  hzvem->Fill(zv);
302  double sq = sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2));
303  double zp = hit.position.z()/cm;
304  if (hit.depth == 1) {
305  em_2d_1->Fill(zp, sq);
306  em_lateral_1->Fill(s);
307  em_long_1->Fill(zp);
308  } else if (hit.depth == 2) {
309  em_2d_2->Fill(zp, sq);
310  em_lateral_2->Fill(sq);
311  em_long_2->Fill(zp);
312  }
313  }
314 #endif
315 #ifdef DebugLog
316  edm::LogInfo("HFShower") << "HFShowerParam: Hit at depth "
317  << hit.depth << " with edep "
318  << hit.edep << " Time " << hit.time;
319 #endif
320  }
321  }
322  }
323  }
324  } else {
325  path = "Rest";
326  edep *= pePerGeV;
327  double tSlice = (aStep->GetPostStepPoint()->GetGlobalTime());
328  double time = fibre->tShift(localPoint,1,0); // remaining part
329  bool ok = true;
330  if (applyFidCut) { // @@ For showerlibrary no z-cut for Short (no z)
331  int npmt = HFFibreFiducial:: PMTNumber(hitPoint);
332  if (npmt <= 0) ok = false;
333  }
334 #ifdef DebugLog
335  edm::LogInfo("HFShower") << "HFShowerParam: getHits hitPoint " << hitPoint << " flag " << ok;
336 #endif
337  if (ok) {
338  hit.depth = 1;
339  hit.time = tSlice+time;
340  hit.edep = edep;
341  hits.push_back(hit);
342 #ifdef DebugLog
343  edm::LogInfo("HFShower") << "HFShowerParam: Hit at depth 1 with edep "
344  << edep << " Time " << tSlice << ":" << time
345  << ":" << hit.time;
346 #endif
347 #ifdef plotDebug
348  double zv = std::abs(hitPoint.z()) - gpar[4];
349  if (fillHisto) {
350  hzvhad->Fill(zv);
351  }
352 #endif
353  if (zz >= gpar[0]) {
354  time = fibre->tShift(localPoint,2,0);
355  hit.depth = 2;
356  hit.time = tSlice+time;
357  hits.push_back(hit);
358 #ifdef DebugLog
359  edm::LogInfo("HFShower") <<"HFShowerParam: Hit at depth 2 with edep "
360  << edep << " Time " << tSlice << ":" << time
361  << hit.time;
362 #endif
363 #ifdef plotDebug
364  if (fillHisto) {
365  hzvhad->Fill(zv);
366  }
367 #endif
368  }
369  }
370  }
371 #ifdef DebugLog
372  for (unsigned int ii=0; ii<hits.size(); ++ii) {
373  double zv = std::abs(hits[ii].position.z());
374  if (zv > 12790) edm::LogInfo("HFShower")<< "HFShowerParam: Abnormal hit along "
375  << path << " in "
376  << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName()
377  << " at " << hits[ii].position << " zz "
378  << zv << " Edep " << edep << " due to "
379  <<track->GetDefinition()->GetParticleName()
380  << " time " << hit.time;
381  }
382 #endif
383  if (kill) {
384  track->SetTrackStatus(fStopAndKill);
385  G4TrackVector tv = *(aStep->GetSecondary());
386  for (unsigned int kk=0; kk<tv.size(); ++kk) {
387  if (tv[kk]->GetVolume() == preStepPoint->GetPhysicalVolume())
388  tv[kk]->SetTrackStatus(fStopAndKill);
389  }
390  }
391 #ifdef DebugLog
392  edm::LogInfo("HFShower") << "HFShowerParam: getHits kill (" << kill
393  << ") track " << track->GetTrackID()
394  << " at " << hitPoint
395  << " and deposit " << edep << " " << hits.size()
396  << " times" << " ZZ " << zz << " " << gpar[0];
397 #endif
398  }
399  }
400  return hits;
401 }
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:112
std::vector< Hit > gfParameterization(G4Step *aStep, bool &ok, bool onlyLong=false)
Definition: HFGflash.cc:91
double attLMeanInv
Definition: HFShowerParam.h:55
bool parametrizeLast
Definition: HFShowerParam.h:56
const double GeV
Definition: MathUtil.h:16
std::vector< Hit > getHits(G4Step *aStep, bool &ok, double weight, bool onlyLong=false)
std::vector< double > gpar
Definition: HFShowerParam.h:58
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:124
Definition: weight.py:1
HFShowerLibrary * showerLibrary
Definition: HFShowerParam.h:52
TH1F * em_lateral_2
Definition: HFShowerParam.h:60
HFFibre * fibre
Definition: HFShowerParam.h:53
static int PMTNumber(const G4ThreeVector &pe_effect)
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HFGflash * gflash
Definition: HFShowerParam.h:54
ii
Definition: cuy.py:588
G4ThreeVector position
Definition: HFShowerParam.h:39
static const G4LogicalVolume * GetVolume(const std::string &name)
double ref_index
Definition: HFShowerParam.h:55
TH1F * em_long_gflash
Definition: HFShowerParam.h:61
TH1F * em_long_sl
Definition: HFShowerParam.h:62
TH1F * em_long_1
Definition: HFShowerParam.h:60
static int position[264][3]
Definition: ReadPGInfo.cc:509
TH1F * em_lateral_1
Definition: HFShowerParam.h:60
TH1F * em_long_2
Definition: HFShowerParam.h:60
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void HFShowerParam::initRun ( G4ParticleTable *  theParticleTable,
HcalDDDSimConstants hcons 
)

Definition at line 97 of file HFShowerParam.cc.

References emPDG, epPDG, fibre, gammaPDG, HcalDDDSimConstants::getGparHF(), gpar, HFFibre::initRun(), HFShowerLibrary::initRun(), and showerLibrary.

Referenced by HCalSD::initRun().

98  {
99  emPDG = theParticleTable->FindParticle("e-")->GetPDGEncoding();
100  epPDG = theParticleTable->FindParticle("e+")->GetPDGEncoding();
101  gammaPDG = theParticleTable->FindParticle("gamma")->GetPDGEncoding();
102 #ifdef DebugLog
103  edm::LogInfo("HFShower") << "HFShowerParam: Particle code for e- = " << emPDG
104  << " for e+ = " << epPDG << " for gamma = " << gammaPDG;
105 #endif
106  if (showerLibrary) showerLibrary->initRun(theParticleTable, hcons);
107  if (fibre) fibre->initRun(hcons);
108 
109  //Special Geometry parameters
110  gpar = hcons->getGparHF();
111  edm::LogInfo("HFShower") << "HFShowerParam: " << gpar.size() <<" gpar (cm)";
112  for (unsigned int ig=0; ig<gpar.size(); ig++)
113  edm::LogInfo("HFShower") << "HFShowerParam: gpar[" << ig << "] = "
114  << gpar[ig]/cm << " cm";
115 }
void initRun(HcalDDDSimConstants *)
Definition: HFFibre.cc:77
std::vector< double > gpar
Definition: HFShowerParam.h:58
HFShowerLibrary * showerLibrary
Definition: HFShowerParam.h:52
HFFibre * fibre
Definition: HFShowerParam.h:53
const std::vector< double > & getGparHF() const
void initRun(G4ParticleTable *, HcalDDDSimConstants *)

Member Data Documentation

double HFShowerParam::aperture
private

Definition at line 55 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

bool HFShowerParam::applyFidCut
private

Definition at line 56 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

double HFShowerParam::attLMeanInv
private

Definition at line 55 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

double HFShowerParam::edMin
private

Definition at line 55 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH2F* HFShowerParam::em_2d_1
private

Definition at line 63 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH2F * HFShowerParam::em_2d_2
private

Definition at line 63 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F * HFShowerParam::em_lateral_1
private

Definition at line 60 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F * HFShowerParam::em_lateral_2
private

Definition at line 60 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F* HFShowerParam::em_long_1
private

Definition at line 60 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F * HFShowerParam::em_long_1_tuned
private

Definition at line 61 of file HFShowerParam.h.

Referenced by HFShowerParam().

TH1F * HFShowerParam::em_long_2
private

Definition at line 60 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F * HFShowerParam::em_long_gflash
private

Definition at line 61 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F* HFShowerParam::em_long_sl
private

Definition at line 62 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

G4int HFShowerParam::emPDG
private

Definition at line 57 of file HFShowerParam.h.

Referenced by getHits(), and initRun().

G4int HFShowerParam::epPDG
private

Definition at line 57 of file HFShowerParam.h.

Referenced by getHits(), and initRun().

HFFibre* HFShowerParam::fibre
private

Definition at line 53 of file HFShowerParam.h.

Referenced by getHits(), HFShowerParam(), initRun(), and ~HFShowerParam().

bool HFShowerParam::fillHisto
private

Definition at line 59 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

G4int HFShowerParam::gammaPDG
private

Definition at line 57 of file HFShowerParam.h.

Referenced by getHits(), and initRun().

HFGflash* HFShowerParam::gflash
private

Definition at line 54 of file HFShowerParam.h.

Referenced by getHits(), HFShowerParam(), and ~HFShowerParam().

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

Definition at line 58 of file HFShowerParam.h.

Referenced by getHits(), and initRun().

TH1F* HFShowerParam::hzvem
private

Definition at line 61 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

TH1F * HFShowerParam::hzvhad
private

Definition at line 61 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

bool HFShowerParam::onlyLong
private

Definition at line 56 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

bool HFShowerParam::parametrizeLast
private

Definition at line 56 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

double HFShowerParam::pePerGeV
private

Definition at line 55 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

double HFShowerParam::ref_index
private

Definition at line 55 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().

HFShowerLibrary* HFShowerParam::showerLibrary
private

Definition at line 52 of file HFShowerParam.h.

Referenced by getHits(), HFShowerParam(), initRun(), and ~HFShowerParam().

bool HFShowerParam::trackEM
private

Definition at line 56 of file HFShowerParam.h.

Referenced by getHits(), and HFShowerParam().