CMS 3D CMS Logo

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

#include <HFShowerLibrary.h>

Classes

struct  BranchCache
 
struct  BranchReader
 
struct  FileParams
 
struct  Hit
 
struct  Params
 
struct  VersionInfo
 

Public Member Functions

std::vector< HitfillHits (const G4ThreeVector &p, const G4ThreeVector &v, int parCode, double parEnergy, bool &ok, double weight, double time, bool onlyLong=false)
 
std::vector< HitgetHits (const G4Step *aStep, bool &ok, double weight, bool onlyLong=false)
 
 HFShowerLibrary (const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
 
 HFShowerLibrary (Params const &, FileParams const &, HFFibre::Params)
 
 ~HFShowerLibrary ()
 

Private Types

enum  FileFormat { FileFormat::kOld, FileFormat::kNew, FileFormat::kNewV3 }
 

Private Member Functions

HFShowerPhotonCollection extrapolate (int, double)
 
HFShowerPhotonCollection getRecord (int, int) const
 
 HFShowerLibrary (const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &hfShower, edm::ParameterSet const &hfShowerLibrary)
 
HFShowerPhotonCollection interpolate (int, double)
 
VersionInfo loadEventInfo (TBranch *, int fileVersion)
 
bool rInside (double r) const
 
void storePhoton (HFShowerPhoton const &iPhoton, HFShowerPhotonCollection &iPhotons) const
 

Private Attributes

bool applyFidCut_
 
double backProb_
 
double dphi_
 
BranchReader emBranch_
 
bool equalizeTimeShift_
 
int evtPerBin_
 
HFFibre fibre_
 
std::vector< double > gpar_
 
BranchReader hadBranch_
 
std::unique_ptr< TFile > hf_
 
int nMomBin_
 
std::vector< double > pmom_
 
double probMax_
 
double rMax_
 
double rMin_
 
int totEvents_
 
bool verbose_
 

Detailed Description

Definition at line 28 of file HFShowerLibrary.h.

Member Enumeration Documentation

◆ FileFormat

enum HFShowerLibrary::FileFormat
strongprivate
Enumerator
kOld 
kNew 
kNewV3 

Definition at line 90 of file HFShowerLibrary.h.

90 { kOld, kNew, kNewV3 };

Constructor & Destructor Documentation

◆ HFShowerLibrary() [1/3]

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

Definition at line 77 of file HFShowerLibrary.cc.

82  hcons,
83  hps,
84  p.getParameter<edm::ParameterSet>("HFShower").getParameter<edm::ParameterSet>("HFShowerBlock"),
85  p.getParameter<edm::ParameterSet>("HFShowerLibrary").getParameter<edm::ParameterSet>("HFLibraryFileBlock")) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HFShowerLibrary(const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)

◆ ~HFShowerLibrary()

HFShowerLibrary::~HFShowerLibrary ( )

Definition at line 197 of file HFShowerLibrary.cc.

References hf_.

197  {
198  if (hf_)
199  hf_->Close();
200 }
std::unique_ptr< TFile > hf_

◆ HFShowerLibrary() [2/3]

HFShowerLibrary::HFShowerLibrary ( Params const &  iParams,
FileParams const &  iFileParams,
HFFibre::Params  iFibreParams 
)

Definition at line 87 of file HFShowerLibrary.cc.

References cms::cuda::assert(), backProb_, HFShowerLibrary::FileParams::branchEvInfo_, HFShowerLibrary::FileParams::cacheBranches_, visDQMUpload::context, dphi_, emBranch_, HFShowerLibrary::FileParams::emBranchName_, equalizeTimeShift_, edmPickEvents::event, evtPerBin_, Exception, DQMFileSaver_cfi::fileFormat, HFShowerLibrary::FileParams::fileName_, HFShowerLibrary::FileParams::fileVersion_, hadBranch_, HFShowerLibrary::FileParams::hadBranchName_, hf_, mps_fire::i, info(), kNew, kNewV3, kOld, HFShowerLibrary::VersionInfo::libVers_, HFShowerLibrary::VersionInfo::listVersion_, loadEventInfo(), nMomBin_, hltrates_dqm_sourceclient-live_cfg::offset, or, pmom_, probMax_, rMax_, rMin_, AlCaHLTBitMon_QueryRunRegistry::string, totEvents_, and verbose_.

