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