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  string path( getenv( "CMSSW_SEARCH_PATH" ) );
100  // TODO make configurable
101  string colliFile = "SimMuon/CSCDigitizer/data/collisions.dat";
102  FileInPath f1( path, colliFile );
103  if (f1() == 0 ) {
104  string errorMessage = "Input file " + colliFile + "not found";
105  edm::LogError("CSCGasCollisions") << errorMessage << " in path " << path
106  << "\nSet Muon:Endcap:CollisionsFile in .orcarc to the "
107  " location of the file relative to ORCA_DATA_PATH." ;
108  throw cms::Exception( " Endcap Muon gas collisions data file not found.");
109  }
110  else {
111  edm::LogInfo(me) << ": reading " << f1.name();
112  }
113 
114  ifstream & fin = *f1();
115 
116  if (fin.fail()) {
117  string errorMessage = "Cannot open input file " + path + colliFile;
118  edm::LogError("CSCGasCollisions") << errorMessage;
119  throw cms::Exception(errorMessage);
120  }
121 
122  fin.clear( ); // Clear eof read status
123  fin.seekg( 0, ios::beg ); // Position at start of file
124 
125  // @@ We had better have the right sizes everywhere or all
126  // hell will break loose. There's no trapping.
127  LogTrace(me) << "Reading gamma bins";
128  int j = 0;
129  for(int i = 0; i<N_GAMMA; ++i) {
130  fin >> theGammaBins[i];
131  LogTrace(me) << ++j << " " << theGammaBins[i];
132  }
133 
134  j = 0;
135  LogTrace(me) << "Reading energy bins \n";
136  for(int i = 0; i<N_ENERGY; ++i) {
137  fin >> theEnergyBins[i];
138  LogTrace(me) << ++j << " " << theEnergyBins[i];
139  }
140 
141  j = 0;
142  LogTrace(me) << "Reading collisions table \n";
143 
144  for(int i = 0; i<N_ENTRIES; ++i) {
145  fin >> theCollisionTable[i];
146  LogTrace(me) << ++j << " " << theCollisionTable[i];
147  }
148 
149  fin.close();
150 }
151 
152 
154 {
155  theParticleDataTable = pdt;
156 }
157 
158 
159 void CSCGasCollisions::simulate( const PSimHit& simhit,
160  std::vector<LocalPoint>& positions, std::vector<int>& electrons, CLHEP::HepRandomEngine* engine ) {
161 
162  const float epsilonL = 0.01; // Shortness of simhit 'length'
163  // const float max_gap_z = 1.5; // Gas gaps are 0.5 or 1.0 cm
164 
165  // Note that what I call the 'gap' may in fact be the 'length' of a PSimHit which
166  // does not start and end on the gap edges. This confuses the nomenclature at least.
167 
168  double mom = simhit.pabs(); // in GeV/c - see MuonSensitiveDetector.cc
169  // int iam = simhit.particleType(); // PDG type
170  delete theCrossGap; // before building new one
171  assert(theParticleDataTable != 0);
172  ParticleData const * particle = theParticleDataTable->particle( simhit.particleType() );
173  double mass = 0.105658; // assume a muon
174  if(particle == 0)
175  {
176  edm::LogError("CSCGasCollisions") << "Cannot find particle of type " << simhit.particleType()
177  << " in the PDT";
178  }
179  else
180  {
181  mass = particle->mass();
182  }
183 
184  theCrossGap = new CSCCrossGap( mass, mom, simhit.exitPoint() - simhit.entryPoint() );
185  float gapSize = theCrossGap->length();
186 
187  // Test the simhit 'length' (beware of angular effects)
188  // if ( gapSize <= epsilonL || gapSize > max_gap_z ) {
189  if ( gapSize <= epsilonL ) {
190  edm::LogVerbatim("CSCDigitizer")
191  << "[CSCGasCollisions] WARNING! simhit entry and exit are too close - skipping simhit:"
192  << "\n entry = " << simhit.entryPoint()
193  << ": exit = " << simhit.exitPoint()
194  << "\n particle type = " << simhit.particleType() << " : momentum = " << simhit.pabs()
195  << " GeV/c : energy loss = " << simhit.energyLoss()*1.E06 << " keV"
196  << ", gapSize = " << gapSize << " cm (< epsilonL = " << epsilonL << " cm)";
197  return; //@@ Just skip this PSimHit
198  }
199 
200  // Interpolate the table for current gamma value
201  // Extract collisions binned by energy loss values, for this gamma
202  std::vector<float> collisions(N_ENERGY);
203  double loggam = theCrossGap->logGamma();
204  fillCollisionsForThisGamma( static_cast<float>(loggam), collisions );
205 
206  double anmin = exp(collisions[N_ENERGY-1]);
207  double anmax = exp(collisions[0]);
208  double amu = anmax - anmin;
209 
210  LogTrace(me) << "collisions extremes = " << collisions[N_ENERGY-1]
211  << ", " << collisions[0] << "\n"
212  << "anmin = " << anmin << ", anmax = " << anmax << "\n"
213  << "amu = " << amu << "\n";
214 
215  float dedx = 0.; // total energy loss
216  double sum_steps = 0.; // total distance across gap (along simhit direction)
217  int n_steps = 0; // no. of steps/primary collisions
218  int n_try = 0; // no. of tries to generate steps
219  double step = -1.; // Sentinel for start
220 
221  LocalPoint layerLocalPoint( simhit.entryPoint() );
222 
223  // step/primary collision loop
224  while ( sum_steps < gapSize) {
225  ++n_try;
226  if ( n_try > MAX_STEPS ) {
227  int maxst = MAX_STEPS;
228  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! n_try = " << n_try
229  << " is > MAX_STEPS = " << maxst << " - skipping simhit:"
230  << "\n entry = " << simhit.entryPoint()
231  << ": exit = " << simhit.exitPoint()
232  << "\n particle type = " << simhit.particleType() << " : momentum = " << simhit.pabs()
233  << " GeV/c : energy loss = " << simhit.energyLoss()*1.E06 << " keV"
234  << "\n gapSize = " << gapSize << " cm, last step = " << step
235  << " cm, sum_steps = " << sum_steps << " cm, n_steps = " << n_steps;
236  break;
237  }
238  step = generateStep( amu, engine );
239  if ( sum_steps + step > gapSize ) break; // this step goes too far
240 
241  float eloss = generateEnergyLoss( amu, anmin, anmax, collisions, engine );
242 
243  // Is the eloss too large? (then GEANT should have produced hits!)
244  if ( eloss > deCut ) {
245  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! eloss = " << eloss
246  << " eV is too large (> " << deCut << " eV) - trying another collision";
247  continue; // to generate another collision/step
248  }
249 
250  dedx += eloss; // the energy lost from the ionizing particle
251  sum_steps += step; // the position of the ionizing particle
252  ++n_steps; // the number of primary collisions
253 
254  if (n_steps > MAX_STEPS ) { // Extra-careful trap for bizarreness
255  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! " << n_steps
256  << " is too many steps - skipping simhit:"
257  << "\n entry = " << simhit.entryPoint()
258  << ": exit = " << simhit.exitPoint()
259  << "\n particle type = " << simhit.particleType() << " : momentum = " << simhit.pabs()
260  << " GeV/c : energy loss = " << simhit.energyLoss()*1.E06 << " keV"
261  << "\ngapSize=" << gapSize << " cm, last step=" << step << " cm, sum_steps=" << sum_steps << " cm";
262  break;
263  }
264  LogTrace(me) << "sum_steps = " << sum_steps << " cm , dedx = " << dedx << " eV";
265 
266  // Generate ionization.
267  // eion is the minimum energy at which ionization can occur in the gas
268  if ( eloss > eion ) {
269  layerLocalPoint += step * theCrossGap->unitVector(); // local point where the collision occurs
270  ionize( eloss, layerLocalPoint );
271  }
272  else {
273  LogTrace(me) << "Energy available = " << eloss << " eV is too low for ionization (< eion = " << eion << " eV)";
274  }
275 
276  } // step/collision loop
277 
278  if ( dumpGasCollisions() ) writeSummary( n_try, n_steps, sum_steps, dedx, simhit );
279 
280  // Return values in two container arguments
281  positions = theCrossGap->ionClusters();
282  electrons = theCrossGap->electrons();
283  return;
284 }
285 
286 double CSCGasCollisions::generateStep( double avCollisions, CLHEP::HepRandomEngine* engine ) const
287 {
288 // Generate a m.f.p. (1/avCollisions = cm/collision)
289  double step = (CLHEP::RandExponential::shoot(engine))/avCollisions;
290 
291 // Without using CLHEP: approx random exponential by...
292 // double da = double(rand())/double(RAND_MAX);
293 // double step = -log(1.-da)/avCollisions;
294 
295  LogTrace(me) << " step = " << step << " cm";
296  // Next line only used to fill a container of 'step's for later diagnostic dumps
297  if ( dumpGasCollisions() ) theCrossGap->addStep( step );
298  return step;
299 }
300 
301 float CSCGasCollisions::generateEnergyLoss( double avCollisions, double anmin, double anmax,
302  const std::vector<float>& collisions, CLHEP::HepRandomEngine* engine ) const
303 {
304 // Generate a no. of collisions between collisions[0] and [N_ENERGY-1]
305  float lnColl = log(CLHEP::RandFlat::shoot(engine, anmin, anmax));
306 
307 // Without using CLHEP: approx random between anmin and anmax
308 // double ra = double(rand())/double(RAND_MAX)*avCollisions;
309 // cout << "ra = " << ra << std::endl;
310 // float lnColl = static_cast<float>( log( ra ) );
311 
312  // Find energy loss for that number
313  float lnE = lnEnergyLoss( lnColl, collisions );
314  float eloss = exp(lnE);
315  // Compensate if gamma was actually below 1.1
316  if ( theCrossGap->gamma() < 1.1 ) eloss = eloss * 0.173554/theCrossGap->beta2();
317  LogTrace(me) << "eloss = " << eloss << " eV";
318  // Next line only used to fill container of eloss's for later diagnostic dumps
319  if ( dumpGasCollisions() ) theCrossGap->addEloss( eloss );
320  return eloss;
321 }
322 
323 void CSCGasCollisions::ionize( double energyAvailable, LocalPoint startHere ) const
324 {
325  while ( energyAvailable > eion ) {
326  LogTrace(me) << " NEW CLUSTER " << theCrossGap->noOfClusters() + 1 <<
327  " AT " << startHere;
328  LocalPoint newCluster( startHere );
329  theCrossGap->addCluster(newCluster);
330 
331  //@@ I consider NOT subtracting eion before calculating range to be a bug.
332  //@@ But this changes tuning of the algorithm so leave it until after the big rush to 7_5_0
333  //@@ energyAvailable -= eion;
334 
335  // Sauli CERN 77-09: delta e range with E in MeV (Sauli references Kobetich & Katz 1968
336  // but that has nothing to do with this expression! He seems to have made a mistake.)
337  // I found the Sauli-quoted relationship in R. Glocker, Z. Naturforsch. Ba, 129 (1948):
338  // delta e range R = aE^n with a=710, n=1.72 for E in MeV and R in mg/cm^2
339  // applicable over the range E = 0.001 to 0.3 MeV.
340 
341  // Take HALF that range. //@@ Why? Why not...
342  double range = 0.5 * (0.71/gasDensity)*pow( energyAvailable*1.E-6, 1.72);
343  LogTrace(me) << " range = " << range << " cm";
344  if ( range < clusterExtent ) {
345 
346  // short-range delta e
347  // How many electrons can we make? Now use *average* energy for ionization (not *minimum*)
348  int nelec = static_cast<int>(energyAvailable/ework);
349  LogTrace(me) << "short-range delta energy in = " << energyAvailable << " eV";
350  //energyAvailable -= nelec*(energyAvailable/ework);
351  energyAvailable -= nelec*ework;
352  // If still above eion (minimum, not average) add one more e
353  if ( energyAvailable > eion ) {
354  ++nelec;
355  energyAvailable -= eion;
356  }
357  LogTrace(me) << "short-range delta energy out = " << energyAvailable << " eV, nelec = " << nelec;
358  theCrossGap->addElectrons( nelec );
359  break;
360 
361  }
362  else {
363  // long-range delta e
364  LogTrace(me) << "long-range delta \n"
365  << "no. of electrons in cluster now = " << theCrossGap->noOfElectrons();
366  theCrossGap->addElectrons( 1 ); // Position is at startHere still
367 
368  bool new_range = false;
369  while ( !new_range && (energyAvailable>ework) ) {
370  energyAvailable -= ework;
371  while ( energyAvailable > eion ) {
372  double range2 = 0.5 * 0.71/gasDensity*pow( 1.E-6*energyAvailable, 1.72);
373  double drange = range - range2;
374  LogTrace(me) << " energy left = " << energyAvailable <<
375  " eV, range2 = " << range2 << " cm, drange = " << drange << " cm";
376  if ( drange < clusterExtent ) {
377  theCrossGap->addElectronToBack(); // increment last element
378  }
379  else {
380  startHere += drange*theCrossGap->unitVector(); // update delta e start position
381  range = range2; // update range
382  new_range = true; // Test range again
383  LogTrace(me) << "reset range to range2 = " << range
384  << " from startHere = " << startHere << " and iterate";
385  }
386  break; // out of inner while energyAvailable>eion
387 
388  } // inner while energyAvailable>eion
389 
390  } // while !new_range && energyAvailable>ework
391 
392  // energyAvailable now less than ework, but still may be over eion...add an e
393  if ( energyAvailable > eion ) {
394  energyAvailable -= ework; // yes, it may go negative
395  theCrossGap->addElectronToBack(); // add one more e
396  }
397 
398  } // if range
399 
400  } // outer while energyAvailable>eion
401 }
402 
403 void CSCGasCollisions::writeSummary( int n_try, int n_steps, double sum_steps, float dedx,
404  const PSimHit& simhit ) const
405 {
406  // Switched from std::cout to LogVerbatim, Mar 2015
407 
408  std::vector<LocalPoint> ion_clusters = theCrossGap->ionClusters();
409  std::vector<int> electrons = theCrossGap->electrons();
410  std::vector<float> elosses = theCrossGap->eLossPerStep();
411  std::vector<double> steps = theCrossGap->stepLengths();
412 
413  edm::LogVerbatim("CSCDigitizer") << "------------------"
414  << "\nAFTER CROSSING GAP";
415  /*
416  edm::LogVerbatim("CSCDigitizer") << "no. of steps tried = " << n_try
417  << "\nno. of steps from theCrossGap = " << theCrossGap->noOfSteps()
418  << "\nsize of steps vector = " << steps.size();
419 
420  edm::LogVerbatim("CSCDigitizer") << "no. of collisions (steps) = " << n_steps
421  << "\nsize of elosses vector = " << elosses.size()
422  << "\nsize of ion clusters vector = " << ion_clusters.size()
423  << "\nsize of electrons vector = " << 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 gap - unless we reach MAX_STEPS
428  if ( (nsteps != MAX_STEPS) && (nsteps != mstep) ) {
429  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! no. of steps = " << nsteps
430  << " .ne. steps.size()-1 = " << mstep;
431  }
432  size_t meloss = elosses.size();
433 
434  if ( nsteps != meloss ){
435  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] WARNING! no. of steps = " << nsteps
436  << " .ne. no. of elosses = " << meloss;
437  }
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  }
452  else {
453  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] # / postion of 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(), bind2nd(greater<float>(), 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 ) edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] average step length = " <<
464  sum_steps/float(nsteps) << " cm";
465  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] total energy loss across gap = " <<
466  dedx << " eV = " << dedx/1000. << " keV";
467  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] no. of primary ionizing collisions across gap = no. with eloss > eion = "
468  << eion << " eV = " << n_ic;
469  if ( nsteps > 0 ) edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] average energy loss/collision = " <<
470  dedx/float(nsteps) << " eV";
471 
472  std::vector<int>::const_iterator bigger = find(electrons.begin(), electrons.end(), 0 );
473  if ( bigger != electrons.end() ) {
474  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] TROUBLE! There is a cluster with 0 electrons.";
475  }
476 
477  int n_e = accumulate(electrons.begin(), electrons.end(), 0 );
478 
479  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] SUMMARY: simhit"
480  << "\n entry = " << simhit.entryPoint()
481  << ": exit = " << simhit.exitPoint()
482  << "\n particle type = " << simhit.particleType() << " : momentum = " << simhit.pabs()
483  << " GeV/c : energy loss = " << simhit.energyLoss()*1.E06 << " keV";
484 
485  if ( nsteps > 0 ) {
486  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisions] SUMMARY: ionization"
487  << " : steps= " << nsteps << " : sum(steps)= " << sum_steps << " cm : <step>= "
488  << sum_steps/float(nsteps) << " cm"
489  << " : ionizing= " << n_ic << " : ionclus= " << mclus << " : total e= " << n_e
490  << " : <dedx/step>= " << dedx/float(nsteps)
491  << " eV : <e/ionclus>= " << float(n_e)/float(mclus)
492  << " : dedx= " << dedx/1000. << " keV";
493  }
494  else {
495  edm::LogVerbatim("CSCDigitizer") << "[CSCGasCollisons] ERROR? no collision steps!";
496  }
497 
498  // Turn off output file - used for initial development
499  // if ( saveGasCollisions() ) {
500  // ofstream of0("osteplen.dat",ios::app);
501  // std::copy( steps.begin(), steps.end(), std::ostream_iterator<float>(of0,"\n"));
502 
503  // ofstream of1("olperc.dat",ios::app);
504  // std::copy( elosses.begin(), elosses.end(), std::ostream_iterator<float>(of1,"\n"));
505 
506  // ofstream of2("oclpos.dat",ios::app);
507  // std::copy( ion_clusters.begin(), ion_clusters.end(), std::ostream_iterator<LocalPoint>(of2,"\n"));
508 
509  // ofstream of3("oepercl.dat",ios::app);
510  // std::copy( electrons.begin(), electrons.end(), std::ostream_iterator<int>(of3,"\n"));
511  // }
512 
513 }
514 
515 float CSCGasCollisions::lnEnergyLoss( float lnCollisions,
516  const std::vector<float>& collisions ) const {
517 
518  float lnE = -1.;
519 
520  // Find collision[] bin in which lnCollisions falls
521  std::vector<float>::const_iterator it = find(collisions.begin(),
522  collisions.end(), lnCollisions );
523 
524  if ( it != collisions.end() ) {
525  // found the value
526  std::vector<float>::difference_type ihi = it - collisions.begin();
527  LogTrace(me) << ": using one energy bin " << ihi << " = "
528  << theEnergyBins[ihi]
529  << " for lnCollisions = " << lnCollisions;
530  lnE = theEnergyBins[ihi];
531  }
532  else {
533  // interpolate the value
534  std::vector<float>::const_iterator loside = find_if(collisions.begin(),
535  collisions.end(), bind2nd(less<float>(), lnCollisions));
536  std::vector<float>::difference_type ilo = loside - collisions.begin();
537  if ( ilo > 0 ) {
538  LogTrace(me) << ": using energy bin "
539  << ilo-1 << " and " << ilo;
540  lnE = theEnergyBins[ilo-1] + (lnCollisions-collisions[ilo-1])*
541  (theEnergyBins[ilo]-theEnergyBins[ilo-1]) /
542  (collisions[ilo]-collisions[ilo-1]);
543  }
544  else {
545  LogTrace(me) << ": using one energy bin 0 = "
546  << theEnergyBins[0]
547  << " for lnCollisions = " << lnCollisions;
548  lnE = theEnergyBins[0]; //@@ WHAT ELSE TO DO?
549  }
550  }
551 
552  return lnE;
553 }
554 
556  std::vector<float>& collisions ) const
557 {
558  std::vector<float>::const_iterator bigger = find_if(theGammaBins.begin(),
559  theGammaBins.end(), bind2nd(greater<float>(), logGamma));
560 
561  if ( bigger == theGammaBins.end() ) {
562  // use highest bin
563  LogTrace(me) << ": using highest gamma bin"
564  << " for logGamma = " << logGamma;
565  for (int i=0; i<N_ENERGY; ++i)
566  collisions[i] = theCollisionTable[i*N_GAMMA];
567  }
568  else {
569  // use bigger and its lower neighbour
570  std::vector<float>::difference_type ihi = bigger - theGammaBins.begin();
571  if ( ihi > 0 ) {
572  double dlg2 = *bigger--; // and decrement after deref
573  //LogTrace(me) << ": using gamma bins "
574  // << ihi-1 << " and " << ihi;
575  double dlg1 = *bigger; // now the preceding element
576  double dlg = (logGamma-dlg1)/(dlg2-dlg1);
577  double omdlg = 1. - dlg;
578  for (int i=0; i<N_ENERGY; ++i)
579  collisions[i] = theCollisionTable[i*N_GAMMA+ihi-1]*omdlg +
580  theCollisionTable[i*N_GAMMA+ihi]*dlg;
581  }
582  else {
583  // bigger has no lower neighbour
584  LogTrace(me) << ": using lowest gamma bin"
585  << " for logGamma = " << logGamma;
586 
587  for (int i=0; i<N_ENERGY; ++i)
588  collisions[i] = theCollisionTable[i*N_GAMMA];
589  }
590  }
591 }
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
const String & name() const
return full name
Definition: FileInPath.h:34
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