88  : fibre_(iFibreParams),
89  hf_(),
90  emBranch_(),
91  hadBranch_(),
92  verbose_(iParams.verbose_),
93  applyFidCut_(iParams.applyFidCut_),
94  equalizeTimeShift_(iParams.equalizeTimeShift_),
95  probMax_(iParams.probMax_),
96  backProb_(iParams.backProb_),
97  dphi_(iParams.dphi_),
98  rMin_(iFibreParams.rTableHF_.front()),
99  rMax_(iFibreParams.rTableHF_.back()),
100  gpar_(iFibreParams.gParHF_) {
101  std::string pTreeName = iFileParams.fileName_;
102 
103  const char* nTree = pTreeName.c_str();
104  {
105  //It is not safe to open a Tfile in one thread and close in another without adding the following:
106  TDirectory::TContext context;
107  hf_ = std::unique_ptr<TFile>(TFile::Open(nTree));
108  }
109  if (!hf_->IsOpen()) {
110  edm::LogError("HFShower") << "HFShowerLibrary: opening " << nTree << " failed";
111  throw cms::Exception("Unknown", "HFShowerLibrary") << "Opening of " << pTreeName << " fails\n";
112  } else {
113  edm::LogVerbatim("HFShower") << "HFShowerLibrary: opening " << nTree << " successfully";
114  }
115 
117  const int fileVersion = iFileParams.fileVersion_;
118 
119  auto newForm = iFileParams.branchEvInfo_.empty();
120  TTree* event(nullptr);
121  if (newForm) {
123  event = (TTree*)hf_->Get("HFSimHits");
124  } else {
125  event = (TTree*)hf_->Get("Events");
126  }
127  VersionInfo versionInfo;
128  if (event) {
129  TBranch* evtInfo(nullptr);
130  if (!newForm) {
131  std::string info = iFileParams.branchEvInfo_;
132  evtInfo = event->GetBranch(info.c_str());
133  }
134  if (evtInfo || newForm) {
135  versionInfo = loadEventInfo(evtInfo, fileVersion);
136  } else {
137  edm::LogError("HFShower") << "HFShowerLibrary: HFShowerLibrayEventInfo"
138  << " Branch does not exist in Event";
139  throw cms::Exception("Unknown", "HFShowerLibrary") << "Event information absent\n";
140  }
141  } else {
142  edm::LogError("HFShower") << "HFShowerLibrary: Events Tree does not "
143  << "exist";
144  throw cms::Exception("Unknown", "HFShowerLibrary") << "Events tree absent\n";
145  }
146 
147  edm::LogVerbatim("HFShower").log([&](auto& logger) {
148  logger << "HFShowerLibrary: Library " << versionInfo.libVers_ << " ListVersion " << versionInfo.listVersion_
149  << " File version " << fileVersion << " Events Total " << totEvents_ << " and " << evtPerBin_
150  << " per bin\n";
151  logger << "HFShowerLibrary: Energies (GeV) with " << nMomBin_ << " bins\n";
152  for (int i = 0; i < nMomBin_; ++i) {
153  if (i / 10 * 10 == i && i > 0) {
154  logger << "\n";
155  }
156  logger << " " << pmom_[i] / CLHEP::GeV;
157  }
158  });
159 
160  std::string nameBr = iFileParams.emBranchName_;
161  auto emBranch = event->GetBranch(nameBr.c_str());
162  if (verbose_)
163  emBranch->Print();
164  nameBr = iFileParams.hadBranchName_;
165  auto hadBranch = event->GetBranch(nameBr.c_str());
166  if (verbose_)
167  hadBranch->Print();
168 
169  if (emBranch->GetClassName() == std::string("vector<float>")) {
172  }
173  emBranch_ = BranchReader(emBranch, fileFormat, 0, iFileParams.cacheBranches_ ? totEvents_ : 0);
174  size_t offset = 0;
175  if (fileFormat == FileFormat::kNewV3 or (fileFormat == FileFormat::kNew and fileVersion < 2)) {
176  //NOTE: for this format, the hadBranch is all empty up to
177  // totEvents_ (which is more like 1/2*GenEntries())
178  offset = totEvents_;
179  }
180  hadBranch_ = BranchReader(hadBranch, fileFormat, offset, iFileParams.cacheBranches_ ? totEvents_ : 0);
181 
182  edm::LogVerbatim("HFShower") << " HFShowerLibrary:Branch " << iFileParams.emBranchName_ << " has "
183  << emBranch->GetEntries() << " entries and Branch " << iFileParams.hadBranchName_
184  << " has " << hadBranch->GetEntries()
185  << " entries\n HFShowerLibrary::No packing information - Assume x, y, z are not in "
186  "packed form\n Maximum probability cut off "
187  << probMax_ << " Back propagation of light probability " << backProb_
188  << " Flag for equalizing Time Shift for different eta " << equalizeTimeShift_;
189 
190  edm::LogVerbatim("HFShower") << "HFShowerLibrary: rMIN " << rMin_ / CLHEP::cm << " cm and rMax " << rMax_ / CLHEP::cm
191  << " (Half) Phi Width of wedge " << dphi_ / CLHEP::deg;
192  if (iFileParams.cacheBranches_) {
193  hf_.reset();
194  }
195 }
Log< level::Info, true > LogVerbatim
static const TGPicture * info(bool iBackgroundIsBlack)
std::vector< double > gParHF_
Definition: HFFibre.h:34
Log< level::Error, false > LogError
std::vector< double > gpar_
assert(be >=bs)
BranchReader hadBranch_
std::unique_ptr< TFile > hf_
VersionInfo loadEventInfo(TBranch *, int fileVersion)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
Definition: logger.py:1
std::vector< double > pmom_
BranchReader emBranch_
std::vector< double > rTableHF_
Definition: HFFibre.h:35
Definition: event.py:1

