CMS 3D CMS Logo

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

#include <HFShowerLibrary.h>

Classes

struct  Hit
 

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 DDCompactView &cpv, edm::ParameterSet const &p)
 
void initRun (G4ParticleTable *, const HcalDDDSimConstants *)
 
 ~HFShowerLibrary ()
 

Protected Member Functions

void extrapolate (int, double)
 
std::vector< double > getDDDArray (const std::string &, const DDsvalues_type &, int &)
 
void getRecord (int, int)
 
void interpolate (int, double)
 
void loadEventInfo (TBranch *)
 
bool rInside (double r)
 
void storePhoton (int j)
 

Private Attributes

bool applyFidCut
 
double backProb
 
double dphi
 
TBranch * emBranch
 
int evtPerBin
 
HFFibrefibre
 
std::vector< double > gpar
 
TBranch * hadBranch
 
TFile * hf
 
float libVers
 
float listVersion
 
bool newForm
 
int nMomBin
 
int npe
 
HFShowerPhotonCollection pe
 
HFShowerPhotonCollectionphoto
 
HFShowerPhotonCollection photon
 
std::vector< double > pmom
 
double probMax
 
double rMax
 
double rMin
 
int totEvents
 
bool v3version
 
bool verbose
 

Detailed Description

Definition at line 29 of file HFShowerLibrary.h.

Constructor & Destructor Documentation

HFShowerLibrary::HFShowerLibrary ( const std::string &  name,
const DDCompactView cpv,
edm::ParameterSet const &  p 
)

Definition at line 26 of file HFShowerLibrary.cc.

References applyFidCut, backProb, emBranch, event(), evtPerBin, Exception, fibre, edm::ParameterSet::getParameter(), GeV, hadBranch, hf, mps_fire::i, info(), libVers, listVersion, loadEventInfo(), newForm, nMomBin, photo, pmom, probMax, AlCaHLTBitMon_QueryRunRegistry::string, totEvents, and v3version.

