CMS 3D CMS Logo

HFShowerLibrary.cc
Go to the documentation of this file.
1 // File: HFShowerLibrary.cc
3 // Description: Shower library for Very forward hadron calorimeter
5 
12 
14 
15 #include "G4VPhysicalVolume.hh"
16 #include "G4NavigationHistory.hh"
17 #include "G4Step.hh"
18 #include "G4Track.hh"
19 #include "G4ParticleTable.hh"
20 #include "Randomize.hh"
21 #include "CLHEP/Units/SystemOfUnits.h"
22 #include "CLHEP/Units/PhysicalConstants.h"
23 
24 //#define DebugLog
25 
27  edm::ParameterSet const & p) : fibre(nullptr),hf(nullptr),
28  emBranch(nullptr),
29  hadBranch(nullptr),
30  npe(0) {
31 
32 
33  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
34  probMax = m_HF.getParameter<double>("ProbMax");
35 
36  edm::ParameterSet m_HS= p.getParameter<edm::ParameterSet>("HFShowerLibrary");
37  edm::FileInPath fp = m_HS.getParameter<edm::FileInPath>("FileName");
38  std::string pTreeName = fp.fullPath();
39  backProb = m_HS.getParameter<double>("BackProbability");
40  std::string emName = m_HS.getParameter<std::string>("TreeEMID");
41  std::string hadName = m_HS.getParameter<std::string>("TreeHadID");
42  std::string branchEvInfo = m_HS.getUntrackedParameter<std::string>("BranchEvt","HFShowerLibraryEventInfos_hfshowerlib_HFShowerLibraryEventInfo");
43  std::string branchPre = m_HS.getUntrackedParameter<std::string>("BranchPre","HFShowerPhotons_hfshowerlib_");
44  std::string branchPost = m_HS.getUntrackedParameter<std::string>("BranchPost","_R.obj");
45  verbose = m_HS.getUntrackedParameter<bool>("Verbosity",false);
46  applyFidCut = m_HS.getParameter<bool>("ApplyFiducialCut");
47 
48  if (pTreeName.find(".") == 0) pTreeName.erase(0,2);
49  const char* nTree = pTreeName.c_str();
50  hf = TFile::Open(nTree);
51 
52  if (!hf->IsOpen()) {
53  edm::LogError("HFShower") << "HFShowerLibrary: opening " << nTree
54  << " failed";
55  throw cms::Exception("Unknown", "HFShowerLibrary")
56  << "Opening of " << pTreeName << " fails\n";
57  } else {
58  edm::LogInfo("HFShower") << "HFShowerLibrary: opening " << nTree
59  << " successfully";
60  }
61 
62  newForm = (branchEvInfo.empty());
63  TTree* event(nullptr);
64  if (newForm) event = (TTree *) hf ->Get("HFSimHits");
65  else event = (TTree *) hf ->Get("Events");
66  if (event) {
67  TBranch *evtInfo(nullptr);
68  if (!newForm) {
69  std::string info = branchEvInfo + branchPost;
70  evtInfo = event->GetBranch(info.c_str());
71  }
72  if (evtInfo || newForm) {
73  loadEventInfo(evtInfo);
74  } else {
75  edm::LogError("HFShower") << "HFShowerLibrary: HFShowerLibrayEventInfo"
76  << " Branch does not exist in Event";
77  throw cms::Exception("Unknown", "HFShowerLibrary")
78  << "Event information absent\n";
79  }
80  } else {
81  edm::LogError("HFShower") << "HFShowerLibrary: Events Tree does not "
82  << "exist";
83  throw cms::Exception("Unknown", "HFShowerLibrary")
84  << "Events tree absent\n";
85  }
86 
87  std::stringstream ss;
88  ss << "HFShowerLibrary: Library " << libVers << " ListVersion " << listVersion
89  << " Events Total " << totEvents << " and " << evtPerBin << " per bin\n";
90  ss << "HFShowerLibrary: Energies (GeV) with " << nMomBin << " bins\n";
91  for (int i=0; i<nMomBin; ++i) {
92  if(i/10*10 == i && i > 0) { ss << "\n"; }
93  ss << " " << pmom[i]/CLHEP::GeV;
94  }
95  edm::LogInfo("HFShower") << ss.str();
96 
97  std::string nameBr = branchPre + emName + branchPost;
98  emBranch = event->GetBranch(nameBr.c_str());
99  if (verbose) emBranch->Print();
100  nameBr = branchPre + hadName + branchPost;
101  hadBranch = event->GetBranch(nameBr.c_str());
102  if (verbose) hadBranch->Print();
103 
104  v3version=false;
105  if ( emBranch->GetClassName() == std::string("vector<float>") ) {
106  v3version=true;
107  }
108 
109  edm::LogInfo("HFShower") << " HFShowerLibrary:Branch " << emName
110  << " has " << emBranch->GetEntries()
111  << " entries and Branch " << hadName
112  << " has " << hadBranch->GetEntries()
113  << " entries"
114  << "\n HFShowerLibrary::No packing information -"
115  << " Assume x, y, z are not in packed form"
116  << "\n Maximum probability cut off "
117  << probMax << " Back propagation of light prob. "
118  << backProb;
119 
120  fibre = new HFFibre(name, cpv, p);
122 }
123 
125  if (hf) hf->Close();
126  delete fibre;
127  delete photo;
128 }
129 
130 void HFShowerLibrary::initRun(G4ParticleTable*, const HcalDDDSimConstants* hcons) {
131 
132  if (fibre) fibre->initRun(hcons);
133 
134  //Radius (minimum and maximum)
135  std::vector<double> rTable = hcons->getRTableHF();
136  rMin = rTable[0];
137  rMax = rTable[rTable.size()-1];
138 
139  //Delta phi
140  std::vector<double> phibin = hcons->getPhiTableHF();
141  dphi = phibin[0];
142  edm::LogInfo("HFShower") << "HFShowerLibrary: rMIN " << rMin/cm
143  << " cm and rMax " << rMax/cm
144  << " (Half) Phi Width of wedge "
145  << dphi/deg;
146 
147  //Special Geometry parameters
148  gpar = hcons->getGparHF();
149 }
150 
151 std::vector<HFShowerLibrary::Hit> HFShowerLibrary::getHits(const G4Step * aStep,
152  bool & isKilled,
153  double weight,
154  bool onlyLong) {
155 
156  auto const preStepPoint = aStep->GetPreStepPoint();
157  auto const postStepPoint = aStep->GetPostStepPoint();
158  auto const track = aStep->GetTrack();
159  // Get Z-direction
160  auto const aParticle = track->GetDynamicParticle();
161  const G4ThreeVector& momDir = aParticle->GetMomentumDirection();
162  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
163  int parCode = track->GetDefinition()->GetPDGEncoding();
164 
165  // VI: for ions use internally pdg code of alpha in order to keep
166  // consistency with previous simulation
167  if(track->GetDefinition()->IsGeneralIon()) { parCode = 1000020040; }
168 
169 #ifdef DebugLog
170  G4String partType = track->GetDefinition()->GetParticleName();
171  const G4ThreeVector localPos =
172  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
173  double zoff = localPos.z() + 0.5*gpar[1];
174 
175  edm::LogInfo("HFShower") << "HFShowerLibrary::getHits " << partType
176  << " of energy " << pin/GeV << " GeV"
177  << " weight= " << weight << " onlyLong: " << onlyLong
178  << " dir.orts " << momDir.x() << ", " <<momDir.y()
179  << ", " << momDir.z() << " Pos x,y,z = "
180  << hitPoint.x() << "," << hitPoint.y() << ","
181  << hitPoint.z() << " (" << zoff
182  << ") sphi,cphi,stheta,ctheta = " << sin(momDir.phi())
183  << "," << cos(momDir.phi()) << ", " << sin(momDir.theta())
184  << "," << cos(momDir.theta());
185 #endif
186 
187  double tSlice = (postStepPoint->GetGlobalTime())/CLHEP::nanosecond;
188 
189  // use kinetic energy for protons and ions
190  double pin = (track->GetDefinition()->GetBaryonNumber() > 0)
191  ? preStepPoint->GetKineticEnergy() : preStepPoint->GetTotalEnergy();
192 
193  return fillHits(hitPoint,momDir,parCode,pin,isKilled,weight,tSlice,onlyLong);
194 }
195 
196 std::vector<HFShowerLibrary::Hit> HFShowerLibrary::fillHits(const G4ThreeVector & hitPoint,
197  const G4ThreeVector & momDir,
198  int parCode, double pin, bool & ok,
199  double weight, double tSlice,bool onlyLong) {
200 
201  std::vector<HFShowerLibrary::Hit> hit;
202  ok = false;
204  // shower is built only for gamma, e+- and stable hadrons
205  if (!isEM && !G4TrackToParticleID::isStableHadron(parCode)) {
206  return hit;
207  }
208  ok = true;
209 
210  // remove low-energy component
211  const double threshold = 50*MeV;
212  if(pin < threshold) { return hit; }
213 
214  double pz = momDir.z();
215  double zint = hitPoint.z();
216 
217  // if particle moves from interaction point or "backwards (halo)
218  bool backward = (pz * zint < 0.) ? true : false;
219 
220  double sphi = sin(momDir.phi());
221  double cphi = cos(momDir.phi());
222  double ctheta = cos(momDir.theta());
223  double stheta = sin(momDir.theta());
224 
225  if(isEM) {
226  if (pin<pmom[nMomBin-1]) {
227  interpolate(0, pin);
228  } else {
229  extrapolate(0, pin);
230  }
231  } else {
232  if (pin<pmom[nMomBin-1]) {
233  interpolate(1, pin);
234  } else {
235  extrapolate(1, pin);
236  }
237  }
238 
239  int nHit = 0;
240  HFShowerLibrary::Hit oneHit;
241  for (int i = 0; i < npe; ++i) {
242  double zv = std::abs(pe[i].z()); // abs local z
243 #ifdef DebugLog
244  edm::LogInfo("HFShower") << "HFShowerLibrary: Hit " << i << " " << pe[i] << " zv " << zv;
245 #endif
246  if (zv <= gpar[1] && pe[i].lambda() > 0 &&
247  (pe[i].z() >= 0 || (zv > gpar[0] && (!onlyLong)))) {
248  int depth = 1;
249  if (onlyLong) {
250  } else if (!backward) { // fully valid only for "front" particles
251  if (pe[i].z() < 0) depth = 2;// with "front"-simulated shower lib.
252  } else { // for "backward" particles - almost equal
253  double r = G4UniformRand(); // share between L and S fibers
254  if (r > 0.5) depth = 2;
255  }
256 
257 
258  // Updated coordinate transformation from local
259  // back to global using two Euler angles: phi and theta
260  double pex = pe[i].x();
261  double pey = pe[i].y();
262 
263  double xx = pex*ctheta*cphi - pey*sphi + zv*stheta*cphi;
264  double yy = pex*ctheta*sphi + pey*cphi + zv*stheta*sphi;
265  double zz = -pex*stheta + zv*ctheta;
266 
267  G4ThreeVector pos = hitPoint + G4ThreeVector(xx,yy,zz);
268  zv = std::abs(pos.z()) - gpar[4] - 0.5*gpar[1];
269  G4ThreeVector lpos = G4ThreeVector(pos.x(),pos.y(),zv);
270 
271  zv = fibre->zShift(lpos,depth,0); // distance to PMT !
272 
273  double r = pos.perp();
274  double p = fibre->attLength(pe[i].lambda());
275  double fi = pos.phi();
276  if (fi < 0) fi += CLHEP::twopi;
277  int isect = int(fi/dphi) + 1;
278  isect = (isect + 1) / 2;
279  double dfi = ((isect*2-1)*dphi - fi);
280  if (dfi < 0) dfi = -dfi;
281  double dfir = r * sin(dfi);
282 #ifdef DebugLog
283  edm::LogInfo("HFShower") << "HFShowerLibrary: Position shift " << xx
284  << ", " << yy << ", " << zz << ": " << pos
285  << " R " << r << " Phi " << fi << " Section "
286  << isect << " R*Dfi " << dfir << " Dist " << zv;
287 #endif
288  zz = std::abs(pos.z());
289  double r1 = G4UniformRand();
290  double r2 = G4UniformRand();
291  double r3 = backward ? G4UniformRand() : -9999.;
292  if (!applyFidCut) dfir += gpar[5];
293 
294 #ifdef DebugLog
295  edm::LogInfo("HFShower") << "HFShowerLibrary: rLimits " << rInside(r)
296  << " attenuation " << r1 <<":" << exp(-p*zv)
297  << " r2 " << r2 << " r3 " << r3 << " rDfi "
298  << gpar[5] << " zz "
299  << zz << " zLim " << gpar[4] << ":"
300  << gpar[4]+gpar[1] << "\n"
301  << " rInside(r) :" << rInside(r)
302  << " r1 <= exp(-p*zv) :" << (r1 <= exp(-p*zv))
303  << " r2 <= probMax :" << (r2 <= probMax*weight)
304  << " r3 <= backProb :" << (r3 <= backProb)
305  << " dfir > gpar[5] :" << (dfir > gpar[5])
306  << " zz >= gpar[4] :" << (zz >= gpar[4])
307  << " zz <= gpar[4]+gpar[1] :"
308  << (zz <= gpar[4]+gpar[1]);
309 #endif
310  if (rInside(r) && r1 <= exp(-p*zv) && r2 <= probMax*weight &&
311  dfir > gpar[5] && zz >= gpar[4] && zz <= gpar[4]+gpar[1] &&
312  r3 <= backProb && (depth != 2 || zz >= gpar[4]+gpar[0])) {
313  oneHit.position = pos;
314  oneHit.depth = depth;
315  oneHit.time = (tSlice+(pe[i].t())+(fibre->tShift(lpos,depth,1)));
316  hit.push_back(oneHit);
317 #ifdef DebugLog
318  edm::LogInfo("HFShower") << "HFShowerLibrary: Final Hit " << nHit
319  <<" position " << (hit[nHit].position)
320  << " Depth " << (hit[nHit].depth) <<" Time "
321  << tSlice << ":" << pe[i].t() << ":"
322  << fibre->tShift(lpos,depth,1) << ":"
323  << (hit[nHit].time);
324 #endif
325  ++nHit;
326  }
327 #ifdef DebugLog
328  else LogDebug("HFShower") << "HFShowerLibrary: REJECTED !!!";
329 #endif
330  if (onlyLong && zz >= gpar[4]+gpar[0] && zz <= gpar[4]+gpar[1]) {
331  r1 = G4UniformRand();
332  r2 = G4UniformRand();
333  if (rInside(r) && r1 <= exp(-p*zv) && r2 <= probMax && dfir > gpar[5]){
334  oneHit.position = pos;
335  oneHit.depth = 2;
336  oneHit.time = (tSlice+(pe[i].t())+(fibre->tShift(lpos,2,1)));
337  hit.push_back(oneHit);
338 #ifdef DebugLog
339  edm::LogInfo("HFShower") << "HFShowerLibrary: Final Hit " << nHit
340  << " position " << (hit[nHit].position)
341  << " Depth " << (hit[nHit].depth) <<" Time "
342  << tSlice << ":" << pe[i].t() << ":"
343  << fibre->tShift(lpos,2,1) << ":"
344  << (hit[nHit].time);
345 #endif
346  ++nHit;
347  }
348  }
349  }
350  }
351 
352 #ifdef DebugLog
353  edm::LogInfo("HFShower") << "HFShowerLibrary: Total Hits " << nHit
354  << " out of " << npe << " PE";
355 #endif
356  if (nHit > npe && !onlyLong) {
357  edm::LogWarning("HFShower") << "HFShowerLibrary: Hit buffer " << npe
358  << " smaller than " << nHit << " Hits";
359  }
360  return hit;
361 }
362 
364 
365  return (r >= rMin && r <= rMax);
366 }
367 
369 
370  int nrc = record-1;
371  photon.clear();
372  photo->clear();
373  if (type > 0) {
374  if (newForm) {
375  if ( !v3version ) {
376  hadBranch->SetAddress(&photo);
377  hadBranch->GetEntry(nrc+totEvents);
378  }
379  else{
380  std::vector<float> t;
381  std::vector<float> *tp=&t;
382  hadBranch->SetAddress(&tp);
383  hadBranch->GetEntry(nrc+totEvents);
384  unsigned int tSize=t.size()/5;
385  photo->reserve(tSize);
386  for ( unsigned int i=0; i<tSize; i++ ) {
387  photo->push_back( HFShowerPhoton( t[i], t[1*tSize+i], t[2*tSize+i], t[3*tSize+i], t[4*tSize+i] ) );
388  }
389  }
390  } else {
391  hadBranch->SetAddress(&photon);
392  hadBranch->GetEntry(nrc);
393  }
394  } else {
395  if (newForm) {
396  if (!v3version) {
397  emBranch->SetAddress(&photo);
398  emBranch->GetEntry(nrc);
399  }
400  else{
401  std::vector<float> t;
402  std::vector<float> *tp=&t;
403  emBranch->SetAddress(&tp);
404  emBranch->GetEntry(nrc);
405  unsigned int tSize=t.size()/5;
406  photo->reserve(tSize);
407  for ( unsigned int i=0; i<tSize; i++ ) {
408  photo->push_back( HFShowerPhoton( t[i], t[1*tSize+i], t[2*tSize+i], t[3*tSize+i], t[4*tSize+i] ) );
409  }
410  }
411  } else {
412  emBranch->SetAddress(&photon);
413  emBranch->GetEntry(nrc);
414  }
415  }
416 #ifdef DebugLog
417  int nPhoton = (newForm) ? photo->size() : photon.size();
418  LogDebug("HFShower") << "HFShowerLibrary::getRecord: Record " << record
419  << " of type " << type << " with " << nPhoton
420  << " photons";
421  for (int j = 0; j < nPhoton; j++)
422  if (newForm) LogDebug("HFShower") << "Photon " << j << " " << photo->at(j);
423  else LogDebug("HFShower") << "Photon " << j << " " << photon[j];
424 #endif
425 }
426 
428 
429  if (branch) {
430  std::vector<HFShowerLibraryEventInfo> eventInfoCollection;
431  branch->SetAddress(&eventInfoCollection);
432  branch->GetEntry(0);
433  edm::LogInfo("HFShower") << "HFShowerLibrary::loadEventInfo loads "
434  << " EventInfo Collection of size "
435  << eventInfoCollection.size() << " records";
436  totEvents = eventInfoCollection[0].totalEvents();
437  nMomBin = eventInfoCollection[0].numberOfBins();
438  evtPerBin = eventInfoCollection[0].eventsPerBin();
439  libVers = eventInfoCollection[0].showerLibraryVersion();
440  listVersion = eventInfoCollection[0].physListVersion();
441  pmom = eventInfoCollection[0].energyBins();
442  } else {
443  edm::LogInfo("HFShower") << "HFShowerLibrary::loadEventInfo loads "
444  << " EventInfo from hardwired numbers";
445  nMomBin = 16;
446  evtPerBin = 5000;
448  libVers = 1.1;
449  listVersion = 3.6;
450  pmom = {2,3,5,7,10,15,20,30,50,75,100,150,250,350,500,1000};
451  }
452  for (int i=0; i<nMomBin; i++)
453  pmom[i] *= GeV;
454 }
455 
456 void HFShowerLibrary::interpolate(int type, double pin) {
457 
458 #ifdef DebugLog
459  LogDebug("HFShower") << "HFShowerLibrary:: Interpolate for Energy " <<pin/GeV
460  << " GeV with " << nMomBin << " momentum bins and "
461  << evtPerBin << " entries/bin -- total " << totEvents;
462 #endif
463  int irc[2]={0,0};
464  double w = 0.;
465  double r = G4UniformRand();
466 
467  if (pin<pmom[0]) {
468  w = pin/pmom[0];
469  irc[1] = int(evtPerBin*r) + 1;
470  irc[0] = 0;
471  } else {
472  for (int j=0; j<nMomBin-1; j++) {
473  if (pin >= pmom[j] && pin < pmom[j+1]) {
474  w = (pin-pmom[j])/(pmom[j+1]-pmom[j]);
475  if (j == nMomBin-2) {
476  irc[1] = int(evtPerBin*0.5*r);
477  } else {
478  irc[1] = int(evtPerBin*r);
479  }
480  irc[1] += (j+1)*evtPerBin + 1;
481  r = G4UniformRand();
482  irc[0] = int(evtPerBin*r) + 1 + j*evtPerBin;
483  if (irc[0]<0) {
484  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[0] = "
485  << irc[0] << " now set to 0";
486  irc[0] = 0;
487  } else if (irc[0] > totEvents) {
488  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[0] = "
489  << irc[0] << " now set to "<< totEvents;
490  irc[0] = totEvents;
491  }
492  }
493  }
494  }
495  if (irc[1]<1) {
496  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[1] = "
497  << irc[1] << " now set to 1";
498  irc[1] = 1;
499  } else if (irc[1] > totEvents) {
500  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[1] = "
501  << irc[1] << " now set to "<< totEvents;
502  irc[1] = totEvents;
503  }
504 
505 #ifdef DebugLog
506  LogDebug("HFShower") << "HFShowerLibrary:: Select records " << irc[0]
507  << " and " << irc[1] << " with weights " << 1-w
508  << " and " << w;
509 #endif
510  pe.clear();
511  npe = 0;
512  int npold = 0;
513  for (int ir=0; ir < 2; ir++) {
514  if (irc[ir]>0) {
515  getRecord (type, irc[ir]);
516  int nPhoton = (newForm) ? photo->size() : photon.size();
517  npold += nPhoton;
518  for (int j=0; j<nPhoton; j++) {
519  r = G4UniformRand();
520  if ((ir==0 && r > w) || (ir > 0 && r < w)) {
521  storePhoton (j);
522  }
523  }
524  }
525  }
526 
527  if ((npe > npold || (npold == 0 && irc[0] > 0)) && !(npe == 0 && npold == 0))
528  edm::LogWarning("HFShower") << "HFShowerLibrary: Interpolation Warning =="
529  << " records " << irc[0] << " and " << irc[1]
530  << " gives a buffer of " << npold
531  << " photons and fills " << npe << " *****";
532 #ifdef DebugLog
533  else
534  LogDebug("HFShower") << "HFShowerLibrary: Interpolation == records "
535  << irc[0] << " and " << irc[1] << " gives a "
536  << "buffer of " << npold << " photons and fills "
537  << npe << " PE";
538  for (int j=0; j<npe; j++)
539  LogDebug("HFShower") << "Photon " << j << " " << pe[j];
540 #endif
541 }
542 
543 void HFShowerLibrary::extrapolate(int type, double pin) {
544 
545  int nrec = int(pin/pmom[nMomBin-1]);
546  double w = (pin - pmom[nMomBin-1]*nrec)/pmom[nMomBin-1];
547  nrec++;
548 #ifdef DebugLog
549  LogDebug("HFShower") << "HFShowerLibrary:: Extrapolate for Energy " << pin
550  << " GeV with " << nMomBin << " momentum bins and "
551  << evtPerBin << " entries/bin -- total " << totEvents
552  << " using " << nrec << " records";
553 #endif
554  std::vector<int> irc(nrec);
555 
556  for (int ir=0; ir<nrec; ir++) {
557  double r = G4UniformRand();
558  irc[ir] = int(evtPerBin*0.5*r) +(nMomBin-1)*evtPerBin + 1;
559  if (irc[ir]<1) {
560  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[" << ir
561  << "] = " << irc[ir] << " now set to 1";
562  irc[ir] = 1;
563  } else if (irc[ir] > totEvents) {
564  edm::LogWarning("HFShower") << "HFShowerLibrary:: Illegal irc[" << ir
565  << "] = " << irc[ir] << " now set to "
566  << totEvents;
567  irc[ir] = totEvents;
568 #ifdef DebugLog
569  } else {
570  LogDebug("HFShower") << "HFShowerLibrary::Extrapolation use irc["
571  << ir << "] = " << irc[ir];
572 #endif
573  }
574  }
575 
576  pe.clear();
577  npe = 0;
578  int npold = 0;
579  for (int ir=0; ir<nrec; ir++) {
580  if (irc[ir]>0) {
581  getRecord (type, irc[ir]);
582  int nPhoton = (newForm) ? photo->size() : photon.size();
583  npold += nPhoton;
584  for (int j=0; j<nPhoton; j++) {
585  double r = G4UniformRand();
586  if (ir != nrec-1 || r < w) {
587  storePhoton (j);
588  }
589  }
590 #ifdef DebugLog
591  LogDebug("HFShower") << "HFShowerLibrary: Record [" << ir << "] = "
592  << irc[ir] << " npold = " << npold;
593 #endif
594  }
595  }
596 #ifdef DebugLog
597  LogDebug("HFShower") << "HFShowerLibrary:: uses " << npold << " photons";
598 #endif
599 
600  if (npe > npold || npold == 0)
601  edm::LogWarning("HFShower") << "HFShowerLibrary: Extrapolation Warning == "
602  << nrec << " records " << irc[0] << ", "
603  << irc[1] << ", ... gives a buffer of " <<npold
604  << " photons and fills " << npe
605  << " *****";
606 #ifdef DebugLog
607  else
608  LogDebug("HFShower") << "HFShowerLibrary: Extrapolation == " << nrec
609  << " records " << irc[0] << ", " << irc[1]
610  << ", ... gives a buffer of " << npold
611  << " photons and fills " << npe << " PE";
612  for (int j=0; j<npe; j++)
613  LogDebug("HFShower") << "Photon " << j << " " << pe[j];
614 #endif
615 }
616 
618 
619  if (newForm) pe.push_back(photo->at(j));
620  else pe.push_back(photon[j]);
621 #ifdef DebugLog
622  LogDebug("HFShower") << "HFShowerLibrary: storePhoton " << j << " npe "
623  << npe << " " << pe[npe];
624 #endif
625  npe++;
626 }
627 
628 std::vector<double> HFShowerLibrary::getDDDArray(const std::string & str,
629  const DDsvalues_type & sv,
630  int & nmin) {
631 
632 #ifdef DebugLog
633  LogDebug("HFShower") << "HFShowerLibrary:getDDDArray called for " << str
634  << " with nMin " << nmin;
635 #endif
636  DDValue value(str);
637  if (DDfetch(&sv,value)) {
638 #ifdef DebugLog
639  LogDebug("HFShower") << value;
640 #endif
641  const std::vector<double> & fvec = value.doubles();
642  int nval = fvec.size();
643  if (nmin > 0) {
644  if (nval < nmin) {
645  edm::LogError("HFShower") << "HFShowerLibrary : # of " << str
646  << " bins " << nval << " < " << nmin
647  << " ==> illegal";
648  throw cms::Exception("Unknown", "HFShowerLibrary")
649  << "nval < nmin for array " << str << "\n";
650  }
651  } else {
652  if (nval < 2) {
653  edm::LogError("HFShower") << "HFShowerLibrary : # of " << str
654  << " bins " << nval << " < 2 ==> illegal"
655  << " (nmin=" << nmin << ")";
656  throw cms::Exception("Unknown", "HFShowerLibrary")
657  << "nval < 2 for array " << str << "\n";
658  }
659  }
660  nmin = nval;
661 
662  return fvec;
663  } else {
664  edm::LogError("HFShower") << "HFShowerLibrary : cannot get array " << str;
665  throw cms::Exception("Unknown", "HFShowerLibrary")
666  << "cannot get array " << str << "\n";
667  }
668 }
#define LogDebug(id)
double tShift(const G4ThreeVector &point, int depth, int fromEndAbs=0)
Definition: HFFibre.cc:115
type
Definition: HCALResponse.h:21
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
Definition: DDValue.cc:140
static const TGPicture * info(bool iBackgroundIsBlack)
const double GeV
Definition: MathUtil.h:16
void getRecord(int, int)
void initRun(const HcalDDDSimConstants *)
Definition: HFFibre.cc:82
const double w
Definition: UKUtility.cc:23
JetCorrectorParameters::Record record
Definition: classes.h:7
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:127
Definition: weight.py:1
#define nullptr
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:83
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
static bool isStableHadron(int pdgCode)
const double MeV
std::vector< HFShowerPhoton > HFShowerPhotonCollection
const std::vector< double > & getRTableHF() const
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
std::maps an index to a DDValue. The index corresponds to the index assigned to the name of the std::...
Definition: DDsvalues.h:20
void loadEventInfo(TBranch *)
const std::vector< double > & getPhiTableHF() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
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
TBranch * emBranch
const std::vector< double > & getGparHF() const
std::vector< double > pmom
G4ThreeVector position
void extrapolate(int, double)
HFShowerLibrary(const std::string &name, const DDCompactView &cpv, edm::ParameterSet const &p)
void storePhoton(int j)
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)
HFShowerPhotonCollection photon
bool rInside(double r)
#define str(s)
void interpolate(int, double)
void initRun(G4ParticleTable *, const HcalDDDSimConstants *)
TBranch * hadBranch