◆ HFShowerLibrary() [3/3]

HFShowerLibrary::HFShowerLibrary ( const HcalDDDSimConstants hcons,
const HcalSimulationParameters hps,
edm::ParameterSet const &  hfShower,
edm::ParameterSet const &  hfShowerLibrary 
)
private

Definition at line 69 of file HFShowerLibrary.cc.

73  : HFShowerLibrary(paramsFrom(hfShower, hfShowerLibrary, hcons->getPhiTableHF().front()),
74  fileParamsFrom(hfShowerLibrary),
75  HFFibre::Params(hfShower.getParameter<double>("CFibre"), hcons, hps)) {}
const std::vector< double > & getPhiTableHF() const
HFShowerLibrary(const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)

Member Function Documentation

◆ extrapolate()

HFShowerPhotonCollection HFShowerLibrary::extrapolate ( int  type,
double  pin 
)
private

Definition at line 656 of file HFShowerLibrary.cc.

References evtPerBin_, getRecord(), createfilelist::int, dqmiolumiharvest::j, nMomBin_, l1ctLayer1_cff::nPhoton, displacedMuons_cfi::photon, BPHMonitor_cfi::photons, pmom_, alignCSCRings::r, storePhoton(), totEvents_, and w().

Referenced by fillHits().

656  {
657  int nrec = int(pin / pmom_[nMomBin_ - 1]);
658  double w = (pin - pmom_[nMomBin_ - 1] * nrec) / pmom_[nMomBin_ - 1];
659  nrec++;
660 #ifdef EDM_ML_DEBUG
661  edm::LogVerbatim("HFShower") << "HFShowerLibrary:: Extrapolate for Energy " << pin / CLHEP::GeV << " GeV with "
662  << nMomBin_ << " momentum bins and " << evtPerBin_ << " entries/bin -- "
663  << "total " << totEvents_ << " using " << nrec << " records";
664 #endif
665  std::vector<int> irc(nrec);
666 
667  for (int ir = 0; ir < nrec; ir++) {
668  double r = G4UniformRand();
669  irc[ir] = int(evtPerBin_ * 0.5 * r) + (nMomBin_ - 1) * evtPerBin_ + 1;
670  if (irc[ir] < 1) {
671  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[" << ir << "] = " << irc[ir] << " now set to 1";
672  irc[ir] = 1;
673  } else if (irc[ir] > totEvents_) {
674  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[" << ir << "] = " << irc[ir] << " now set to "
675  << totEvents_;
676  irc[ir] = totEvents_;
677 #ifdef EDM_ML_DEBUG
678  } else {
679  edm::LogVerbatim("HFShower") << "HFShowerLibrary::Extrapolation use irc[" << ir << "] = " << irc[ir];
680 #endif
681  }
682  }
683 
685  std::size_t npold = 0;
686  for (int ir = 0; ir < nrec; ir++) {
687  if (irc[ir] > 0) {
688  auto const photons = getRecord(type, irc[ir]);
689  int nPhoton = photons.size();
690  npold += nPhoton;
691  pe.reserve(pe.size() + nPhoton);
692  for (auto const& photon : photons) {
693  double r = G4UniformRand();
694  if (ir != nrec - 1 || r < w) {
695  storePhoton(photon, pe);
696  }
697  }
698 #ifdef EDM_ML_DEBUG
699  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Record [" << ir << "] = " << irc[ir] << " npold = " << npold;
700 #endif
701  }
702  }
703 #ifdef EDM_ML_DEBUG
704  edm::LogVerbatim("HFShower") << "HFShowerLibrary:: uses " << npold << " photons";
705 #endif
706 
707  if (pe.size() > npold || npold == 0)
708  edm::LogWarning("HFShower") << "HFShowerLibrary: Extrapolation Warning == " << nrec << " records " << irc[0] << ", "
709  << irc[1] << ", ... gives a buffer of " << npold << " photons and fills " << pe.size()
710  << " *****";
711 #ifdef EDM_ML_DEBUG
712  else
713  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Extrapolation == " << nrec << " records " << irc[0] << ", "
714  << irc[1] << ", ... gives a buffer of " << npold << " photons and fills " << pe.size()
715  << " PE";
716  for (std::size_t j = 0; j < pe.size(); j++)
717  edm::LogVerbatim("HFShower") << "Photon " << j << " " << pe[j];
718 #endif
719  return pe;
720 }
Log< level::Info, true > LogVerbatim
T w() const
void storePhoton(HFShowerPhoton const &iPhoton, HFShowerPhotonCollection &iPhotons) const
HFShowerPhotonCollection getRecord(int, int) const
std::vector< HFShowerPhoton > HFShowerPhotonCollection
std::vector< double > pmom_
Log< level::Warning, false > LogWarning