27  : fibre(nullptr), hf(nullptr), emBranch(nullptr), hadBranch(nullptr), npe(0) {
28  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
29  probMax = m_HF.getParameter<double>("ProbMax");
30 
31  edm::ParameterSet m_HS = p.getParameter<edm::ParameterSet>("HFShowerLibrary");
32  edm::FileInPath fp = m_HS.getParameter<edm::FileInPath>("FileName");
33  std::string pTreeName = fp.fullPath();
34  backProb = m_HS.getParameter<double>("BackProbability");
35  std::string emName = m_HS.getParameter<std::string>("TreeEMID");
36  std::string hadName = m_HS.getParameter<std::string>("TreeHadID");
37  std::string branchEvInfo = m_HS.getUntrackedParameter<std::string>(
38  "BranchEvt", "HFShowerLibraryEventInfos_hfshowerlib_HFShowerLibraryEventInfo");
39  std::string branchPre = m_HS.getUntrackedParameter<std::string>("BranchPre", "HFShowerPhotons_hfshowerlib_");
40  std::string branchPost = m_HS.getUntrackedParameter<std::string>("BranchPost", "_R.obj");
41  verbose = m_HS.getUntrackedParameter<bool>("Verbosity", false);
42  applyFidCut = m_HS.getParameter<bool>("ApplyFiducialCut");
43 
44  if (pTreeName.find(".") == 0)
45  pTreeName.erase(0, 2);
46  const char* nTree = pTreeName.c_str();
47  hf = TFile::Open(nTree);
48 
49  if (!hf->IsOpen()) {
50  edm::LogError("HFShower") << "HFShowerLibrary: opening " << nTree << " failed";
51  throw cms::Exception("Unknown", "HFShowerLibrary") << "Opening of " << pTreeName << " fails\n";
52  } else {
53  edm::LogVerbatim("HFShower") << "HFShowerLibrary: opening " << nTree << " successfully";
54  }
55 
56  newForm = (branchEvInfo.empty());
57  TTree* event(nullptr);
58  if (newForm)
59  event = (TTree*)hf->Get("HFSimHits");
60  else
61  event = (TTree*)hf->Get("Events");
62  if (event) {
63  TBranch* evtInfo(nullptr);
64  if (!newForm) {
65  std::string info = branchEvInfo + branchPost;
66  evtInfo = event->GetBranch(info.c_str());
67  }
68  if (evtInfo || newForm) {
69  loadEventInfo(evtInfo);
70  } else {
71  edm::LogError("HFShower") << "HFShowerLibrary: HFShowerLibrayEventInfo"
72  << " Branch does not exist in Event";
73  throw cms::Exception("Unknown", "HFShowerLibrary") << "Event information absent\n";
74  }
75  } else {
76  edm::LogError("HFShower") << "HFShowerLibrary: Events Tree does not "
77  << "exist";
78  throw cms::Exception("Unknown", "HFShowerLibrary") << "Events tree absent\n";
79  }
80 
81  std::stringstream ss;
82  ss << "HFShowerLibrary: Library " << libVers << " ListVersion " << listVersion << " Events Total " << totEvents
83  << " and " << evtPerBin << " per bin\n";
84  ss << "HFShowerLibrary: Energies (GeV) with " << nMomBin << " bins\n";
85  for (int i = 0; i < nMomBin; ++i) {
86  if (i / 10 * 10 == i && i > 0) {
87  ss << "\n";
88  }
89  ss << " " << pmom[i] / CLHEP::GeV;
90  }
91  edm::LogVerbatim("HFShower") << ss.str();
92 
93  std::string nameBr = branchPre + emName + branchPost;
94  emBranch = event->GetBranch(nameBr.c_str());
95  if (verbose)
96  emBranch->Print();
97  nameBr = branchPre + hadName + branchPost;
98  hadBranch = event->GetBranch(nameBr.c_str());
99  if (verbose)
100  hadBranch->Print();
101 
102  v3version = false;
103  if (emBranch->GetClassName() == std::string("vector<float>")) {
104  v3version = true;
105  }
106 
107  edm::LogVerbatim("HFShower") << " HFShowerLibrary:Branch " << emName << " has " << emBranch->GetEntries()
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. "
112  << backProb;
113 
114  fibre = new HFFibre(name, cpv, p);
116 }
T getParameter(std::string const &) const
static const TGPicture * info(bool iBackgroundIsBlack)
const double GeV
Definition: MathUtil.h:16
std::vector< HFShowerPhoton > HFShowerPhotonCollection
void loadEventInfo(TBranch *)
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
TBranch * emBranch
std::vector< double > pmom
Definition: event.py:1
TBranch * hadBranch
HFShowerLibrary::~HFShowerLibrary ( )

Definition at line 118 of file HFShowerLibrary.cc.

References fibre, hf, and photo.

118  {
119  if (hf)
120  hf->Close();
121  delete fibre;
122  delete photo;
123 }
HFShowerPhotonCollection * photo

Member Function Documentation

void HFShowerLibrary::extrapolate ( int  type,
double  pin 
)
protected

Definition at line 514 of file HFShowerLibrary.cc.

References evtPerBin, getRecord(), createfilelist::int, newForm, nMomBin, npe, pe, photo, photon, pmom, alignCSCRings::r, storePhoton(), totEvents, and w.

Referenced by fillHits().

