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 
11 
15 
16 #include "G4VPhysicalVolume.hh"
17 #include "G4Step.hh"
18 #include "G4Track.hh"
19 #include "G4NavigationHistory.hh"
20 #include "Randomize.hh"
21 
22 #include "CLHEP/Units/GlobalPhysicalConstants.h"
23 #include "CLHEP/Units/GlobalSystemOfUnits.h"
24 
25 #include <iostream>
26 
27 //#define EDM_ML_DEBUG
28 //#define plotDebug
29 //#define mkdebug
30 
32  : showerLibrary(nullptr), fibre(nullptr), gflash(nullptr), fillHisto(false) {
33  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
34  pePerGeV = m_HF.getParameter<double>("PEPerGeV");
35  trackEM = m_HF.getParameter<bool>("TrackEM");
36  bool useShowerLibrary = m_HF.getParameter<bool>("UseShowerLibrary");
37  bool useGflash = m_HF.getParameter<bool>("UseHFGflash");
38  edMin = m_HF.getParameter<double>("EminLibrary");
39  onlyLong = m_HF.getParameter<bool>("OnlyLong");
40  ref_index = m_HF.getParameter<double>("RefIndex");
41  double lambdaMean = m_HF.getParameter<double>("LambdaMean");
42  aperture = cos(asin(m_HF.getParameter<double>("Aperture")));
43  applyFidCut = m_HF.getParameter<bool>("ApplyFiducialCut");
44  parametrizeLast = m_HF.getUntrackedParameter<bool>("ParametrizeLast", false);
45  edm::LogVerbatim("HFShower") << "HFShowerParam::Use of shower library is set to " << useShowerLibrary
46  << " Use of Gflash is set to " << useGflash << " P.E. per GeV " << pePerGeV
47  << ", ref. index of fibre " << ref_index << ", Track EM Flag " << trackEM << ", edMin "
48  << edMin << " GeV, use of Short fibre info in"
49  << " shower library set to " << !(onlyLong)
50  << ", use of parametrization for last part set to " << parametrizeLast
51  << ", Mean lambda " << lambdaMean << ", aperture (cutoff) " << aperture
52  << ", Application of Fiducial Cut " << applyFidCut;
53 
54 #ifdef plotDebug
56  if (tfile.isAvailable()) {
57  fillHisto = true;
58  edm::LogVerbatim("HFShower") << "HFShowerParam::Save histos in directory "
59  << "ProfileFromParam";
60  TFileDirectory showerDir = tfile->mkdir("ProfileFromParam");
61  hzvem = showerDir.make<TH1F>("hzvem", "Longitudinal Profile (EM Part);Number of PE", 330, 0.0, 1650.0);
62  hzvhad = showerDir.make<TH1F>("hzvhad", "Longitudinal Profile (Had Part);Number of PE", 330, 0.0, 1650.0);
63  em_2d_1 = showerDir.make<TH2F>(
64  "em_2d_1", "Lateral Profile vs. Shower Depth;cm;Events", 800, 800.0, 1600.0, 100, 50.0, 150.0);
65  em_long_1 =
66  showerDir.make<TH1F>("em_long_1", "Longitudinal Profile;Radiation Length;Number of Spots", 800, 800.0, 1600.0);
67  em_long_1_tuned = showerDir.make<TH1F>(
68  "em_long_1_tuned", "Longitudinal Profile;Radiation Length;Number of Spots", 800, 800.0, 1600.0);
69  em_lateral_1 = showerDir.make<TH1F>("em_lateral_1", "Lateral Profile;cm;Events", 100, 50.0, 150.0);
70  em_2d_2 = showerDir.make<TH2F>(
71  "em_2d_2", "Lateral Profile vs. Shower Depth;cm;Events", 800, 800.0, 1600.0, 100, 50.0, 150.0);
72  em_long_2 =
73  showerDir.make<TH1F>("em_long_2", "Longitudinal Profile;Radiation Length;Number of Spots", 800, 800.0, 1600.0);
74  em_lateral_2 = showerDir.make<TH1F>("em_lateral_2", "Lateral Profile;cm;Events", 100, 50.0, 150.0);
75  em_long_gflash = showerDir.make<TH1F>(
76  "em_long_gflash", "Longitudinal Profile From GFlash;cm;Number of Spots", 800, 800.0, 1600.0);
77  em_long_sl = showerDir.make<TH1F>(
78  "em_long_sl", "Longitudinal Profile From Shower Library;cm;Number of Spots", 800, 800.0, 1600.0);
79  } else {
80  fillHisto = false;
81  edm::LogVerbatim("HFShower") << "HFShowerParam::No file is available for saving histos so the "
82  << "flag is set to false";
83  }
84 #endif
85 
86  if (useShowerLibrary)
87  showerLibrary = new HFShowerLibrary(name, cpv, p);
88  if (useGflash)
89  gflash = new HFGflash(p);
90  fibre = new HFFibre(name, cpv, p);
91  attLMeanInv = fibre->attLength(lambdaMean);
92  edm::LogVerbatim("HFShower") << "att. length used for (lambda=" << lambdaMean << ") = " << 1 / (attLMeanInv * cm)
93  << " cm";
94 }
95 
97  delete fibre;
98  delete gflash;
99  delete showerLibrary;
100 }
101 
103  if (showerLibrary)
104  showerLibrary->initRun(nullptr, hcons);
105  if (fibre)
106  fibre->initRun(hcons);
107 }
108 
109 std::vector<HFShowerParam::Hit> HFShowerParam::getHits(const G4Step* aStep, double weight, bool& isKilled) {
110  auto const preStepPoint = aStep->GetPreStepPoint();
111  auto const track = aStep->GetTrack();
113  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
114  double zv = std::abs(hitPoint.z()) - gpar[4] - 0.5 * gpar[1];
115  G4ThreeVector localPoint = G4ThreeVector(hitPoint.x(), hitPoint.y(), zv);
116 
117  double pin = (preStepPoint->GetTotalEnergy()) / GeV;
118  double zint = hitPoint.z();
119  double zz = std::abs(zint) - gpar[4];
120 
121 #ifdef EDM_ML_DEBUG
122  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits " << track->GetDefinition()->GetParticleName()
123  << " of energy " << pin << " GeV Pos x,y,z = " << hitPoint.x() << "," << hitPoint.y()
124  << "," << zint << " (" << zz << "," << localPoint.z() << ", "
125  << (localPoint.z() + 0.5 * gpar[1]) << ") Local " << localPoint;
126 #endif
127  std::vector<HFShowerParam::Hit> hits;
129  hit.position = hitPoint;
130 
131  // look for other charged particles
132  bool other = false;
133  double pBeta = track->GetDynamicParticle()->GetTotalMomentum() / track->GetDynamicParticle()->GetTotalEnergy();
134  double dirz = (track->GetDynamicParticle()->GetMomentumDirection()).z();
135  if (hitPoint.z() < 0)
136  dirz *= -1.;
137 #ifdef EDM_ML_DEBUG
138  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits Momentum "
139  << track->GetDynamicParticle()->GetMomentumDirection() << " HitPoint " << hitPoint
140  << " dirz " << dirz;
141 #endif
142  if (!isEM && track->GetDefinition()->GetPDGCharge() != 0 && pBeta > (1 / ref_index) &&
143  aStep->GetTotalEnergyDeposit() > 0.) {
144  other = true;
145  }
146 
147  // take only e+-/gamma/or special particles
148  if (isEM || other) {
149  // Leave out the last part
150  double edep = 0.;
151  if ((!trackEM) && ((zz < (gpar[1] - gpar[2])) || parametrizeLast) && (!other)) {
152  edep = pin;
153  isKilled = true;
154  } else if ((track->GetDefinition()->GetPDGCharge() != 0) && (pBeta > (1 / ref_index)) && (dirz > aperture)) {
155  edep = (aStep->GetTotalEnergyDeposit()) / GeV;
156  }
157  std::string path = "ShowerLibrary";
158 #ifdef EDM_ML_DEBUG
159  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits edep = " << edep << " weight " << weight << " final "
160  << edep * weight << ", Kill = " << isKilled << ", pin = " << pin
161  << ", edMin = " << edMin << " Other " << other;
162 #endif
163  edep *= weight;
164  if (edep > 0) {
165  if ((showerLibrary || gflash) && isKilled && pin > edMin && (!other)) {
166  if (showerLibrary) {
167  std::vector<HFShowerLibrary::Hit> hitSL = showerLibrary->getHits(aStep, isKilled, weight, onlyLong);
168  for (unsigned int i = 0; i < hitSL.size(); i++) {
169  bool ok = true;
170 #ifdef EDM_ML_DEBUG
171  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits applyFidCut = " << applyFidCut;
172 #endif
173  if (applyFidCut) { // @@ For showerlibrary no z-cut for Short (no z)
174  int npmt = HFFibreFiducial::PMTNumber(hitSL[i].position);
175  if (npmt <= 0)
176  ok = false;
177  }
178  if (ok) {
179  hit.position = hitSL[i].position;
180  hit.depth = hitSL[i].depth;
181  hit.time = hitSL[i].time;
182  hit.edep = 1;
183  hits.push_back(hit);
184 #ifdef plotDebug
185  if (fillHisto) {
186  double zv = std::abs(hit.position.z()) - gpar[4];
187  hzvem->Fill(zv);
188  em_long_sl->Fill(hit.position.z() / cm);
189  double sq = sqrt(pow(hit.position.x() / cm, 2) + pow(hit.position.y() / cm, 2));
190  double zp = hit.position.z() / cm;
191  if (hit.depth == 1) {
192  em_2d_1->Fill(zp, sq);
193  em_lateral_1->Fill(sq);
194  em_long_1->Fill(zp);
195  } else if (hit.depth == 2) {
196  em_2d_2->Fill(zp, sq);
197  em_lateral_2->Fill(sq);
198  em_long_2->Fill(zp);
199  }
200  }
201 #endif
202 #ifdef EDM_ML_DEBUG
203  edm::LogVerbatim("HFShower")
204  << "HFShowerParam: Hit at depth " << hit.depth << " with edep " << hit.edep << " Time " << hit.time;
205 #endif
206  }
207  }
208  } else { // GFlash clusters with known z
209  std::vector<HFGflash::Hit> hitSL = gflash->gfParameterization(aStep, onlyLong);
210  for (unsigned int i = 0; i < hitSL.size(); ++i) {
211  bool ok = true;
212  G4ThreeVector pe_effect(hitSL[i].position.x(), hitSL[i].position.y(), hitSL[i].position.z());
213  double zv = std::abs(pe_effect.z()) - gpar[4];
214  //depth
215  int depth = 1;
216  int npmt = 0;
217  if (zv < 0. || zv > gpar[1]) {
218 #ifdef mkdebug
219  std::cout << "-#Zcut-HFShowerParam::getHits:z=" << zv << ",m=" << gpar[1] << std::endl;
220 #endif
221  ok = false;
222  }
223  if (ok && applyFidCut) {
224  npmt = HFFibreFiducial::PMTNumber(pe_effect);
225 #ifdef EDM_ML_DEBUG
226  edm::LogVerbatim("HFShower") << "HFShowerParam::getHits:#PMT= " << npmt << ",z = " << zv;
227 #endif
228  if (npmt <= 0) {
229 #ifdef EDM_ML_DEBUG
230  edm::LogVerbatim("HFShower") << "-#PMT=0 cut-HFShowerParam::getHits: npmt = " << npmt;
231 #endif
232  ok = false;
233  } else if (npmt > 24) { // a short fibre
234  if (zv > gpar[0]) {
235  depth = 2;
236  } else {
237 #ifdef EDM_ML_DEBUG
238  edm::LogVerbatim("HFShower") << "-SHORT cut-HFShowerParam::getHits:zMin=" << gpar[0];
239 #endif
240  ok = false;
241  }
242  }
243 #ifdef EDM_ML_DEBUG
244  edm::LogVerbatim("HFShower")
245  << "HFShowerParam: npmt " << npmt << " zv " << std::abs(pe_effect.z()) << ":" << gpar[4] << ":" << zv
246  << ":" << gpar[0] << " ok " << ok << " depth " << depth;
247 #endif
248  } else {
249  if (G4UniformRand() > 0.5)
250  depth = 2;
251  if (depth == 2 && zv < gpar[0])
252  ok = false;
253  }
254  //attenuation
255  double dist = fibre->zShift(localPoint, depth, 0); // distance to PMT
256  double r1 = G4UniformRand();
257 #ifdef EDM_ML_DEBUG
258  edm::LogVerbatim("HFShower") << "HFShowerParam:Distance to PMT (" << npmt << ") " << dist
259  << ", exclusion flag " << (r1 > exp(-attLMeanInv * zv));
260 #endif
261  if (r1 > exp(-attLMeanInv * dist))
262  ok = false;
263  if (ok) {
264  double r2 = G4UniformRand();
265 #ifdef EDM_ML_DEBUG
266  edm::LogVerbatim("HFShower")
267  << "HFShowerParam:Extra exclusion " << r2 << ">" << weight << " " << (r2 > weight);
268 #endif
269  if (r2 < weight) {
270  double time = 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() / cm, hitSL[i].edep);
280  hzvem->Fill(zv);
281  double sq = sqrt(pow(hit.position.x() / cm, 2) + pow(hit.position.y() / cm, 2));
282  double zp = hit.position.z() / 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 = fibre->tShift(localPoint, 1, 0); // remaining part
307  bool ok = true;
308  if (applyFidCut) { // @@ For showerlibrary no z-cut for Short (no z)
309  int npmt = HFFibreFiducial::PMTNumber(hitPoint);
310  if (npmt <= 0)
311  ok = false;
312  }
313 #ifdef EDM_ML_DEBUG
314  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits hitPoint " << hitPoint << " flag " << ok;
315 #endif
316  if (ok) {
317  hit.depth = 1;
318  hit.time = tSlice + time;
319  hit.edep = edep;
320  hits.push_back(hit);
321 #ifdef EDM_ML_DEBUG
322  edm::LogVerbatim("HFShower") << "HFShowerParam: Hit at depth 1 with edep " << edep << " Time " << tSlice
323  << ":" << time << ":" << hit.time;
324 #endif
325 #ifdef plotDebug
326  double zv = std::abs(hitPoint.z()) - gpar[4];
327  if (fillHisto) {
328  hzvhad->Fill(zv);
329  }
330 #endif
331  if (zz >= gpar[0]) {
332  time = fibre->tShift(localPoint, 2, 0);
333  hit.depth = 2;
334  hit.time = tSlice + time;
335  hits.push_back(hit);
336 #ifdef EDM_ML_DEBUG
337  edm::LogVerbatim("HFShower") << "HFShowerParam: Hit at depth 2 with edep " << edep << " Time " << tSlice
338  << ":" << time << hit.time;
339 #endif
340 #ifdef plotDebug
341  if (fillHisto) {
342  hzvhad->Fill(zv);
343  }
344 #endif
345  }
346  }
347  }
348 #ifdef EDM_ML_DEBUG
349  for (unsigned int ii = 0; ii < hits.size(); ++ii) {
350  double zv = std::abs(hits[ii].position.z());
351  if (zv > 12790)
352  edm::LogVerbatim("HFShower") << "HFShowerParam: Abnormal hit along " << path << " in "
353  << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName() << " at "
354  << hits[ii].position << " zz " << zv << " Edep " << edep << " due to "
355  << track->GetDefinition()->GetParticleName() << " time " << hit.time;
356  }
357  edm::LogVerbatim("HFShower") << "HFShowerParam: getHits kill (" << isKilled << ") track " << track->GetTrackID()
358  << " at " << hitPoint << " and deposit " << edep << " " << hits.size() << " times"
359  << " ZZ " << zz << " " << gpar[0];
360 #endif
361  }
362  }
363  return hits;
364 }
365 
366 std::vector<double> HFShowerParam::getDDDArray(const std::string& str, const DDsvalues_type& sv) {
367 #ifdef EDM_ML_DEBUG
368  edm::LogVerbatim("HFShower") << "HFShowerParam:getDDDArray called for " << str;
369 #endif
370  DDValue value(str);
371  if (DDfetch(&sv, value)) {
372 #ifdef EDM_ML_DEBUG
373  edm::LogVerbatim("HFShower") << value;
374 #endif
375  const std::vector<double>& fvec = value.doubles();
376  int nval = fvec.size();
377  if (nval < 2) {
378  edm::LogError("HFShower") << "HFShowerParam : # of " << str << " bins " << nval << " < 2 ==> illegal";
379  throw cms::Exception("Unknown", "HFShowerParam") << "nval < 2 for array " << str << "\n";
380  }
381  return fvec;
382  } else {
383  edm::LogError("HFShower") << "HFShowerParam : cannot get array " << str;
384  throw cms::Exception("Unknown", "HFShowerParam") << "cannot get array " << str << "\n";
385  }
386 }
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:113
T getParameter(std::string const &) const
HFShowerParam(const std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
double attLMeanInv
Definition: HFShowerParam.h:49
T getUntrackedParameter(std::string const &, T const &) const
std::vector< Hit > gfParameterization(const G4Step *aStep, bool onlyLong=false)
Definition: HFGflash.cc:128
std::vector< Hit > getHits(const G4Step *aStep, bool &ok, double weight, bool onlyLong=false)
const std::vector< double > & doubles() const
a reference to the double-valued values stored in the given instance of DDValue
Definition: DDValue.cc:140
bool parametrizeLast
Definition: HFShowerParam.h:50
const double GeV
Definition: MathUtil.h:16
void initRun(const HcalDDDSimConstants *)
Definition: HFFibre.cc:81
virtual ~HFShowerParam()
double attLength(double lambda)
Definition: HFFibre.cc:97
#define nullptr
void initRun(const HcalDDDSimConstants *)
std::vector< double > gpar
Definition: HFShowerParam.h:52
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:123
Definition: weight.py:1
HFShowerLibrary * showerLibrary
Definition: HFShowerParam.h:46
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
TH1F * em_lateral_2
Definition: HFShowerParam.h:54
HFFibre * fibre
Definition: HFShowerParam.h:47
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
bool isAvailable() const
Definition: Service.h:40
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< Hit > getHits(const G4Step *aStep, double weight, bool &isKilled)
T * make(const Args &...args) const
make new ROOT object
HFGflash * gflash
Definition: HFShowerParam.h:48
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
ii
Definition: cuy.py:590
G4ThreeVector position
Definition: HFShowerParam.h:34
int PMTNumber(const G4ThreeVector &pe_effect)
double ref_index
Definition: HFShowerParam.h:49
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:55
TH1F * em_long_sl
Definition: HFShowerParam.h:56
TH1F * em_long_1
Definition: HFShowerParam.h:54
TH1F * em_long_1_tuned
Definition: HFShowerParam.h:55
static int position[264][3]
Definition: ReadPGInfo.cc:509
static bool isGammaElectronPositron(int pdgCode)
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
Definition: DDsvalues.h:12
TH1F * em_lateral_1
Definition: HFShowerParam.h:54
TH1F * em_long_2
Definition: HFShowerParam.h:54
#define str(s)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void initRun(G4ParticleTable *, const HcalDDDSimConstants *)