◆ fillHits()

std::vector< HFShowerLibrary::Hit > HFShowerLibrary::fillHits ( const G4ThreeVector &  p,
const G4ThreeVector &  v,
int  parCode,
double  parEnergy,
bool &  ok,
double  weight,
double  time,
bool  onlyLong = false 
)

Definition at line 243 of file HFShowerLibrary.cc.

References funct::abs(), applyFidCut_, HFFibre::attLength(), backProb_, funct::cos(), HFShowerLibrary::Hit::depth, LEDCalibrationChannels::depth, dphi_, equalizeTimeShift_, JetChargeProducer_cfi::exp, extrapolate(), fibre_, gpar_, mps_fire::i, createfilelist::int, interpolate(), G4TrackToParticleID::isGammaElectronPositron(), G4TrackToParticleID::isStableHadron(), nMomBin_, convertSQLiteXML::ok, AlCaHLTBitMon_ParallelJobs::p, displacedMuons_cfi::photon, pmom_, HFShowerLibrary::Hit::position, probMax_, alignCSCRings::r, diffTwoXMLs::r1, diffTwoXMLs::r2, rInside(), funct::sin(), remoteMonitoring_LASER_era2018_cfg::threshold, HFShowerLibrary::Hit::time, funct::true, HFFibre::tShift(), geometryCSVtoXML::xx, geometryCSVtoXML::yy, hit::z, HFFibre::zShift(), gpuVertexFinder::zv, and geometryCSVtoXML::zz.

Referenced by getHits().