514  {
515  int nrec = int(pin / pmom[nMomBin - 1]);
516  double w = (pin - pmom[nMomBin - 1] * nrec) / pmom[nMomBin - 1];
517  nrec++;
518 #ifdef EDM_ML_DEBUG
519  edm::LogVerbatim("HFShower") << "HFShowerLibrary:: Extrapolate for Energy " << pin << " GeV with " << nMomBin
520  << " momentum bins and " << evtPerBin << " entries/bin -- "
521  << "total " << totEvents << " using " << nrec << " records";
522 #endif
523  std::vector<int> irc(nrec);
524 
525  for (int ir = 0; ir < nrec; ir++) {
526  double r = G4UniformRand();
527  irc[ir] = int(evtPerBin * 0.5 * r) + (nMomBin - 1) * evtPerBin + 1;
528  if (irc[ir] < 1) {
529  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[" << ir << "] = " << irc[ir] << " now set to 1";
530  irc[ir] = 1;
531  } else if (irc[ir] > totEvents) {
532  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[" << ir << "] = " << irc[ir] << " now set to "
533  << totEvents;
534  irc[ir] = totEvents;
535 #ifdef EDM_ML_DEBUG
536  } else {
537  edm::LogVerbatim("HFShower") << "HFShowerLibrary::Extrapolation use irc[" << ir << "] = " << irc[ir];
538 #endif
539  }
540  }
541 
542  pe.clear();
543  npe = 0;
544  int npold = 0;
545  for (int ir = 0; ir < nrec; ir++) {
546  if (irc[ir] > 0) {
547  getRecord(type, irc[ir]);
548  int nPhoton = (newForm) ? photo->size() : photon.size();
549  npold += nPhoton;
550  for (int j = 0; j < nPhoton; j++) {
551  double r = G4UniformRand();
552  if (ir != nrec - 1 || r < w) {
553  storePhoton(j);
554  }
555  }
556 #ifdef EDM_ML_DEBUG
557  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Record [" << ir << "] = " << irc[ir] << " npold = " << npold;
558 #endif
559  }
560  }
561 #ifdef EDM_ML_DEBUG
562  edm::LogVerbatim("HFShower") << "HFShowerLibrary:: uses " << npold << " photons";
563 #endif
564 
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
568  << " *****";
569 #ifdef EDM_ML_DEBUG
570  else
571  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Extrapolation == " << nrec << " records " << irc[0] << ", "
572  << irc[1] << ", ... gives a buffer of " << npold << " photons and fills " << npe
573  << " PE";
574  for (int j = 0; j < npe; j++)
575  edm::LogVerbatim("HFShower") << "Photon " << j << " " << pe[j];
576 #endif
577 }
type
Definition: HCALResponse.h:21
void getRecord(int, int)
const double w
Definition: UKUtility.cc:23
HFShowerPhotonCollection pe
HFShowerPhotonCollection * photo
std::vector< double > pmom
void storePhoton(int j)
HFShowerPhotonCollection photon
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 185 of file HFShowerLibrary.cc.

References funct::abs(), applyFidCut, HFFibre::attLength(), backProb, funct::cos(), egammaForCoreTracking_cff::depth, HFShowerLibrary::Hit::depth, dphi, JetChargeProducer_cfi::exp, extrapolate(), fibre, gpar, mps_fire::i, createfilelist::int, interpolate(), G4TrackToParticleID::isGammaElectronPositron(), G4TrackToParticleID::isStableHadron(), MeV, nMomBin, npe, AlCaHLTBitMon_ParallelJobs::p, pe, pmom, HFShowerLibrary::Hit::position, probMax, alignCSCRings::r, diffTwoXMLs::r1, diffTwoXMLs::r2, rInside(), funct::sin(), electronIdCutBased_cfi::threshold, HFShowerLibrary::Hit::time, HFFibre::tShift(), geometryCSVtoXML::xx, geometryCSVtoXML::yy, z, HFFibre::zShift(), and geometryCSVtoXML::zz.

Referenced by getHits().

