CMS 3D CMS Logo

CSCGasCollisions.cc
Go to the documentation of this file.
1 // This implements the CSC gas ionization simulation.
2 // It requires the GEANT3-generated tables in the input file
3 // (default) MuonData/EndcapData/collisions.dat
4 
5 // Outstanding problems to fix 15-Oct-2003
6 // - ework must be re-tuned once eion or ework are removed before delta e range
7 // estimation.
8 // - logic of long-range delta e's must be revisited.
9 // - The code here needs to have diagnostic output removed, or reduced.
10 // PARTICULARLY filling of std::vectors!
11 // - 'gap' is a misnomer, when simhit entry and exit don't coincide with gap
12 // edges
13 // so the CSCCrossGap might better be named as a simhit-related thing.
14 // 22-Jan-2004 Corrected position of trap for 'infinite' loop while
15 // generating steps. Output files (debugV-flagged only) require SC flag.
16 // Increase deCut from 10 keV to 100 keV to accomodate protons!
17 // Mar-2015: cout to LogVerbatim. Change superseded debugV-flagged code to flag
18 // from config.
19 
25 
26 #include "CLHEP/Random/RandExponential.h"
27 #include "CLHEP/Random/RandFlat.h"
28 
29 // Only needed if file writing code is activated again
30 // #include <iostream>
31 #include <cassert>
32 #include <fstream>
33 
34 // stdlib math functions
35 // for 'abs':
36 #include <cstdlib>
37 // for usual math functions:
38 #include <cmath>
39 
40 // stdlib container trickery
41 // for 'for_each':
42 #include <algorithm>
43 // for 'greater' and 'less':
44 #include <functional>
45 // for 'accumulate':
46 #include <iterator>
47 #include <numeric>
48 
49 using namespace std;
50 
51 /* Gas mixture is Ar/CO2/CF4 = 40/50/10
52 We'll use the ionization energies
53 Ar 15.8 eV
54 CO2 13.7 eV
55 CF4 17.8
56 to arrive at a weighted average of eion = 14.95 */
57 
59  : me("CSCGasCollisions"),
60  gasDensity(2.1416e-03),
61  deCut(1.e05),
62  eion(14.95),
63  ework(34.0),
64  clusterExtent(0.001),
65  theGammaBins(N_GAMMA, 0.),
66  theEnergyBins(N_ENERGY, 0.),
67  theCollisionTable(N_ENTRIES, 0.),
68  theCrossGap(nullptr),
69  theParticleDataTable(nullptr),
70  saveGasCollisions_(false),
71  dumpGasCollisions_(false) {
72  dumpGasCollisions_ = pset.getUntrackedParameter<bool>("dumpGasCollisions");
73 
74  edm::LogInfo(me) << "Constructing a " << me << ":";
75  edm::LogInfo(me) << "gas density = " << gasDensity << " g/cm3";
76  edm::LogInfo(me) << "max eloss per collision allowed = " << deCut / 1000.
77  << " keV (for higher elosses, hits should have been simulated.)";
78  edm::LogInfo(me) << "ionization threshold = " << eion << " eV";
79  edm::LogInfo(me) << "effective work function = " << ework << " eV";
80  edm::LogInfo(me) << "cluster extent = " << clusterExtent * 1.e04 << " micrometres";
81  edm::LogInfo(me) << "dump gas collision and simhit information? " << dumpGasCollisions();
82  edm::LogInfo(me) << "save gas collision information? -NOT YET IMPLEMENTED- " << saveGasCollisions();
83 
85 }
86 
88  edm::LogInfo(me) << "Destructing a " << me;
89  delete theCrossGap;
90 }
91 
93  // I'd prefer to allow comments in the data file which means
94  // complications of reading line-by-line and then item-by-item.
95 
96  // Is float OK? Or do I need double?
97 
98  // Do I need any sort of error trapping?
99 
100  // We use the default CMSSW data file path
101  // We use the default file name 'collisions.dat'
102 
103  // This can be reset in .orcarc by SimpleConfigurable
104  // Muon:Endcap:CollisionsFile
105 
106  // TODO make configurable
107  string colliFile = "SimMuon/CSCDigitizer/data/collisions.dat";
108  edm::FileInPath f1{colliFile};
109  edm::LogInfo(me) << ": reading " << f1.fullPath();
110 
111  ifstream fin{f1.fullPath()};
112 
113  if (fin.fail()) {
114  string errorMessage = "Cannot open input file " + f1.fullPath();
115  edm::LogError("CSCGasCollisions") << errorMessage;
116  throw cms::Exception(errorMessage);
117  }
118 
119  fin.clear(); // Clear eof read status
120  fin.seekg(0, ios::beg); // Position at start of file
121 
122  // @@ We had better have the right sizes everywhere or all
123  // hell will break loose. There's no trapping.
124  LogTrace(me) << "Reading gamma bins";
125  int j = 0;
126  for (int i = 0; i < N_GAMMA; ++i) {
127  fin >> theGammaBins[i];
128  LogTrace(me) << ++j << " " << theGammaBins[i];
129  }
130 
131  j = 0;
132  LogTrace(me) << "Reading energy bins \n";
133  for (int i = 0; i < N_ENERGY; ++i) {
134  fin >> theEnergyBins[i];
135  LogTrace(me) << ++j << " " << theEnergyBins[i];
136  }
137 
138  j = 0;
139  LogTrace(me) << "Reading collisions table \n";
140 
141  for (int i = 0; i < N_ENTRIES; ++i) {
142  fin >> theCollisionTable[i];
143  LogTrace(me) << ++j << " " << theCollisionTable[i];
144  }
145 
146  fin.close();
147 }
148 
150 
152  std::vector<LocalPoint> &positions,
153  std::vector<int> &electrons,
154  CLHEP::HepRandomEngine *engine) {
155  const float epsilonL = 0.01; // Shortness of simhit 'length'
156  // const float max_gap_z = 1.5; // Gas gaps are 0.5
157  // or 1.0 cm
158 
159  // Note that what I call the 'gap' may in fact be the 'length' of a PSimHit
160  // which does not start and end on the gap edges. This confuses the
161  // nomenclature at least.
162 
163  double mom = simhit.pabs(); // in GeV/c - see MuonSensitiveDetector.cc
164  // int iam = simhit.particleType(); // PDG type
165  delete theCrossGap; // before building new one
166  assert(theParticleDataTable != nullptr);
167  ParticleData const *particle = theParticleDataTable->particle(simhit.particleType());
168  double mass = 0.105658; // assume a muon
169  if (particle == nullptr) {
170  edm::LogError("CSCGasCollisions") << "Cannot find particle of type " << simhit.particleType() << " in the PDT";
171  } else {
172  mass = particle->mass();
173  }
174 
175  theCrossGap = new CSCCrossGap(mass, mom, simhit.exitPoint() - simhit.entryPoint());
176  float gapSize = theCrossGap->length();
177 
178  // Test the simhit 'length' (beware of angular effects)
179  // if ( gapSize <= epsilonL || gapSize > max_gap_z ) {
180  if (gapSize <= epsilonL) {
181  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! simhit entry and exit are too close - "
182  "skipping simhit:"
183  << "\n entry = " << simhit.entryPoint() << ": exit = " << simhit.exitPoint()
184  << "\n particle type = " << simhit.particleType()
185  << " : momentum = " << simhit.pabs()
186  << " GeV/c : energy loss = " << simhit.energyLoss() * 1.E06 << " keV"
187  << ", gapSize = " << gapSize << " cm (< epsilonL = " << epsilonL << " cm)";
188  return; //@@ Just skip this PSimHit
189  }
190 
191  // Interpolate the table for current gamma value
192  // Extract collisions binned by energy loss values, for this gamma
193  std::vector<float> collisions(N_ENERGY);
194  double loggam = theCrossGap->logGamma();
195  fillCollisionsForThisGamma(static_cast<float>(loggam), collisions);
196 
197  double anmin = exp(collisions[N_ENERGY - 1]);
198  double anmax = exp(collisions[0]);
199  double amu = anmax - anmin;
200 
201  LogTrace(me) << "collisions extremes = " << collisions[N_ENERGY - 1] << ", " << collisions[0] << "\n"
202  << "anmin = " << anmin << ", anmax = " << anmax << "\n"
203  << "amu = " << amu << "\n";
204 
205  float dedx = 0.; // total energy loss
206  double sum_steps = 0.; // total distance across gap (along simhit direction)
207  int n_steps = 0; // no. of steps/primary collisions
208  int n_try = 0; // no. of tries to generate steps
209  double step = -1.; // Sentinel for start
210 
211  LocalPoint layerLocalPoint(simhit.entryPoint());
212 
213  // step/primary collision loop
214  while (sum_steps < gapSize) {
215  ++n_try;
216  if (n_try > MAX_STEPS) {
217  int maxst = MAX_STEPS;
218  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! n_try = " << n_try
219  << " is > MAX_STEPS = " << maxst << " - skipping simhit:"
220  << "\n entry = " << simhit.entryPoint() << ": exit = " << simhit.exitPoint()
221  << "\n particle type = " << simhit.particleType()
222  << " : momentum = " << simhit.pabs()
223  << " GeV/c : energy loss = " << simhit.energyLoss() * 1.E06 << " keV"
224  << "\n gapSize = " << gapSize << " cm, last step = " << step
225  << " cm, sum_steps = " << sum_steps << " cm, n_steps = " << n_steps;
226  break;
227  }
228  step = generateStep(amu, engine);
229  if (sum_steps + step > gapSize)
230  break; // this step goes too far
231 
232  float eloss = generateEnergyLoss(amu, anmin, anmax, collisions, engine);
233 
234  // Is the eloss too large? (then GEANT should have produced hits!)
235  if (eloss > deCut) {
236  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! eloss = " << eloss << " eV is too large (> "
237  << deCut << " eV) - trying another collision";
238  continue; // to generate another collision/step
239  }
240 
241  dedx += eloss; // the energy lost from the ionizing particle
242  sum_steps += step; // the position of the ionizing particle
243  ++n_steps; // the number of primary collisions
244 
245  if (n_steps > MAX_STEPS) { // Extra-careful trap for bizarreness
246  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! " << n_steps
247  << " is too many steps - skipping simhit:"
248  << "\n entry = " << simhit.entryPoint() << ": exit = " << simhit.exitPoint()
249  << "\n particle type = " << simhit.particleType()
250  << " : momentum = " << simhit.pabs()
251  << " GeV/c : energy loss = " << simhit.energyLoss() * 1.E06 << " keV"
252  << "\ngapSize=" << gapSize << " cm, last step=" << step
253  << " cm, sum_steps=" << sum_steps << " cm";
254  break;
255  }
256  LogTrace(me) << "sum_steps = " << sum_steps << " cm , dedx = " << dedx << " eV";
257 
258  // Generate ionization.
259  // eion is the minimum energy at which ionization can occur in the gas
260  if (eloss > eion) {
261  layerLocalPoint += step * theCrossGap->unitVector(); // local point where the collision occurs
262  ionize(eloss, layerLocalPoint);
263  } else {
264  LogTrace(me) << "Energy available = " << eloss << " eV is too low for ionization (< eion = " << eion << " eV)";
265  }
266 
267  } // step/collision loop
268 
269  if (dumpGasCollisions())
270  writeSummary(n_try, n_steps, sum_steps, dedx, simhit);
271 
272  // Return values in two container arguments
273  positions = theCrossGap->ionClusters();
274  electrons = theCrossGap->electrons();
275  return;
276 }
277 
278 double CSCGasCollisions::generateStep(double avCollisions, CLHEP::HepRandomEngine *engine) const {
279  // Generate a m.f.p. (1/avCollisions = cm/collision)
280  double step = (CLHEP::RandExponential::shoot(engine)) / avCollisions;
281 
282  // Without using CLHEP: approx random exponential by...
283  // double da = double(rand())/double(RAND_MAX);
284  // double step = -log(1.-da)/avCollisions;
285 
286  LogTrace(me) << " step = " << step << " cm";
287  // Next line only used to fill a container of 'step's for later diagnostic
288  // dumps
289  if (dumpGasCollisions())
290  theCrossGap->addStep(step);
291  return step;
292 }
293 
294 float CSCGasCollisions::generateEnergyLoss(double avCollisions,
295  double anmin,
296  double anmax,
297  const std::vector<float> &collisions,
298  CLHEP::HepRandomEngine *engine) const {
299  // Generate a no. of collisions between collisions[0] and [N_ENERGY-1]
300  float lnColl = log(CLHEP::RandFlat::shoot(engine, anmin, anmax));
301 
302  // Without using CLHEP: approx random between anmin and anmax
303  // double ra = double(rand())/double(RAND_MAX)*avCollisions;
304  // cout << "ra = " << ra << std::endl;
305  // float lnColl = static_cast<float>( log( ra ) );
306 
307  // Find energy loss for that number
308  float lnE = lnEnergyLoss(lnColl, collisions);
309  float eloss = exp(lnE);
310  // Compensate if gamma was actually below 1.1
311  if (theCrossGap->gamma() < 1.1)
312  eloss = eloss * 0.173554 / theCrossGap->beta2();
313  LogTrace(me) << "eloss = " << eloss << " eV";
314  // Next line only used to fill container of eloss's for later diagnostic dumps
315  if (dumpGasCollisions())
316  theCrossGap->addEloss(eloss);
317  return eloss;
318 }
319 
320 void CSCGasCollisions::ionize(double energyAvailable, LocalPoint startHere) const {
321  while (energyAvailable > eion) {
322  LogTrace(me) << " NEW CLUSTER " << theCrossGap->noOfClusters() + 1 << " AT " << startHere;
323  LocalPoint newCluster(startHere);
324  theCrossGap->addCluster(newCluster);
325 
326  //@@ I consider NOT subtracting eion before calculating range to be a bug.
327  //@@ But this changes tuning of the algorithm so leave it until after the
328  // big rush to 7_5_0
329  //@@ energyAvailable -= eion;
330 
331  // Sauli CERN 77-09: delta e range with E in MeV (Sauli references Kobetich
332  // & Katz 1968 but that has nothing to do with this expression! He seems to
333  // have made a mistake.) I found the Sauli-quoted relationship in R.
334  // Glocker, Z. Naturforsch. Ba, 129 (1948): delta e range R = aE^n with
335  // a=710, n=1.72 for E in MeV and R in mg/cm^2 applicable over the range E =
336  // 0.001 to 0.3 MeV.
337 
338  // Take HALF that range. //@@ Why? Why not...
339  double range = 0.5 * (0.71 / gasDensity) * pow(energyAvailable * 1.E-6, 1.72);
340  LogTrace(me) << " range = " << range << " cm";
341  if (range < clusterExtent) {
342  // short-range delta e
343  // How many electrons can we make? Now use *average* energy for ionization
344  // (not *minimum*)
345  int nelec = static_cast<int>(energyAvailable / ework);
346  LogTrace(me) << "short-range delta energy in = " << energyAvailable << " eV";
347  // energyAvailable -= nelec*(energyAvailable/ework);
348  energyAvailable -= nelec * ework;
349  // If still above eion (minimum, not average) add one more e
350  if (energyAvailable > eion) {
351  ++nelec;
352  energyAvailable -= eion;
353  }
354  LogTrace(me) << "short-range delta energy out = " << energyAvailable << " eV, nelec = " << nelec;
355  theCrossGap->addElectrons(nelec);
356  break;
357 
358  } else {
359  // long-range delta e
360  LogTrace(me) << "long-range delta \n"
361  << "no. of electrons in cluster now = " << theCrossGap->noOfElectrons();
362  theCrossGap->addElectrons(1); // Position is at startHere still
363 
364  bool new_range = false;
365  while (!new_range && (energyAvailable > ework)) {
366  energyAvailable -= ework;
367  while (energyAvailable > eion) {
368  double range2 = 0.5 * 0.71 / gasDensity * pow(1.E-6 * energyAvailable, 1.72);
369  double drange = range - range2;
370  LogTrace(me) << " energy left = " << energyAvailable << " eV, range2 = " << range2
371  << " cm, drange = " << drange << " cm";
372  if (drange < clusterExtent) {
373  theCrossGap->addElectronToBack(); // increment last element
374  } else {
375  startHere += drange * theCrossGap->unitVector(); // update delta e start position
376  range = range2; // update range
377  new_range = true; // Test range again
378  LogTrace(me) << "reset range to range2 = " << range << " from startHere = " << startHere << " and iterate";
379  }
380  break; // out of inner while energyAvailable>eion
381 
382  } // inner while energyAvailable>eion
383 
384  } // while !new_range && energyAvailable>ework
385 
386  // energyAvailable now less than ework, but still may be over eion...add
387  // an e
388  if (energyAvailable > eion) {
389  energyAvailable -= ework; // yes, it may go negative
390  theCrossGap->addElectronToBack(); // add one more e
391  }
392 
393  } // if range
394 
395  } // outer while energyAvailable>eion
396 }
397 
398 void CSCGasCollisions::writeSummary(int n_try, int n_steps, double sum_steps, float dedx, const PSimHit &simhit) const {
399  // Switched from std::cout to LogVerbatim, Mar 2015
400 
401  std::vector<LocalPoint> ion_clusters = theCrossGap->ionClusters();
402  std::vector<int> electrons = theCrossGap->electrons();
403  std::vector<float> elosses = theCrossGap->eLossPerStep();
404  std::vector<double> steps = theCrossGap->stepLengths();
405 
406  edm::LogVerbatim("CSCDigitizer") << "------------------"
407  << "\nAFTER CROSSING GAP";
408  /*
409  edm::LogVerbatim("CSCDigitizer") << "no. of steps tried = " <<
410  n_try
411  << "\nno. of steps from theCrossGap = " <<
412  theCrossGap->noOfSteps()
413  << "\nsize of steps vector = " <<
414  steps.size();
415 
416  edm::LogVerbatim("CSCDigitizer") << "no. of collisions (steps) = " <<
417  n_steps
418  << "\nsize of elosses vector = " <<
419  elosses.size()
420  << "\nsize of ion clusters vector = " <<
421  ion_clusters.size()
422  << "\nsize of electrons vector = " <<
423  electrons.size();
424  */
425 
426  size_t nsteps = n_steps; // force ridiculous type conversion
427  size_t mstep = steps.size() - 1; // final step gets filled but is outside gas
428  // gap - unless we reach MAX_STEPS
429  if ((nsteps != MAX_STEPS) && (nsteps != mstep)) {
430  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! no. of steps = " << nsteps
431  << " .ne. steps.size()-1 = " << mstep;
432  }
433  size_t meloss = elosses.size();
434 
435  if (nsteps != meloss) {
436  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! no. of steps = " << nsteps
437  << " .ne. no. of elosses = " << meloss;
438  } else {
439  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] # / length of step / energy loss per collision:";
440  for (size_t i = 0; i != nsteps; ++i) {
441  edm::LogVerbatim("CSCDigitizer") << i + 1 << " / S: " << steps[i] << " / E: " << elosses[i];
442  }
443  }
444 
445  size_t mclus = ion_clusters.size();
446  size_t melec = electrons.size();
447 
448  if (mclus != melec) {
449  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! size of cluster vector = " << mclus
450  << " .ne. size of electrons vector = " << melec;
451  } else {
452  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] # / postion of "
453  "cluster / electrons per cluster: ";
454  for (size_t i = 0; i != mclus; ++i) {
455  edm::LogVerbatim("CSCDigitizer") << i + 1 << " / I: " << ion_clusters[i] << " / E: " << electrons[i];
456  }
457  }
458 
459  int n_ic = count_if(elosses.begin(), elosses.end(), [&](auto c) { return c > this->eion; });
460 
461  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] total no. of collision steps = " << n_steps;
462  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] total sum of steps = " << sum_steps << " cm";
463  if (nsteps > 0)
464  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] average step length = " << sum_steps / float(nsteps)
465  << " cm";
466  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] total energy loss across gap = " << dedx
467  << " eV = " << dedx / 1000. << " keV";
468  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] no. of primary ionizing collisions across gap = "
469  "no. with eloss > eion = "
470  << eion << " eV = " << n_ic;
471  if (nsteps > 0)
472  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] average energy loss/collision = " << dedx / float(nsteps)
473  << " eV";
474 
475  std::vector<int>::const_iterator bigger = find(electrons.begin(), electrons.end(), 0);
476  if (bigger != electrons.end()) {
477  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] TROUBLE! There is a cluster with 0 electrons.";
478  }
479 
480  int n_e = accumulate(electrons.begin(), electrons.end(), 0);
481 
482  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] SUMMARY: simhit"
483  << "\n entry = " << simhit.entryPoint() << ": exit = " << simhit.exitPoint()
484  << "\n particle type = " << simhit.particleType()
485  << " : momentum = " << simhit.pabs()
486  << " GeV/c : energy loss = " << simhit.energyLoss() * 1.E06 << " keV";
487 
488  if (nsteps > 0) {
489  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] SUMMARY: ionization"
490  << " : steps= " << nsteps << " : sum(steps)= " << sum_steps
491  << " cm : <step>= " << sum_steps / float(nsteps) << " cm"
492  << " : ionizing= " << n_ic << " : ionclus= " << mclus << " : total e= " << n_e
493  << " : <dedx/step>= " << dedx / float(nsteps)
494  << " eV : <e/ionclus>= " << float(n_e) / float(mclus)
495  << " : dedx= " << dedx / 1000. << " keV";
496  } else {
497  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisons] ERROR? no collision steps!";
498  }
499 
500  // Turn off output file - used for initial development
501  // if ( saveGasCollisions() ) {
502  // ofstream of0("osteplen.dat",ios::app);
503  // std::copy( steps.begin(), steps.end(),
504  // std::ostream_iterator<float>(of0,"\n"));
505 
506  // ofstream of1("olperc.dat",ios::app);
507  // std::copy( elosses.begin(), elosses.end(),
508  // std::ostream_iterator<float>(of1,"\n"));
509 
510  // ofstream of2("oclpos.dat",ios::app);
511  // std::copy( ion_clusters.begin(), ion_clusters.end(),
512  // std::ostream_iterator<LocalPoint>(of2,"\n"));
513 
514  // ofstream of3("oepercl.dat",ios::app);
515  // std::copy( electrons.begin(), electrons.end(),
516  // std::ostream_iterator<int>(of3,"\n"));
517  // }
518 }
519 
520 float CSCGasCollisions::lnEnergyLoss(float lnCollisions, const std::vector<float> &collisions) const {
521  float lnE = -1.;
522 
523  // Find collision[] bin in which lnCollisions falls
524  std::vector<float>::const_iterator it = find(collisions.begin(), collisions.end(), lnCollisions);
525 
526  if (it != collisions.end()) {
527  // found the value
528  std::vector<float>::difference_type ihi = it - collisions.begin();
529  LogTrace(me) << ": using one energy bin " << ihi << " = " << theEnergyBins[ihi]
530  << " for lnCollisions = " << lnCollisions;
531  lnE = theEnergyBins[ihi];
532  } else {
533  // interpolate the value
534  std::vector<float>::const_iterator loside =
535  find_if(collisions.begin(), collisions.end(), [&lnCollisions](auto c) { return c < lnCollisions; });
536  std::vector<float>::difference_type ilo = loside - collisions.begin();
537  if (ilo > 0) {
538  LogTrace(me) << ": using energy bin " << ilo - 1 << " and " << ilo;
539  lnE = theEnergyBins[ilo - 1] + (lnCollisions - collisions[ilo - 1]) *
540  (theEnergyBins[ilo] - theEnergyBins[ilo - 1]) /
541  (collisions[ilo] - collisions[ilo - 1]);
542  } else {
543  LogTrace(me) << ": using one energy bin 0 = " << theEnergyBins[0] << " for lnCollisions = " << lnCollisions;
544  lnE = theEnergyBins[0]; //@@ WHAT ELSE TO DO?
545  }
546  }
547 
548  return lnE;
549 }
550 
551 void CSCGasCollisions::fillCollisionsForThisGamma(float logGamma, std::vector<float> &collisions) const {
552  std::vector<float>::const_iterator bigger =
553  find_if(theGammaBins.begin(), theGammaBins.end(), [&logGamma](auto c) { return c > logGamma; });
554 
555  if (bigger == theGammaBins.end()) {
556  // use highest bin
557  LogTrace(me) << ": using highest gamma bin"
558  << " for logGamma = " << logGamma;
559  for (int i = 0; i < N_ENERGY; ++i)
560  collisions[i] = theCollisionTable[i * N_GAMMA];
561  } else {
562  // use bigger and its lower neighbour
563  std::vector<float>::difference_type ihi = bigger - theGammaBins.begin();
564  if (ihi > 0) {
565  double dlg2 = *bigger--; // and decrement after deref
566  // LogTrace(me) << ": using gamma bins "
567  // << ihi-1 << " and " << ihi;
568  double dlg1 = *bigger; // now the preceding element
569  double dlg = (logGamma - dlg1) / (dlg2 - dlg1);
570  double omdlg = 1. - dlg;
571  for (int i = 0; i < N_ENERGY; ++i)
572  collisions[i] = theCollisionTable[i * N_GAMMA + ihi - 1] * omdlg + theCollisionTable[i * N_GAMMA + ihi] * dlg;
573  } else {
574  // bigger has no lower neighbour
575  LogTrace(me) << ": using lowest gamma bin"
576  << " for logGamma = " << logGamma;
577 
578  for (int i = 0; i < N_ENERGY; ++i)
579  collisions[i] = theCollisionTable[i * N_GAMMA];
580  }
581  }
582 }
T getUntrackedParameter(std::string const &, T const &) const
std::vector< int > electrons() const
Definition: CSCCrossGap.h:34
HepPDT::ParticleDataTable ParticleDataTable
bool saveGasCollisions(void) const
void addEloss(float eloss)
Definition: CSCCrossGap.h:46
int noOfClusters() const
Definition: CSCCrossGap.h:33
static const int N_ENERGY
#define nullptr
std::vector< LocalPoint > ionClusters() const
Definition: CSCCrossGap.h:32
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
void addElectronToBack()
Definition: CSCCrossGap.h:43
CSCGasCollisions(const edm::ParameterSet &pset)
std::vector< float > eLossPerStep() const
Definition: CSCCrossGap.h:38
std::vector< double > stepLengths() const
Definition: CSCCrossGap.h:36
std::vector< float > theGammaBins
const std::string me
float lnEnergyLoss(float, const std::vector< float > &) const
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
bool dumpGasCollisions(void) const
const ParticleDataTable * theParticleDataTable
double gamma() const
Definition: CSCCrossGap.h:51
LocalVector unitVector() const
Definition: CSCCrossGap.h:53
int noOfElectrons() const
Definition: CSCCrossGap.h:35
static const int N_GAMMA
static const int MAX_STEPS
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:67
HepPDT::ParticleData ParticleData
void setParticleDataTable(const ParticleDataTable *pdt)
void ionize(double energyTransferred, LocalPoint startHere) const
#define LogTrace(id)
static const int N_ENTRIES
void writeSummary(int n_try, int n_steps, double sum_steps, float dedx, const PSimHit &simhit) const
double logGamma(double mass, float momentum)
Definition: CSCCrossGap.cc:14
std::vector< float > theCollisionTable
double generateStep(double avCollisions, CLHEP::HepRandomEngine *) const
double beta2() const
Definition: CSCCrossGap.h:50
std::vector< float > theEnergyBins
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:79
int particleType() const
Definition: PSimHit.h:89
float generateEnergyLoss(double avCollisions, double anmin, double anmax, const std::vector< float > &collisions, CLHEP::HepRandomEngine *) const
float length() const
Definition: CSCCrossGap.h:54
void addStep(double step)
Definition: CSCCrossGap.h:45
void fillCollisionsForThisGamma(float, std::vector< float > &) const
CSCCrossGap * theCrossGap
step
Definition: StallMonitor.cc:94
void simulate(const PSimHit &, std::vector< LocalPoint > &clusters, std::vector< int > &electrons, CLHEP::HepRandomEngine *)
void addCluster(LocalPoint here)
Definition: CSCCrossGap.h:41
virtual ~CSCGasCollisions()
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void addElectrons(int nelec=1)
Definition: CSCCrossGap.h:42