250  {
251  std::vector<HFShowerLibrary::Hit> hit;
252  ok = false;
254  // shower is built only for gamma, e+- and stable hadrons
255  if (!isEM && !G4TrackToParticleID::isStableHadron(parCode)) {
256  return hit;
257  }
258  ok = true;
259 
260  // remove low-energy component
261  const double threshold = 50 * MeV;
262  if (pin < threshold) {
263  return hit;
264  }
265 
266  double pz = momDir.z();
267  double zint = hitPoint.z();
268 
269  // if particle moves from interaction point or "backwards (halo)
270  bool backward = (pz * zint < 0.) ? true : false;
271 
272  double sphi = sin(momDir.phi());
273  double cphi = cos(momDir.phi());
274  double ctheta = cos(momDir.theta());
275  double stheta = sin(momDir.theta());
276 
278  if (isEM) {
279  if (pin < pmom_[nMomBin_ - 1]) {
280  pe = interpolate(0, pin);
281  } else {
282  pe = extrapolate(0, pin);
283  }
284  } else {
285  if (pin < pmom_[nMomBin_ - 1]) {
286  pe = interpolate(1, pin);
287  } else {
288  pe = extrapolate(1, pin);
289  }
290  }
291 
292  std::size_t nHit = 0;
293  HFShowerLibrary::Hit oneHit;
294 #ifdef EDM_ML_DEBUG
295  int i = 0;
296 #endif
297  for (auto const& photon : pe) {
298  double zv = std::abs(photon.z()); // abs local z
299 #ifdef EDM_ML_DEBUG
300  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Hit " << i++ << " " << photon << " zv " << zv;
301 #endif
302  if (zv <= gpar_[1] && photon.lambda() > 0 && (photon.z() >= 0 || (zv > gpar_[0] && (!onlyLong)))) {
303  int depth = 1;
304  if (onlyLong) {
305  } else if (!backward) { // fully valid only for "front" particles
306  if (photon.z() < 0)
307  depth = 2; // with "front"-simulated shower lib.
308  } else { // for "backward" particles - almost equal
309  double r = G4UniformRand(); // share between L and S fibers
310  if (r > 0.5)
311  depth = 2;
312  }
313 
314  // Updated coordinate transformation from local
315  // back to global using two Euler angles: phi and theta
316  double pex = photon.x();
317  double pey = photon.y();
318 
319  double xx = pex * ctheta * cphi - pey * sphi + zv * stheta * cphi;
320  double yy = pex * ctheta * sphi + pey * cphi + zv * stheta * sphi;
321  double zz = -pex * stheta + zv * ctheta;
322 
323  G4ThreeVector pos = hitPoint + G4ThreeVector(xx, yy, zz);
324  zv = std::abs(pos.z()) - gpar_[4] - 0.5 * gpar_[1];
325  G4ThreeVector lpos = G4ThreeVector(pos.x(), pos.y(), zv);
326 
327  zv = fibre_.zShift(lpos, depth, 0); // distance to PMT !
328 
329  double r = pos.perp();
330  double p = fibre_.attLength(photon.lambda());
331  double fi = pos.phi();
332  if (fi < 0)
333  fi += CLHEP::twopi;
334  int isect = int(fi / dphi_) + 1;
335  isect = (isect + 1) / 2;
336  double dfi = ((isect * 2 - 1) * dphi_ - fi);
337  if (dfi < 0)
338  dfi = -dfi;
339  double dfir = r * sin(dfi);
340 #ifdef EDM_ML_DEBUG
341  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Position shift " << xx << ", " << yy << ", " << zz << ": "
342  << pos << " R " << r << " Phi " << fi << " Section " << isect << " R*Dfi " << dfir
343  << " Dist " << zv;
344 #endif
345  zz = std::abs(pos.z());
346  double r1 = G4UniformRand();
347  double r2 = G4UniformRand();
348  double r3 = backward ? G4UniformRand() : -9999.;
349  if (!applyFidCut_)
350  dfir += gpar_[5];
351 
352 #ifdef EDM_ML_DEBUG
353  edm::LogVerbatim("HFShower") << "HFShowerLibrary: rLimits " << rInside(r) << " attenuation " << r1 << ":"
354  << exp(-p * zv) << " r2 " << r2 << " r3 " << r3 << " rDfi " << gpar_[5] << " zz "
355  << zz << " zLim " << gpar_[4] << ":" << gpar_[4] + gpar_[1] << "\n"
356  << " rInside(r) :" << rInside(r) << " r1 <= exp(-p*zv) :" << (r1 <= exp(-p * zv))
357  << " r2 <= probMax :" << (r2 <= probMax_ * weight)
358  << " r3 <= backProb :" << (r3 <= backProb_)
359  << " dfir > gpar[5] :" << (dfir > gpar_[5])
360  << " zz >= gpar[4] :" << (zz >= gpar_[4])
361  << " zz <= gpar[4]+gpar[1] :" << (zz <= gpar_[4] + gpar_[1]);
362 #endif
363  if (rInside(r) && r1 <= exp(-p * zv) && r2 <= probMax_ * weight && dfir > gpar_[5] && zz >= gpar_[4] &&
364  zz <= gpar_[4] + gpar_[1] && r3 <= backProb_ && (depth != 2 || zz >= gpar_[4] + gpar_[0])) {
365  double tdiff = (equalizeTimeShift_) ? (fibre_.tShift(lpos, depth, -1)) : (fibre_.tShift(lpos, depth, 1));
366  oneHit.position = pos;
367  oneHit.depth = depth;
368  oneHit.time = (tSlice + (photon.t()) + tdiff);
369  hit.push_back(oneHit);
370 
371 #ifdef EDM_ML_DEBUG
372  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Final Hit " << nHit << " position " << (hit[nHit].position)
373  << " Depth " << (hit[nHit].depth) << " Time " << tSlice << ":" << photon.t() << ":"
374  << tdiff << ":" << (hit[nHit].time);
375 #endif
376  ++nHit;
377  }
378 #ifdef EDM_ML_DEBUG
379  else
380  edm::LogVerbatim("HFShower") << "HFShowerLibrary: REJECTED !!!";
381 #endif
382  if (onlyLong && zz >= gpar_[4] + gpar_[0] && zz <= gpar_[4] + gpar_[1]) {
383  r1 = G4UniformRand();
384  r2 = G4UniformRand();
385  if (rInside(r) && r1 <= exp(-p * zv) && r2 <= probMax_ && dfir > gpar_[5]) {
386  double tdiff = (equalizeTimeShift_) ? (fibre_.tShift(lpos, 2, -1)) : (fibre_.tShift(lpos, 2, 1));
387  oneHit.position = pos;
388  oneHit.depth = 2;
389  oneHit.time = (tSlice + (photon.t()) + tdiff);
390  hit.push_back(oneHit);
391 #ifdef EDM_ML_DEBUG
392  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Final Hit " << nHit << " position " << (hit[nHit].position)
393  << " Depth " << (hit[nHit].depth) << " Time " << tSlice << ":" << photon.t()
394  << ":" << tdiff << ":" << (hit[nHit].time);
395 #endif
396  ++nHit;
397  }
398  }
399  }
400  }
401 
402 #ifdef EDM_ML_DEBUG
403  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Total Hits " << nHit << " out of " << pe.size() << " PE";
404 #endif
405  if (nHit > pe.size() && !onlyLong) {
406  edm::LogWarning("HFShower") << "HFShowerLibrary: Hit buffer " << pe.size() << " smaller than " << nHit << " Hits";
407  }
408  return hit;
409 }
Log< level::Info, true > LogVerbatim
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0) const
Definition: HFFibre.cc:125
float *__restrict__ zv
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Definition: weight.py:1
HFShowerPhotonCollection extrapolate(int, double)
std::vector< double > gpar_
static bool isStableHadron(int pdgCode)
std::vector< HFShowerPhoton > HFShowerPhotonCollection
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
HFShowerPhotonCollection interpolate(int, double)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0) const
Definition: HFFibre.cc:114
std::vector< double > pmom_
G4ThreeVector position
static bool isGammaElectronPositron(int pdgCode)
Log< level::Warning, false > LogWarning
bool rInside(double r) const
double attLength(double lambda) const
Definition: HFFibre.cc:98

◆ getHits()