192  {
193  std::vector<HFShowerLibrary::Hit> hit;
194  ok = false;
196  // shower is built only for gamma, e+- and stable hadrons
197  if (!isEM && !G4TrackToParticleID::isStableHadron(parCode)) {
198  return hit;
199  }
200  ok = true;
201 
202  // remove low-energy component
203  const double threshold = 50 * MeV;
204  if (pin < threshold) {
205  return hit;
206  }
207 
208  double pz = momDir.z();
209  double zint = hitPoint.z();
210 
211  // if particle moves from interaction point or "backwards (halo)
212  bool backward = (pz * zint < 0.) ? true : false;
213 
214  double sphi = sin(momDir.phi());
215  double cphi = cos(momDir.phi());
216  double ctheta = cos(momDir.theta());
217  double stheta = sin(momDir.theta());
218 
219  if (isEM) {
220  if (pin < pmom[nMomBin - 1]) {
221  interpolate(0, pin);
222  } else {
223  extrapolate(0, pin);
224  }
225  } else {
226  if (pin < pmom[nMomBin - 1]) {
227  interpolate(1, pin);
228  } else {
229  extrapolate(1, pin);
230  }
231  }
232 
233  int nHit = 0;
234  HFShowerLibrary::Hit oneHit;
235  for (int i = 0; i < npe; ++i) {
236  double zv = std::abs(pe[i].z()); // abs local z
237 #ifdef EDM_ML_DEBUG
238  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Hit " << i << " " << pe[i] << " zv " << zv;
239 #endif
240  if (zv <= gpar[1] && pe[i].lambda() > 0 && (pe[i].z() >= 0 || (zv > gpar[0] && (!onlyLong)))) {
241  int depth = 1;
242  if (onlyLong) {
243  } else if (!backward) { // fully valid only for "front" particles
244  if (pe[i].z() < 0)
245  depth = 2; // with "front"-simulated shower lib.
246  } else { // for "backward" particles - almost equal
247  double r = G4UniformRand(); // share between L and S fibers
248  if (r > 0.5)
249  depth = 2;
250  }
251 
252  // Updated coordinate transformation from local
253  // back to global using two Euler angles: phi and theta
254  double pex = pe[i].x();
255  double pey = pe[i].y();
256 
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;
260 
261  G4ThreeVector pos = hitPoint + G4ThreeVector(xx, yy, zz);
262  zv = std::abs(pos.z()) - gpar[4] - 0.5 * gpar[1];
263  G4ThreeVector lpos = G4ThreeVector(pos.x(), pos.y(), zv);
264 
265  zv = fibre->zShift(lpos, depth, 0); // distance to PMT !
266 
267  double r = pos.perp();
268  double p = fibre->attLength(pe[i].lambda());
269  double fi = pos.phi();
270  if (fi < 0)
271  fi += CLHEP::twopi;
272  int isect = int(fi / dphi) + 1;
273  isect = (isect + 1) / 2;
274  double dfi = ((isect * 2 - 1) * dphi - fi);
275  if (dfi < 0)
276  dfi = -dfi;
277  double dfir = r * sin(dfi);
278 #ifdef EDM_ML_DEBUG
279  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Position shift " << xx << ", " << yy << ", " << zz << ": "
280  << pos << " R " << r << " Phi " << fi << " Section " << isect << " R*Dfi " << dfir
281  << " Dist " << zv;
282 #endif
283  zz = std::abs(pos.z());
284  double r1 = G4UniformRand();
285  double r2 = G4UniformRand();
286  double r3 = backward ? G4UniformRand() : -9999.;
287  if (!applyFidCut)
288  dfir += gpar[5];
289 
290 #ifdef EDM_ML_DEBUG
291  edm::LogVerbatim("HFShower") << "HFShowerLibrary: rLimits " << rInside(r) << " attenuation " << r1 << ":"
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]);
299 #endif
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])) {
302  oneHit.position = pos;
303  oneHit.depth = depth;
304  oneHit.time = (tSlice + (pe[i].t()) + (fibre->tShift(lpos, depth, 1)));
305  hit.push_back(oneHit);
306 #ifdef EDM_ML_DEBUG
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);
310 #endif
311  ++nHit;
312  }
313 #ifdef EDM_ML_DEBUG
314  else
315  edm::LogVerbatim("HFShower") << "HFShowerLibrary: REJECTED !!!";
316 #endif
317  if (onlyLong && zz >= gpar[4] + gpar[0] && zz <= gpar[4] + gpar[1]) {
318  r1 = G4UniformRand();
319  r2 = G4UniformRand();
320  if (rInside(r) && r1 <= exp(-p * zv) && r2 <= probMax && dfir > gpar[5]) {
321  oneHit.position = pos;
322  oneHit.depth = 2;
323  oneHit.time = (tSlice + (pe[i].t()) + (fibre->tShift(lpos, 2, 1)));
324  hit.push_back(oneHit);
325 #ifdef EDM_ML_DEBUG
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);
329 #endif
330  ++nHit;
331  }
332  }
333  }
334  }
335 
336 #ifdef EDM_ML_DEBUG
337  edm::LogVerbatim("HFShower") << "HFShowerLibrary: Total Hits " << nHit << " out of " << npe << " PE";
338 #endif
339  if (nHit > npe && !onlyLong) {
340  edm::LogWarning("HFShower") << "HFShowerLibrary: Hit buffer " << npe << " smaller than " << nHit << " Hits";
341  }
342  return hit;
343 }
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:113
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double attLength(double lambda)
Definition: HFFibre.cc:97
HFShowerPhotonCollection pe
double zShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:123
Definition: weight.py:1
static bool isStableHadron(int pdgCode)
const double MeV
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > gpar
std::vector< double > pmom
G4ThreeVector position
void extrapolate(int, double)
static bool isGammaElectronPositron(int pdgCode)
bool rInside(double r)
void interpolate(int, double)
std::vector< double > HFShowerLibrary::getDDDArray ( const std::string &  str,
const DDsvalues_type sv,
int &  nmin 
)
protected

