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