std::vector< HFShowerLibrary::Hit > HFShowerLibrary::getHits ( const G4Step *  aStep,
bool &  ok,
double  weight,
bool  onlyLong = false 
)

Definition at line 202 of file HFShowerLibrary.cc.

References funct::cos(), fillHits(), gpar_, funct::sin(), and HLT_2023v12_cff::track.

205  {
206  auto const preStepPoint = aStep->GetPreStepPoint();
207  auto const postStepPoint = aStep->GetPostStepPoint();
208  auto const track = aStep->GetTrack();
209  // Get Z-direction
210  auto const aParticle = track->GetDynamicParticle();
211  const G4ThreeVector& momDir = aParticle->GetMomentumDirection();
212  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
213  int parCode = track->GetDefinition()->GetPDGEncoding();
214 
215  // VI: for ions use internally pdg code of alpha in order to keep
216  // consistency with previous simulation
217  if (track->GetDefinition()->IsGeneralIon()) {
218  parCode = 1000020040;
219  }
220 
221 #ifdef EDM_ML_DEBUG
222  G4String partType = track->GetDefinition()->GetParticleName();
223  const G4ThreeVector localPos = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
224  double zoff = localPos.z() + 0.5 * gpar_[1];
225 
226  edm::LogVerbatim("HFShower") << "HFShowerLibrary::getHits " << partType << " of energy "
227  << track->GetKineticEnergy() / CLHEP::GeV << " GeV weight= " << weight
228  << " onlyLong: " << onlyLong << " dir.orts " << momDir.x() << ", " << momDir.y() << ", "
229  << momDir.z() << " Pos x,y,z = " << hitPoint.x() << "," << hitPoint.y() << ","
230  << hitPoint.z() << " (" << zoff << ") sphi,cphi,stheta,ctheta = " << sin(momDir.phi())
231  << "," << cos(momDir.phi()) << ", " << sin(momDir.theta()) << "," << cos(momDir.theta());
232 #endif
233 
234  double tSlice = (postStepPoint->GetGlobalTime()) / CLHEP::nanosecond;
235 
236  // use kinetic energy for protons and ions
237  double pin = (track->GetDefinition()->GetBaryonNumber() > 0) ? preStepPoint->GetKineticEnergy()
238  : preStepPoint->GetTotalEnergy();
239 
240  return fillHits(hitPoint, momDir, parCode, pin, isKilled, weight, tSlice, onlyLong);
241 }
Log< level::Info, true > LogVerbatim
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Definition: weight.py:1
std::vector< double > gpar_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< Hit > fillHits(const G4ThreeVector &p, const G4ThreeVector &v, int parCode, double parEnergy, bool &ok, double weight, double time, bool onlyLong=false)

◆ getRecord()

HFShowerPhotonCollection HFShowerLibrary::getRecord ( int  type,
int  record 
) const
private

Definition at line 527 of file HFShowerLibrary.cc.

References emBranch_, HFShowerLibrary::BranchReader::getRecord(), hadBranch_, dqmiolumiharvest::j, l1ctLayer1_cff::nPhoton, and AlCaHarvesting_cff::record.

Referenced by extrapolate(), and interpolate().

527  {
529  if (type > 0) {
530  photo = hadBranch_.getRecord(record);
531  } else {
532  photo = emBranch_.getRecord(record);
533  }
534 #ifdef EDM_ML_DEBUG
535  int nPhoton = photo.size();
536  edm::LogVerbatim("HFShower") << "HFShowerLibrary::getRecord: Record " << record << " of type " << type << " with "
537  << nPhoton << " photons";
538  for (int j = 0; j < nPhoton; j++)
539  edm::LogVerbatim("HFShower") << "Photon " << j << " " << photo[j];
540 #endif
541  return photo;
542 }
Log< level::Info, true > LogVerbatim
BranchReader hadBranch_
std::vector< HFShowerPhoton > HFShowerPhotonCollection
BranchReader emBranch_
HFShowerPhotonCollection getRecord(int) const

◆ interpolate()

HFShowerPhotonCollection HFShowerLibrary::interpolate ( int  type,
double  pin 
)
private

Definition at line 575 of file HFShowerLibrary.cc.

References evtPerBin_, getRecord(), createfilelist::int, dqmiolumiharvest::j, nMomBin_, l1ctLayer1_cff::nPhoton, displacedMuons_cfi::photon, BPHMonitor_cfi::photons, pmom_, alignCSCRings::r, storePhoton(), totEvents_, and w().

Referenced by fillHits().