Definition at line 590 of file HFShowerLibrary.cc.

References DDfetch(), DDValue::doubles(), Exception, str, and relativeConstraints::value.

590  {
591 #ifdef EDM_ML_DEBUG
592  edm::LogVerbatim("HFShower") << "HFShowerLibrary:getDDDArray called for " << str << " with nMin " << nmin;
593 #endif
594  DDValue value(str);
595  if (DDfetch(&sv, value)) {
596 #ifdef EDM_ML_DEBUG
597  edm::LogVerbatim("HFShower") << value;
598 #endif
599  const std::vector<double>& fvec = value.doubles();
600  int nval = fvec.size();
601  if (nmin > 0) {
602  if (nval < nmin) {
603  edm::LogError("HFShower") << "HFShowerLibrary : # of " << str << " bins " << nval << " < " << nmin
604  << " ==> illegal";
605  throw cms::Exception("Unknown", "HFShowerLibrary") << "nval < nmin for array " << str << "\n";
606  }
607  } else {
608  if (nval < 2) {
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";
612  }
613  }
614  nmin = nval;
615 
616  return fvec;
617  } else {
618  edm::LogError("HFShower") << "HFShowerLibrary : cannot get array " << str;
619  throw cms::Exception("Unknown", "HFShowerLibrary") << "cannot get array " << str << "\n";
620  }
621 }
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
Definition: value.py:1
#define str(s)
std::vector< HFShowerLibrary::Hit > HFShowerLibrary::getHits ( const G4Step *  aStep,
bool &  ok,
double  weight,
bool  onlyLong = false 
)

Definition at line 144 of file HFShowerLibrary.cc.

References funct::cos(), fillHits(), GeV, gpar, funct::sin(), and HiIsolationCommonParameters_cff::track.

Referenced by HFShowerParam::getHits().

147  {
148  auto const preStepPoint = aStep->GetPreStepPoint();
149  auto const postStepPoint = aStep->GetPostStepPoint();
150  auto const track = aStep->GetTrack();
151  // Get Z-direction
152  auto const aParticle = track->GetDynamicParticle();
153  const G4ThreeVector& momDir = aParticle->GetMomentumDirection();
154  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
155  int parCode = track->GetDefinition()->GetPDGEncoding();
156 
157  // VI: for ions use internally pdg code of alpha in order to keep
158  // consistency with previous simulation
159  if (track->GetDefinition()->IsGeneralIon()) {
160  parCode = 1000020040;
161  }
162 
163 #ifdef EDM_ML_DEBUG
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];
167 
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());
174 #endif
175 
176  double tSlice = (postStepPoint->GetGlobalTime()) / CLHEP::nanosecond;
177 
178  // use kinetic energy for protons and ions
179  double pin = (track->GetDefinition()->GetBaryonNumber() > 0) ? preStepPoint->GetKineticEnergy()
180  : preStepPoint->GetTotalEnergy();
181 
182  return fillHits(hitPoint, momDir, parCode, pin, isKilled, weight, tSlice, onlyLong);
183 }
const double GeV
Definition: MathUtil.h:16
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Definition: weight.py:1
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< double > gpar
std::vector< Hit > fillHits(const G4ThreeVector &p, const G4ThreeVector &v, int parCode, double parEnergy, bool &ok, double weight, double time, bool onlyLong=false)
void HFShowerLibrary::getRecord ( int  type,
int  record 
)
protected

Definition at line 347 of file HFShowerLibrary.cc.

References emBranch, hadBranch, mps_fire::i, newForm, photo, photon, lumiQTWidget::t, totEvents, and v3version.

Referenced by extrapolate(), and interpolate().

