15 #include "G4VPhysicalVolume.hh" 16 #include "G4NavigationHistory.hh" 19 #include "G4ParticleTable.hh" 20 #include "Randomize.hh" 21 #include "CLHEP/Units/SystemOfUnits.h" 22 #include "CLHEP/Units/PhysicalConstants.h" 34 backProb = m_HS.getParameter<
double>(
"BackProbability");
38 "BranchEvt",
"HFShowerLibraryEventInfos_hfshowerlib_HFShowerLibraryEventInfo");
39 std::string branchPre = m_HS.getUntrackedParameter<
std::string>(
"BranchPre",
"HFShowerPhotons_hfshowerlib_");
41 verbose = m_HS.getUntrackedParameter<
bool>(
"Verbosity",
false);
42 applyFidCut = m_HS.getParameter<
bool>(
"ApplyFiducialCut");
44 if (pTreeName.find(
".") == 0)
45 pTreeName.erase(0, 2);
46 const char* nTree = pTreeName.c_str();
47 hf = TFile::Open(nTree);
50 edm::LogError(
"HFShower") <<
"HFShowerLibrary: opening " << nTree <<
" failed";
51 throw cms::Exception(
"Unknown",
"HFShowerLibrary") <<
"Opening of " << pTreeName <<
" fails\n";
53 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: opening " << nTree <<
" successfully";
56 newForm = (branchEvInfo.empty());
57 TTree*
event(
nullptr);
59 event = (TTree*)
hf->Get(
"HFSimHits");
61 event = (TTree*)
hf->Get(
"Events");
63 TBranch* evtInfo(
nullptr);
66 evtInfo =
event->GetBranch(info.c_str());
71 edm::LogError(
"HFShower") <<
"HFShowerLibrary: HFShowerLibrayEventInfo" 72 <<
" Branch does not exist in Event";
73 throw cms::Exception(
"Unknown",
"HFShowerLibrary") <<
"Event information absent\n";
76 edm::LogError(
"HFShower") <<
"HFShowerLibrary: Events Tree does not " 78 throw cms::Exception(
"Unknown",
"HFShowerLibrary") <<
"Events tree absent\n";
84 ss <<
"HFShowerLibrary: Energies (GeV) with " <<
nMomBin <<
" bins\n";
86 if (
i / 10 * 10 ==
i &&
i > 0) {
93 std::string nameBr = branchPre + emName + branchPost;
94 emBranch =
event->GetBranch(nameBr.c_str());
97 nameBr = branchPre + hadName + branchPost;
98 hadBranch =
event->GetBranch(nameBr.c_str());
108 <<
" entries and Branch " << hadName <<
" has " <<
hadBranch->GetEntries() <<
" entries" 109 <<
"\n HFShowerLibrary::No packing information -" 110 <<
" Assume x, y, z are not in packed form" 111 <<
"\n Maximum probability cut off " << probMax <<
" Back propagation of light prob. " 132 rMax = rTable[rTable.size() - 1];
138 <<
" (Half) Phi Width of wedge " <<
dphi / deg;
148 auto const preStepPoint = aStep->GetPreStepPoint();
149 auto const postStepPoint = aStep->GetPostStepPoint();
150 auto const track = aStep->GetTrack();
152 auto const aParticle =
track->GetDynamicParticle();
153 const G4ThreeVector& momDir = aParticle->GetMomentumDirection();
154 const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
155 int parCode =
track->GetDefinition()->GetPDGEncoding();
159 if (
track->GetDefinition()->IsGeneralIon()) {
160 parCode = 1000020040;
164 G4String partType =
track->GetDefinition()->GetParticleName();
165 const G4ThreeVector localPos = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
166 double zoff = localPos.z() + 0.5 *
gpar[1];
168 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary::getHits " << partType <<
" of energy " 169 <<
track->GetKineticEnergy() /
GeV <<
" GeV weight= " << weight
170 <<
" onlyLong: " << onlyLong <<
" dir.orts " << momDir.x() <<
", " << momDir.y() <<
", " 171 << momDir.z() <<
" Pos x,y,z = " << hitPoint.x() <<
"," << hitPoint.y() <<
"," 172 << hitPoint.z() <<
" (" << zoff <<
") sphi,cphi,stheta,ctheta = " <<
sin(momDir.phi())
173 <<
"," <<
cos(momDir.phi()) <<
", " <<
sin(momDir.theta()) <<
"," <<
cos(momDir.theta());
176 double tSlice = (postStepPoint->GetGlobalTime()) / CLHEP::nanosecond;
179 double pin = (
track->GetDefinition()->GetBaryonNumber() > 0) ? preStepPoint->GetKineticEnergy()
180 : preStepPoint->GetTotalEnergy();
182 return fillHits(hitPoint, momDir, parCode, pin, isKilled, weight, tSlice, onlyLong);
186 const G4ThreeVector& momDir,
193 std::vector<HFShowerLibrary::Hit>
hit;
204 if (pin < threshold) {
208 double pz = momDir.z();
209 double zint = hitPoint.z();
212 bool backward = (pz * zint < 0.) ?
true :
false;
214 double sphi =
sin(momDir.phi());
215 double cphi =
cos(momDir.phi());
216 double ctheta =
cos(momDir.theta());
217 double stheta =
sin(momDir.theta());
235 for (
int i = 0;
i <
npe; ++
i) {
240 if (zv <=
gpar[1] &&
pe[
i].lambda() > 0 && (
pe[
i].
z() >= 0 || (zv >
gpar[0] && (!onlyLong)))) {
243 }
else if (!backward) {
247 double r = G4UniformRand();
254 double pex =
pe[
i].x();
255 double pey =
pe[
i].y();
257 double xx = pex * ctheta * cphi - pey * sphi + zv * stheta * cphi;
258 double yy = pex * ctheta * sphi + pey * cphi + zv * stheta * sphi;
259 double zz = -pex * stheta + zv * ctheta;
261 G4ThreeVector
pos = hitPoint + G4ThreeVector(xx, yy, zz);
263 G4ThreeVector lpos = G4ThreeVector(pos.x(), pos.y(), zv);
267 double r = pos.perp();
269 double fi = pos.phi();
272 int isect =
int(fi /
dphi) + 1;
273 isect = (isect + 1) / 2;
274 double dfi = ((isect * 2 - 1) *
dphi - fi);
277 double dfir = r *
sin(dfi);
279 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Position shift " << xx <<
", " << yy <<
", " << zz <<
": " 280 << pos <<
" R " << r <<
" Phi " << fi <<
" Section " << isect <<
" R*Dfi " << dfir
284 double r1 = G4UniformRand();
285 double r2 = G4UniformRand();
286 double r3 = backward ? G4UniformRand() : -9999.;
292 <<
exp(-p * zv) <<
" r2 " << r2 <<
" r3 " << r3 <<
" rDfi " <<
gpar[5] <<
" zz " 293 << zz <<
" zLim " <<
gpar[4] <<
":" <<
gpar[4] +
gpar[1] <<
"\n" 294 <<
" rInside(r) :" <<
rInside(r) <<
" r1 <= exp(-p*zv) :" << (r1 <=
exp(-p * zv))
295 <<
" r2 <= probMax :" << (r2 <=
probMax * weight)
296 <<
" r3 <= backProb :" << (r3 <=
backProb)
297 <<
" dfir > gpar[5] :" << (dfir >
gpar[5]) <<
" zz >= gpar[4] :" << (zz >=
gpar[4])
298 <<
" zz <= gpar[4]+gpar[1] :" << (zz <=
gpar[4] +
gpar[1]);
300 if (
rInside(r) && r1 <=
exp(-p * zv) && r2 <= probMax * weight && dfir >
gpar[5] && zz >=
gpar[4] &&
301 zz <= gpar[4] + gpar[1] && r3 <= backProb && (depth != 2 || zz >=
gpar[4] +
gpar[0])) {
305 hit.push_back(oneHit);
307 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Final Hit " << nHit <<
" position " << (hit[nHit].position)
308 <<
" Depth " << (hit[nHit].depth) <<
" Time " << tSlice <<
":" <<
pe[
i].t() <<
":" 309 <<
fibre->
tShift(lpos, depth, 1) <<
":" << (hit[nHit].time);
318 r1 = G4UniformRand();
319 r2 = G4UniformRand();
320 if (
rInside(r) && r1 <=
exp(-p * zv) && r2 <= probMax && dfir >
gpar[5]) {
324 hit.push_back(oneHit);
326 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Final Hit " << nHit <<
" position " << (hit[nHit].position)
327 <<
" Depth " << (hit[nHit].depth) <<
" Time " << tSlice <<
":" <<
pe[
i].t()
328 <<
":" <<
fibre->
tShift(lpos, 2, 1) <<
":" << (hit[nHit].time);
337 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Total Hits " << nHit <<
" out of " << npe <<
" PE";
339 if (nHit > npe && !onlyLong) {
340 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary: Hit buffer " << npe <<
" smaller than " << nHit <<
" Hits";
348 int nrc = record - 1;
357 std::vector<float>
t;
358 std::vector<float>* tp = &
t;
361 unsigned int tSize = t.size() / 5;
362 photo->reserve(tSize);
363 for (
unsigned int i = 0;
i < tSize;
i++) {
365 HFShowerPhoton(t[
i], t[1 * tSize + i], t[2 * tSize + i], t[3 * tSize + i], t[4 * tSize + i]));
378 std::vector<float>
t;
379 std::vector<float>* tp = &
t;
382 unsigned int tSize = t.size() / 5;
383 photo->reserve(tSize);
384 for (
unsigned int i = 0;
i < tSize;
i++) {
386 HFShowerPhoton(t[
i], t[1 * tSize + i], t[2 * tSize + i], t[3 * tSize + i], t[4 * tSize + i]));
396 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary::getRecord: Record " << record <<
" of type " << type <<
" with " 397 << nPhoton <<
" photons";
398 for (
int j = 0; j < nPhoton; j++)
408 std::vector<HFShowerLibraryEventInfo> eventInfoCollection;
409 branch->SetAddress(&eventInfoCollection);
412 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary::loadEventInfo loads EventInfo Collection of size " 413 << eventInfoCollection.size() <<
" records";
415 totEvents = eventInfoCollection[0].totalEvents();
416 nMomBin = eventInfoCollection[0].numberOfBins();
417 evtPerBin = eventInfoCollection[0].eventsPerBin();
418 libVers = eventInfoCollection[0].showerLibraryVersion();
419 listVersion = eventInfoCollection[0].physListVersion();
420 pmom = eventInfoCollection[0].energyBins();
423 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary::loadEventInfo loads EventInfo from hardwired" 431 pmom = {2, 3, 5, 7, 10, 15, 20, 30, 50, 75, 100, 150, 250, 350, 500, 1000};
444 double r = G4UniformRand();
451 for (
int j = 0; j <
nMomBin - 1; j++) {
452 if (pin >=
pmom[j] && pin <
pmom[j + 1]) {
454 if (j == nMomBin - 2) {
463 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary:: Illegal irc[0] = " << irc[0] <<
" now set to 0";
465 }
else if (irc[0] > totEvents) {
473 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary:: Illegal irc[1] = " << irc[1] <<
" now set to 1";
475 }
else if (irc[1] > totEvents) {
481 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary:: Select records " << irc[0] <<
" and " << irc[1] <<
" with weights " 482 << 1 - w <<
" and " <<
w;
487 for (
int ir = 0; ir < 2; ir++) {
492 for (
int j = 0; j < nPhoton; j++) {
494 if ((ir == 0 && r > w) || (ir > 0 && r < w)) {
501 if ((
npe > npold || (npold == 0 && irc[0] > 0)) && !(
npe == 0 && npold == 0))
502 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary: Interpolation Warning ==" 503 <<
" records " << irc[0] <<
" and " << irc[1] <<
" gives a buffer of " << npold
504 <<
" photons and fills " <<
npe <<
" *****";
507 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Interpolation == records " << irc[0] <<
" and " << irc[1]
508 <<
" gives a buffer of " << npold <<
" photons and fills " << npe <<
" PE";
509 for (
int j = 0; j <
npe; j++)
520 <<
" momentum bins and " <<
evtPerBin <<
" entries/bin -- " 521 <<
"total " <<
totEvents <<
" using " << nrec <<
" records";
523 std::vector<int> irc(nrec);
525 for (
int ir = 0; ir < nrec; ir++) {
526 double r = G4UniformRand();
529 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary:: Illegal irc[" << ir <<
"] = " << irc[ir] <<
" now set to 1";
532 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary:: Illegal irc[" << ir <<
"] = " << irc[ir] <<
" now set to " 537 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary::Extrapolation use irc[" << ir <<
"] = " << irc[ir];
545 for (
int ir = 0; ir < nrec; ir++) {
550 for (
int j = 0; j < nPhoton; j++) {
551 double r = G4UniformRand();
552 if (ir != nrec - 1 || r < w) {
557 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Record [" << ir <<
"] = " << irc[ir] <<
" npold = " << npold;
562 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary:: uses " << npold <<
" photons";
565 if (
npe > npold || npold == 0)
566 edm::LogWarning(
"HFShower") <<
"HFShowerLibrary: Extrapolation Warning == " << nrec <<
" records " << irc[0] <<
", " 567 << irc[1] <<
", ... gives a buffer of " << npold <<
" photons and fills " <<
npe 571 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary: Extrapolation == " << nrec <<
" records " << irc[0] <<
", " 572 << irc[1] <<
", ... gives a buffer of " << npold <<
" photons and fills " << npe
574 for (
int j = 0; j <
npe; j++)
592 edm::LogVerbatim(
"HFShower") <<
"HFShowerLibrary:getDDDArray called for " << str <<
" with nMin " << nmin;
599 const std::vector<double>& fvec = value.
doubles();
600 int nval = fvec.size();
603 edm::LogError(
"HFShower") <<
"HFShowerLibrary : # of " << str <<
" bins " << nval <<
" < " << nmin
605 throw cms::Exception(
"Unknown",
"HFShowerLibrary") <<
"nval < nmin for array " << str <<
"\n";
609 edm::LogError(
"HFShower") <<
"HFShowerLibrary : # of " << str <<
" bins " << nval <<
" < 2 ==> illegal" 610 <<
" (nmin=" << nmin <<
")";
611 throw cms::Exception(
"Unknown",
"HFShowerLibrary") <<
"nval < 2 for array " << str <<
"\n";
619 throw cms::Exception(
"Unknown",
"HFShowerLibrary") <<
"cannot get array " << str <<
"\n";
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
T getParameter(std::string const &) const
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
void initRun(const HcalDDDSimConstants *)
Sin< T >::type sin(const T &t)
double attLength(double lambda)
HFShowerPhotonCollection pe
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Compact representation of the geometrical detector hierarchy.
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
static bool isStableHadron(int pdgCode)
std::vector< HFShowerPhoton > HFShowerPhotonCollection
const std::vector< double > & getRTableHF() const
Cos< T >::type cos(const T &t)
void loadEventInfo(TBranch *)
const std::vector< double > & getPhiTableHF() const
Abs< T >::type abs(const T &t)
std::vector< double > gpar
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &, int &)
HFShowerPhotonCollection * photo
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
const std::vector< double > & getGparHF() const
std::vector< double > pmom
void extrapolate(int, double)
HFShowerLibrary(const std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
static bool isGammaElectronPositron(int pdgCode)
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
std::vector< Hit > fillHits(const G4ThreeVector &p, const G4ThreeVector &v, int parCode, double parEnergy, bool &ok, double weight, double time, bool onlyLong=false)
HFShowerPhotonCollection photon
void interpolate(int, double)
void initRun(G4ParticleTable *, const HcalDDDSimConstants *)