575  {
576 #ifdef EDM_ML_DEBUG
577  edm::LogVerbatim("HFShower") << "HFShowerLibrary:: Interpolate for Energy " << pin / CLHEP::GeV << " GeV with "
578  << nMomBin_ << " momentum bins and " << evtPerBin_ << " entries/bin -- total "
579  << totEvents_;
580 #endif
581  int irc[2] = {0, 0};
582  double w = 0.;
583  double r = G4UniformRand();
584 
585  if (pin < pmom_[0]) {
586  w = pin / pmom_[0];
587  irc[1] = int(evtPerBin_ * r) + 1;
588  irc[0] = 0;
589  } else {
590  for (int j = 0; j < nMomBin_ - 1; j++) {
591  if (pin >= pmom_[j] && pin < pmom_[j + 1]) {
592  w = (pin - pmom_[j]) / (pmom_[j + 1] - pmom_[j]);
593  if (j == nMomBin_ - 2) {
594  irc[1] = int(evtPerBin_ * 0.5 * r);
595  } else {
596  irc[1] = int(evtPerBin_ * r);
597  }
598  irc[1] += (j + 1) * evtPerBin_ + 1;
599  r = G4UniformRand();
600  irc[0] = int(evtPerBin_ * r) + 1 + j * evtPerBin_;
601  if (irc[0] < 0) {
602  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[0] = " << irc[0] << " now set to 0";
603  irc[0] = 0;
604  } else if (irc[0] > totEvents_) {
605  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[0] = " << irc[0] << " now set to "
606  << totEvents_;
607  irc[0] = totEvents_;
608  }
609  }
610  }
611  }
612  if (irc[1] < 1) {
613  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[1] = " << irc[1] << " now set to 1";
614  irc[1] = 1;
615  } else if (irc[1] > totEvents_) {
616  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[1] = " << irc[1] << " now set to " << totEvents_;
617  irc[1] = totEvents_;
618  }
619 
620 #ifdef EDM_ML_DEBUG
621  edm::LogVerbatim("HFShower") << "HFShowerLibrary:: Select records " << irc[0] << " and " << irc[1] << " with weights "
622  << 1 - w << " and " << w;
623 #endif
625 
626  std::size_t npold = 0;
627  for (int ir = 0; ir < 2; ir++) {
628  if (irc[ir] > 0) {
629  auto photons = getRecord(type, irc[ir]);
630  int nPhoton = photons.size();
631  npold += nPhoton;
632  pe.reserve(pe.size() + nPhoton);
633  for (auto const& photon : photons) {
634  r = G4UniformRand();
635  if ((ir == 0 && r > w) || (ir > 0 && r < w)) {
636  storePhoton(photon, pe);
637  }
638  }
639  }
640  }
641 
642  if ((pe.size() > npold || (npold == 0 && irc[0] > 0)) && !(pe.empty() && npold == 0))
643  edm::LogWarning("HFShower") << "HFShowerLibrary: Interpolation Warning =="
644  << " records " << irc[0] << " and " << irc[1] << " gives a buffer of " << npold
645  << " photons and fills " << pe.size() << " *****";
646 #ifdef EDM_ML_DEBUG
647  else
648  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Interpolation == records " << irc[0] << " and " << irc[1]
649  << " gives a buffer of " << npold << " photons and fills " << pe.size() << " PE";
650  for (std::size_t j = 0; j < pe.size(); j++)
651  edm::LogVerbatim("HFShower") << "Photon " << j << " " << pe[j];
652 #endif
653  return pe;
654 }
Log< level::Info, true > LogVerbatim
T w() const
void storePhoton(HFShowerPhoton const &iPhoton, HFShowerPhotonCollection &iPhotons) const
HFShowerPhotonCollection getRecord(int, int) const
std::vector< HFShowerPhoton > HFShowerPhotonCollection
std::vector< double > pmom_
Log< level::Warning, false > LogWarning

◆ loadEventInfo()

HFShowerLibrary::VersionInfo HFShowerLibrary::loadEventInfo ( TBranch *  branch,
int  fileVersion 
)
private

Definition at line 544 of file HFShowerLibrary.cc.

References MicroEventContent_cff::branch, evtPerBin_, mps_fire::i, HFShowerLibrary::VersionInfo::libVers_, HFShowerLibrary::VersionInfo::listVersion_, nMomBin_, pmom_, and totEvents_.

Referenced by HFShowerLibrary().

544  {
545  VersionInfo versionInfo;
546  if (branch) {
547  std::vector<HFShowerLibraryEventInfo> eventInfoCollection;
548  branch->SetAddress(&eventInfoCollection);
549  branch->GetEntry(0);
550  edm::LogVerbatim("HFShower") << "HFShowerLibrary::loadEventInfo loads EventInfo Collection of size "
551  << eventInfoCollection.size() << " records";
552 
553  totEvents_ = eventInfoCollection[0].totalEvents();
554  nMomBin_ = eventInfoCollection[0].numberOfBins();
555  evtPerBin_ = eventInfoCollection[0].eventsPerBin();
556  versionInfo.libVers_ = eventInfoCollection[0].showerLibraryVersion();
557  versionInfo.listVersion_ = eventInfoCollection[0].physListVersion();
558  pmom_ = eventInfoCollection[0].energyBins();
559  } else {
560  edm::LogVerbatim("HFShower") << "HFShowerLibrary::loadEventInfo loads EventInfo from hardwired"
561  << " numbers";
562 
563  nMomBin_ = 16;
564  evtPerBin_ = (fileVersion == 0) ? 5000 : 10000;
566  versionInfo.libVers_ = (fileVersion == 0) ? 1.1 : 1.2;
567  versionInfo.listVersion_ = 3.6;
568  pmom_ = {2, 3, 5, 7, 10, 15, 20, 30, 50, 75, 100, 150, 250, 350, 500, 1000};
569  }
570  for (int i = 0; i < nMomBin_; i++)
571  pmom_[i] *= CLHEP::GeV;
572  return versionInfo;
573 }
Log< level::Info, true > LogVerbatim
std::vector< double > pmom_