347  {
348  int nrc = record - 1;
349  photon.clear();
350  photo->clear();
351  if (type > 0) {
352  if (newForm) {
353  if (!v3version) {
354  hadBranch->SetAddress(&photo);
355  hadBranch->GetEntry(nrc + totEvents);
356  } else {
357  std::vector<float> t;
358  std::vector<float>* tp = &t;
359  hadBranch->SetAddress(&tp);
360  hadBranch->GetEntry(nrc + totEvents);
361  unsigned int tSize = t.size() / 5;
362  photo->reserve(tSize);
363  for (unsigned int i = 0; i < tSize; i++) {
364  photo->push_back(
365  HFShowerPhoton(t[i], t[1 * tSize + i], t[2 * tSize + i], t[3 * tSize + i], t[4 * tSize + i]));
366  }
367  }
368  } else {
369  hadBranch->SetAddress(&photon);
370  hadBranch->GetEntry(nrc);
371  }
372  } else {
373  if (newForm) {
374  if (!v3version) {
375  emBranch->SetAddress(&photo);
376  emBranch->GetEntry(nrc);
377  } else {
378  std::vector<float> t;
379  std::vector<float>* tp = &t;
380  emBranch->SetAddress(&tp);
381  emBranch->GetEntry(nrc);
382  unsigned int tSize = t.size() / 5;
383  photo->reserve(tSize);
384  for (unsigned int i = 0; i < tSize; i++) {
385  photo->push_back(
386  HFShowerPhoton(t[i], t[1 * tSize + i], t[2 * tSize + i], t[3 * tSize + i], t[4 * tSize + i]));
387  }
388  }
389  } else {
390  emBranch->SetAddress(&photon);
391  emBranch->GetEntry(nrc);
392  }
393  }
394 #ifdef EDM_ML_DEBUG
395  int nPhoton = (newForm) ? photo->size() : photon.size();
396  edm::LogVerbatim("HFShower") << "HFShowerLibrary::getRecord: Record " << record << " of type " << type << " with "
397  << nPhoton << " photons";
398  for (int j = 0; j < nPhoton; j++)
399  if (newForm)
400  edm::LogVerbatim("HFShower") << "Photon " << j << " " << photo->at(j);
401  else
402  edm::LogVerbatim("HFShower") << "Photon " << j << " " << photon[j];
403 #endif
404 }
type
Definition: HCALResponse.h:21
JetCorrectorParameters::Record record
Definition: classes.h:7
HFShowerPhotonCollection * photo
TBranch * emBranch
HFShowerPhotonCollection photon
TBranch * hadBranch
void HFShowerLibrary::initRun ( G4ParticleTable *  ,
const HcalDDDSimConstants hcons 
)

Definition at line 125 of file HFShowerLibrary.cc.

References dphi, fibre, HcalDDDSimConstants::getGparHF(), HcalDDDSimConstants::getPhiTableHF(), HcalDDDSimConstants::getRTableHF(), gpar, HFFibre::initRun(), rMax, and rMin.

Referenced by HFShowerParam::initRun().

125  {
126  if (fibre)
127  fibre->initRun(hcons);
128 
129  //Radius (minimum and maximum)
130  std::vector<double> rTable = hcons->getRTableHF();
131  rMin = rTable[0];
132  rMax = rTable[rTable.size() - 1];
133 
134  //Delta phi
135  std::vector<double> phibin = hcons->getPhiTableHF();
136  dphi = phibin[0];
137  edm::LogVerbatim("HFShower") << "HFShowerLibrary: rMIN " << rMin / cm << " cm and rMax " << rMax / cm
138  << " (Half) Phi Width of wedge " << dphi / deg;
139 
140  //Special Geometry parameters
141  gpar = hcons->getGparHF();
142 }
void initRun(const HcalDDDSimConstants *)
Definition: HFFibre.cc:81
const std::vector< double > & getRTableHF() const
const std::vector< double > & getPhiTableHF() const
std::vector< double > gpar
const std::vector< double > & getGparHF() const
void HFShowerLibrary::interpolate ( int  type,
double  pin 
)
protected

Definition at line 437 of file HFShowerLibrary.cc.

References evtPerBin, getRecord(), GeV, createfilelist::int, newForm, nMomBin, npe, pe, photo, photon, pmom, alignCSCRings::r, storePhoton(), totEvents, and w.

Referenced by fillHits().

