12 #include "G4VPhysicalVolume.hh" 13 #include "G4NavigationHistory.hh" 16 #include "G4ParticleTable.hh" 17 #include "Randomize.hh" 18 #include "CLHEP/Units/SystemOfUnits.h" 19 #include "CLHEP/Units/PhysicalConstants.h" 45 if (
params.fileName_.find(
'.') == 0)
46 params.fileName_.erase(0, 2);
51 "BranchEvt",
"HFShowerLibraryEventInfos_hfshowerlib_HFShowerLibraryEventInfo");
57 params.emBranchName_ = branchPre + emName + branchPost;
58 params.hadBranchName_ = branchPre + hadName + branchPost;
60 if (not branchEvInfo.empty()) {
61 params.branchEvInfo_ = branchEvInfo + branchPost;
73 :
HFShowerLibrary(paramsFrom(hfShower, hfShowerLibrary, hcons->getPhiTableHF().front()),
74 fileParamsFrom(hfShowerLibrary),
75 HFFibre::
Params(hfShower.getParameter<double>(
"CFibre"), hcons, hps)) {}
88 : fibre_(iFibreParams),
92 verbose_(iParams.verbose_),
93 applyFidCut_(iParams.applyFidCut_),
94 equalizeTimeShift_(iParams.equalizeTimeShift_),
95 probMax_(iParams.probMax_),
96 backProb_(iParams.backProb_),
98 rMin_(iFibreParams.rTableHF_.front()),
99 rMax_(iFibreParams.rTableHF_.back()),
100 gpar_(iFibreParams.gParHF_) {
103 const char* nTree = pTreeName.c_str();
107 hf_ = std::unique_ptr<TFile>(TFile::Open(nTree));
109 if (!
hf_->IsOpen()) {
110 edm::LogError(
"HFShower") <<
"HFShowerLibrary: opening " << nTree <<
" failed";
111 throw cms::Exception(
"Unknown",
"HFShowerLibrary") <<
"Opening of " << pTreeName <<
" fails\n";
113 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: opening " << nTree <<
" successfully";
120 TTree*
event(
nullptr);
123 event = (TTree*)
hf_->Get(
"HFSimHits");
125 event = (TTree*)
hf_->Get(
"Events");
129 TBranch* evtInfo(
nullptr);
132 evtInfo =
event->GetBranch(
info.c_str());
134 if (evtInfo || newForm) {
137 edm::LogError(
"HFShower") <<
"HFShowerLibrary: HFShowerLibrayEventInfo" 138 <<
" Branch does not exist in Event";
139 throw cms::Exception(
"Unknown",
"HFShowerLibrary") <<
"Event information absent\n";
142 edm::LogError(
"HFShower") <<
"HFShowerLibrary: Events Tree does not " 144 throw cms::Exception(
"Unknown",
"HFShowerLibrary") <<
"Events tree absent\n";
151 logger <<
"HFShowerLibrary: Energies (GeV) with " <<
nMomBin_ <<
" bins\n";
153 if (
i / 10 * 10 ==
i &&
i > 0) {
161 auto emBranch =
event->GetBranch(nameBr.c_str());
165 auto hadBranch =
event->GetBranch(nameBr.c_str());
169 if (emBranch->GetClassName() ==
std::string(
"vector<float>")) {
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 " 191 <<
" (Half) Phi Width of wedge " <<
dphi_ / CLHEP::deg;
206 auto const preStepPoint = aStep->GetPreStepPoint();
207 auto const postStepPoint = aStep->GetPostStepPoint();
208 auto const track = aStep->GetTrack();
210 auto const aParticle =
track->GetDynamicParticle();
211 const G4ThreeVector& momDir = aParticle->GetMomentumDirection();
212 const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
213 int parCode =
track->GetDefinition()->GetPDGEncoding();
217 if (
track->GetDefinition()->IsGeneralIon()) {
218 parCode = 1000020040;
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];
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());
234 double tSlice = (postStepPoint->GetGlobalTime()) / CLHEP::nanosecond;
237 double pin = (
track->GetDefinition()->GetBaryonNumber() > 0) ? preStepPoint->GetKineticEnergy()
238 : preStepPoint->GetTotalEnergy();
240 return fillHits(hitPoint, momDir, parCode, pin, isKilled,
weight, tSlice, onlyLong);
244 const G4ThreeVector& momDir,
251 std::vector<HFShowerLibrary::Hit>
hit;
266 double pz = momDir.
z();
267 double zint = hitPoint.z();
270 bool backward = (pz * zint < 0.) ?
true :
false;
272 double sphi =
sin(momDir.phi());
273 double cphi =
cos(momDir.phi());
274 double ctheta =
cos(momDir.theta());
275 double stheta =
sin(momDir.theta());
292 std::size_t nHit = 0;
297 for (
auto const&
photon : pe) {
305 }
else if (!backward) {
309 double r = G4UniformRand();
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;
323 G4ThreeVector
pos = hitPoint + G4ThreeVector(
xx,
yy,
zz);
325 G4ThreeVector lpos = G4ThreeVector(
pos.x(),
pos.y(),
zv);
329 double r =
pos.perp();
331 double fi =
pos.phi();
335 isect = (isect + 1) / 2;
336 double dfi = ((isect * 2 - 1) *
dphi_ - fi);
339 double dfir =
r *
sin(dfi);
341 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Position shift " <<
xx <<
", " <<
yy <<
", " <<
zz <<
": " 342 <<
pos <<
" R " <<
r <<
" Phi " << fi <<
" Section " << isect <<
" R*Dfi " << dfir
346 double r1 = G4UniformRand();
347 double r2 = G4UniformRand();
348 double r3 = backward ? G4UniformRand() : -9999.;
354 <<
exp(-
p *
zv) <<
" r2 " <<
r2 <<
" r3 " << r3 <<
" rDfi " <<
gpar_[5] <<
" zz " 356 <<
" rInside(r) :" <<
rInside(
r) <<
" r1 <= exp(-p*zv) :" << (
r1 <=
exp(-
p *
zv))
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]);
364 zz <= gpar_[4] + gpar_[1] && r3 <= backProb_ && (depth != 2 || zz >=
gpar_[4] +
gpar_[0])) {
368 oneHit.
time = (tSlice + (
photon.t()) + tdiff);
369 hit.push_back(oneHit);
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);
383 r1 = G4UniformRand();
384 r2 = G4UniformRand();
389 oneHit.
time = (tSlice + (
photon.t()) + tdiff);
390 hit.push_back(oneHit);
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);
403 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Total Hits " << nHit <<
" out of " << pe.size() <<
" PE";
405 if (nHit > pe.size() && !onlyLong) {
406 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary: Hit buffer " << pe.size() <<
" smaller than " << nHit <<
" Hits";
415 if (nRecords > maxRecordsToCache) {
416 nRecords = maxRecordsToCache;
422 for (
size_t r = 0;
r < nRecords; ++
r) {
430 for (
size_t r = 0;
r < nRecords; ++
r) {
441 cache_ = makeCache(*
this, maxRecordsToCache, branch_->GetDirectory()->GetFile()->GetName(), branch_->GetName());
448 std::lock_guard<std::mutex> guard(
s_mutex);
449 CMS_SA_ALLOW static std::map<std::pair<std::string, std::string>, std::weak_ptr<BranchCache>> s_map;
450 auto v = s_map[{iFileName, iBranchName}].lock();
454 v = std::make_shared<BranchCache>(iReader, maxRecordsToCache);
455 s_map[{iFileName, iBranchName}] =
v;
461 assert(static_cast<size_t>(iRecord + 1) < offsets_.size());
463 auto start = offsets_[iRecord - 1];
464 auto end = offsets_[iRecord];
471 iBranch->SetAddress(&photo);
472 iBranch->GetEntry(iEntry);
478 auto temp = std::make_unique<HFShowerPhotonCollection>();
479 iBranch->SetAddress(&
temp);
480 iBranch->GetEntry(iEntry);
488 std::vector<float>
t;
489 std::vector<float>*
tp = &
t;
490 iBranch->SetAddress(&
tp);
491 iBranch->GetEntry(iEntry);
492 unsigned int tSize =
t.size() / 5;
493 photo.reserve(tSize);
494 for (
unsigned int i = 0;
i < tSize;
i++) {
495 photo.emplace_back(
t[
i],
t[1 * tSize +
i],
t[2 * tSize +
i],
t[3 * tSize +
i],
t[4 * tSize +
i]);
503 return cache_->getRecord(
record);
510 photo = getRecordNewForm(branch_, nrc + offset_);
514 photo = getRecordNewFormV3(branch_, nrc + offset_);
518 photo = getRecordOldForm(branch_, nrc);
547 std::vector<HFShowerLibraryEventInfo> eventInfoCollection;
548 branch->SetAddress(&eventInfoCollection);
550 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary::loadEventInfo loads EventInfo Collection of size " 551 << eventInfoCollection.size() <<
" records";
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();
560 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary::loadEventInfo loads EventInfo from hardwired" 564 evtPerBin_ = (fileVersion == 0) ? 5000 : 10000;
566 versionInfo.
libVers_ = (fileVersion == 0) ? 1.1 : 1.2;
568 pmom_ = {2, 3, 5, 7, 10, 15, 20, 30, 50, 75, 100, 150, 250, 350, 500, 1000};
577 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary:: Interpolate for Energy " << pin / CLHEP::GeV <<
" GeV with " 583 double r = G4UniformRand();
585 if (pin <
pmom_[0]) {
602 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary:: Illegal irc[0] = " << irc[0] <<
" now set to 0";
605 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary:: Illegal irc[0] = " << irc[0] <<
" now set to " 613 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary:: Illegal irc[1] = " << irc[1] <<
" now set to 1";
621 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary:: Select records " << irc[0] <<
" and " << irc[1] <<
" with weights " 622 << 1 -
w <<
" and " <<
w;
626 std::size_t npold = 0;
627 for (
int ir = 0; ir < 2; ir++) {
632 pe.reserve(pe.size() +
nPhoton);
635 if ((ir == 0 &&
r >
w) || (ir > 0 &&
r <
w)) {
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() <<
" *****";
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++)
661 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary:: Extrapolate for Energy " << pin / CLHEP::GeV <<
" GeV with " 663 <<
"total " <<
totEvents_ <<
" using " << nrec <<
" records";
665 std::vector<int> irc(nrec);
667 for (
int ir = 0; ir < nrec; ir++) {
668 double r = G4UniformRand();
671 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary:: Illegal irc[" << ir <<
"] = " << irc[ir] <<
" now set to 1";
674 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary:: Illegal irc[" << ir <<
"] = " << irc[ir] <<
" now set to " 679 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary::Extrapolation use irc[" << ir <<
"] = " << irc[ir];
685 std::size_t npold = 0;
686 for (
int ir = 0; ir < nrec; ir++) {
691 pe.reserve(pe.size() +
nPhoton);
693 double r = G4UniformRand();
694 if (ir != nrec - 1 ||
r <
w) {
699 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Record [" << ir <<
"] = " << irc[ir] <<
" npold = " << npold;
704 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary:: uses " << npold <<
" photons";
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()
713 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Extrapolation == " << nrec <<
" records " << irc[0] <<
", " 714 << irc[1] <<
", ... gives a buffer of " << npold <<
" photons and fills " << pe.size()
716 for (std::size_t
j = 0;
j < pe.size();
j++)
723 iPhotons.push_back(iPhoton);
725 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: storePhoton " << iPhoton <<
" npe " << iPhotons.size() <<
" "
Log< level::Info, true > LogVerbatim
static HFShowerPhotonCollection getRecordOldForm(TBranch *, int iEntry)
BranchCache(BranchReader &, size_t maxRecordsToCache)
T getParameter(std::string const &) const
std::vector< Hit > getHits(const G4Step *aStep, bool &ok, double weight, bool onlyLong=false)
static std::recursive_mutex s_mutex
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0) const
void storePhoton(HFShowerPhoton const &iPhoton, HFShowerPhotonCollection &iPhotons) const
HFShowerPhotonCollection getRecord(int, int) const
std::string emBranchName_
std::string hadBranchName_
std::size_t numberOfRecords() const
Sin< T >::type sin(const T &t)
HFShowerPhotonCollection extrapolate(int, double)
Log< level::Error, false > LogError
std::vector< double > gpar_
void doCaching(size_t maxRecords)
static bool isStableHadron(int pdgCode)
T getUntrackedParameter(std::string const &, T const &) const
std::unique_ptr< TFile > hf_
VersionInfo loadEventInfo(TBranch *, int fileVersion)
std::vector< HFShowerPhoton > HFShowerPhotonCollection
Cos< T >::type cos(const T &t)
HFShowerPhotonCollection interpolate(int, double)
Abs< T >::type abs(const T &t)
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0) const
std::string branchEvInfo_
static std::shared_ptr< BranchCache const > makeCache(BranchReader &, size_t maxRecordsToCache, std::string const &iFileName, std::string const &iBranchName)
std::vector< double > pmom_
HFShowerPhotonCollection getRecord(int) const
HFShowerLibrary(const std::string &name, const HcalDDDSimConstants *hcons, const HcalSimulationParameters *hps, edm::ParameterSet const &p)
static HFShowerPhotonCollection getRecordNewForm(TBranch *, int iEntry)
std::vector< std::size_t > offsets_
static HFShowerPhotonCollection getRecordNewFormV3(TBranch *, int iEntry)
HFShowerPhotonCollection photons_
static bool isGammaElectronPositron(int pdgCode)
std::vector< Hit > fillHits(const G4ThreeVector &p, const G4ThreeVector &v, int parCode, double parEnergy, bool &ok, double weight, double time, bool onlyLong=false)
Log< level::Warning, false > LogWarning
bool rInside(double r) const
HFShowerPhotonCollection getRecord(int) const
double attLength(double lambda) const