◆ rInside()

bool HFShowerLibrary::rInside ( double  r) const
private

Definition at line 411 of file HFShowerLibrary.cc.

References alignCSCRings::r, rMax_, and rMin_.

Referenced by fillHits().

411 { return (r >= rMin_ && r <= rMax_); }

◆ storePhoton()

void HFShowerLibrary::storePhoton ( HFShowerPhoton const &  iPhoton,
HFShowerPhotonCollection iPhotons 
) const
private

Definition at line 722 of file HFShowerLibrary.cc.

Referenced by extrapolate(), and interpolate().

722  {
723  iPhotons.push_back(iPhoton);
724 #ifdef EDM_ML_DEBUG
725  edm::LogVerbatim("HFShower") << "HFShowerLibrary: storePhoton " << iPhoton << " npe " << iPhotons.size() << " "
726  << iPhotons.back();
727 #endif
728 }
Log< level::Info, true > LogVerbatim

Member Data Documentation

◆ applyFidCut_

bool HFShowerLibrary::applyFidCut_
private

Definition at line 139 of file HFShowerLibrary.h.

Referenced by fillHits().

◆ backProb_

double HFShowerLibrary::backProb_
private

Definition at line 144 of file HFShowerLibrary.h.

Referenced by fillHits(), and HFShowerLibrary().

◆ dphi_

double HFShowerLibrary::dphi_
private

Definition at line 145 of file HFShowerLibrary.h.

Referenced by fillHits(), and HFShowerLibrary().

◆ emBranch_

BranchReader HFShowerLibrary::emBranch_
private

Definition at line 136 of file HFShowerLibrary.h.

Referenced by getRecord(), and HFShowerLibrary().

◆ equalizeTimeShift_

bool HFShowerLibrary::equalizeTimeShift_
private

Definition at line 143 of file HFShowerLibrary.h.

Referenced by fillHits(), and HFShowerLibrary().

◆ evtPerBin_

int HFShowerLibrary::evtPerBin_
private

Definition at line 140 of file HFShowerLibrary.h.

Referenced by extrapolate(), HFShowerLibrary(), interpolate(), and loadEventInfo().

◆ fibre_

HFFibre HFShowerLibrary::fibre_
private

Definition at line 133 of file HFShowerLibrary.h.

Referenced by fillHits().

◆ gpar_

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

Definition at line 146 of file HFShowerLibrary.h.

Referenced by fillHits(), and getHits().

◆ hadBranch_

BranchReader HFShowerLibrary::hadBranch_
private

Definition at line 137 of file HFShowerLibrary.h.

Referenced by getRecord(), and HFShowerLibrary().

◆ hf_

std::unique_ptr<TFile> HFShowerLibrary::hf_
private

Definition at line 134 of file HFShowerLibrary.h.

Referenced by HFShowerLibrary(), and ~HFShowerLibrary().

◆ nMomBin_

int HFShowerLibrary::nMomBin_
private

Definition at line 140 of file HFShowerLibrary.h.

Referenced by extrapolate(), fillHits(), HFShowerLibrary(), interpolate(), and loadEventInfo().

◆ pmom_

std::vector<double> HFShowerLibrary::pmom_
private

Definition at line 141 of file HFShowerLibrary.h.

Referenced by extrapolate(), fillHits(), HFShowerLibrary(), interpolate(), and loadEventInfo().

◆ probMax_

double HFShowerLibrary::probMax_
private

Definition at line 144 of file HFShowerLibrary.h.

Referenced by fillHits(), and HFShowerLibrary().

◆ rMax_

double HFShowerLibrary::rMax_
private

Definition at line 145 of file HFShowerLibrary.h.

Referenced by HFShowerLibrary(), and rInside().

◆ rMin_

double HFShowerLibrary::rMin_
private

Definition at line 145 of file HFShowerLibrary.h.

Referenced by HFShowerLibrary(), and rInside().

◆ totEvents_

int HFShowerLibrary::totEvents_
private

Definition at line 140 of file HFShowerLibrary.h.

Referenced by extrapolate(), HFShowerLibrary(), interpolate(), and loadEventInfo().

◆ verbose_

bool HFShowerLibrary::verbose_
private

Definition at line 139 of file HFShowerLibrary.h.

Referenced by HFShowerLibrary().