437  {
438 #ifdef EDM_ML_DEBUG
439  edm::LogVerbatim("HFShower") << "HFShowerLibrary:: Interpolate for Energy " << pin / GeV << " GeV with " << nMomBin
440  << " momentum bins and " << evtPerBin << " entries/bin -- total " << totEvents;
441 #endif
442  int irc[2] = {0, 0};
443  double w = 0.;
444  double r = G4UniformRand();
445 
446  if (pin < pmom[0]) {
447  w = pin / pmom[0];
448  irc[1] = int(evtPerBin * r) + 1;
449  irc[0] = 0;
450  } else {
451  for (int j = 0; j < nMomBin - 1; j++) {
452  if (pin >= pmom[j] && pin < pmom[j + 1]) {
453  w = (pin - pmom[j]) / (pmom[j + 1] - pmom[j]);
454  if (j == nMomBin - 2) {
455  irc[1] = int(evtPerBin * 0.5 * r);
456  } else {
457  irc[1] = int(evtPerBin * r);
458  }
459  irc[1] += (j + 1) * evtPerBin + 1;
460  r = G4UniformRand();
461  irc[0] = int(evtPerBin * r) + 1 + j * evtPerBin;
462  if (irc[0] < 0) {
463  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[0] = " << irc[0] << " now set to 0";
464  irc[0] = 0;
465  } else if (irc[0] > totEvents) {
466  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[0] = " << irc[0] << " now set to " << totEvents;
467  irc[0] = totEvents;
468  }
469  }
470  }
471  }
472  if (irc[1] < 1) {
473  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[1] = " << irc[1] << " now set to 1";
474  irc[1] = 1;
475  } else if (irc[1] > totEvents) {
476  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[1] = " << irc[1] << " now set to " << totEvents;
477  irc[1] = totEvents;
478  }
479 
480 #ifdef EDM_ML_DEBUG
481  edm::LogVerbatim("HFShower") << "HFShowerLibrary:: Select records " << irc[0] << " and " << irc[1] << " with weights "
482  << 1 - w << " and " << w;
483 #endif
484  pe.clear();
485  npe = 0;
486  int npold = 0;
487  for (int ir = 0; ir < 2; ir++) {
488  if (irc[ir] > 0) {
489  getRecord(type, irc[ir]);
490  int nPhoton = (newForm) ? photo->size() : photon.size();
491  npold += nPhoton;
492  for (int j = 0; j < nPhoton; j++) {
493  r = G4UniformRand();
494  if ((ir == 0 && r > w) || (ir > 0 && r < w)) {
495  storePhoton(j);
496  }
497  }
498  }
499  }
500 
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 << " *****";
505 #ifdef EDM_ML_DEBUG
506  else
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++)
510  edm::LogVerbatim("HFShower") << "Photon " << j << " " << pe[j];
511 #endif
512 }
type
Definition: HCALResponse.h:21
const double GeV
Definition: MathUtil.h:16
void getRecord(int, int)
const double w
Definition: UKUtility.cc:23
HFShowerPhotonCollection pe
HFShowerPhotonCollection * photo
std::vector< double > pmom
void storePhoton(int j)
HFShowerPhotonCollection photon
void HFShowerLibrary::loadEventInfo ( TBranch *  branch)
protected

Definition at line 406 of file HFShowerLibrary.cc.

References evtPerBin, GeV, mps_fire::i, libVers, listVersion, nMomBin, pmom, and totEvents.

Referenced by HFShowerLibrary().

406  {
407  if (branch) {
408  std::vector<HFShowerLibraryEventInfo> eventInfoCollection;
409  branch->SetAddress(&eventInfoCollection);
410  branch->GetEntry(0);
411 #ifdef EDM_ML_DEBUG
412  edm::LogVerbatim("HFShower") << "HFShowerLibrary::loadEventInfo loads EventInfo Collection of size "
413  << eventInfoCollection.size() << " records";
414 #endif
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();
421  } else {
422 #ifdef EDM_ML_DEBUG
423  edm::LogVerbatim("HFShower") << "HFShowerLibrary::loadEventInfo loads EventInfo from hardwired"
424  << " numbers";
425 #endif
426  nMomBin = 16;
427  evtPerBin = 5000;
429  libVers = 1.1;
430  listVersion = 3.6;
431  pmom = {2, 3, 5, 7, 10, 15, 20, 30, 50, 75, 100, 150, 250, 350, 500, 1000};
432  }
433  for (int i = 0; i < nMomBin; i++)
434  pmom[i] *= GeV;
435 }
const double GeV
Definition: MathUtil.h:16
std::vector< double > pmom
bool HFShowerLibrary::rInside ( double  r)
protected

Definition at line 345 of file HFShowerLibrary.cc.

References rMax, and rMin.

Referenced by fillHits().

345 { return (r >= rMin && r <= rMax); }
void HFShowerLibrary::storePhoton ( int  j)
protected

Definition at line 579 of file HFShowerLibrary.cc.

References newForm, npe, pe, photo, and photon.

Referenced by extrapolate(), and interpolate().

579  {
580  if (newForm)
581  pe.push_back(photo->at(j));
582  else
583  pe.push_back(photon[j]);
584 #ifdef EDM_ML_DEBUG
585  edm::LogVerbatim("HFShower") << "HFShowerLibrary: storePhoton " << j << " npe " << npe << " " << pe[npe];
586 #endif
587  npe++;
588 }
HFShowerPhotonCollection pe
HFShowerPhotonCollection * photo
HFShowerPhotonCollection photon

Member Data Documentation

bool HFShowerLibrary::applyFidCut
private

Definition at line 68 of file HFShowerLibrary.h.

Referenced by fillHits(), and HFShowerLibrary().

double HFShowerLibrary::backProb
private

Definition at line 73 of file HFShowerLibrary.h.

Referenced by fillHits(), and HFShowerLibrary().

double HFShowerLibrary::dphi
private

Definition at line 74 of file HFShowerLibrary.h.

Referenced by fillHits(), and initRun().

TBranch* HFShowerLibrary::emBranch
private

Definition at line 66 of file HFShowerLibrary.h.

Referenced by getRecord(), and HFShowerLibrary().

int HFShowerLibrary::evtPerBin
private

Definition at line 69 of file HFShowerLibrary.h.

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

HFFibre* HFShowerLibrary::fibre
private

Definition at line 64 of file HFShowerLibrary.h.

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

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

Definition at line 75 of file HFShowerLibrary.h.

Referenced by fillHits(), getHits(), and initRun().

TBranch * HFShowerLibrary::hadBranch
private

Definition at line 66 of file HFShowerLibrary.h.

Referenced by getRecord(), and HFShowerLibrary().

TFile* HFShowerLibrary::hf
private

Definition at line 65 of file HFShowerLibrary.h.

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

float HFShowerLibrary::libVers
private

Definition at line 70 of file HFShowerLibrary.h.

Referenced by HFShowerLibrary(), and loadEventInfo().

float HFShowerLibrary::listVersion
private

Definition at line 70 of file HFShowerLibrary.h.

Referenced by HFShowerLibrary(), and loadEventInfo().

bool HFShowerLibrary::newForm
private

Definition at line 68 of file HFShowerLibrary.h.

Referenced by extrapolate(), getRecord(), HFShowerLibrary(), interpolate(), and storePhoton().

int HFShowerLibrary::nMomBin
private

Definition at line 69 of file HFShowerLibrary.h.

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

int HFShowerLibrary::npe
private

Definition at line 77 of file HFShowerLibrary.h.

Referenced by extrapolate(), fillHits(), interpolate(), and storePhoton().

HFShowerPhotonCollection HFShowerLibrary::pe
private

Definition at line 78 of file HFShowerLibrary.h.

Referenced by extrapolate(), fillHits(), interpolate(), and storePhoton().

HFShowerPhotonCollection* HFShowerLibrary::photo
private
HFShowerPhotonCollection HFShowerLibrary::photon
private

Definition at line 80 of file HFShowerLibrary.h.

Referenced by extrapolate(), getRecord(), interpolate(), and storePhoton().

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

Definition at line 71 of file HFShowerLibrary.h.

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

double HFShowerLibrary::probMax
private

Definition at line 73 of file HFShowerLibrary.h.

Referenced by fillHits(), and HFShowerLibrary().

double HFShowerLibrary::rMax
private

Definition at line 74 of file HFShowerLibrary.h.

Referenced by initRun(), and rInside().

double HFShowerLibrary::rMin
private

Definition at line 74 of file HFShowerLibrary.h.

Referenced by initRun(), and rInside().

int HFShowerLibrary::totEvents
private

Definition at line 69 of file HFShowerLibrary.h.

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

bool HFShowerLibrary::v3version
private

Definition at line 68 of file HFShowerLibrary.h.

Referenced by getRecord(), and HFShowerLibrary().

bool HFShowerLibrary::verbose
private

Definition at line 68 of file HFShowerLibrary.h.