CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes
CSCGasCollisions Class Reference

#include <CSCGasCollisions.h>

Public Member Functions

 CSCGasCollisions ()
 
void setParticleDataTable (const ParticleDataTable *pdt)
 
void setRandomEngine (CLHEP::HepRandomEngine &engine)
 
void simulate (const PSimHit &, std::vector< LocalPoint > &clusters, std::vector< int > &electrons)
 
virtual ~CSCGasCollisions ()
 

Static Public Attributes

static const int MAX_STEPS = 400
 
static const int N_ENERGY = 63
 
static const int N_ENTRIES = N_GAMMA*N_ENERGY
 
static const int N_GAMMA = 21
 

Private Member Functions

void fillCollisionsForThisGamma (float, std::vector< float > &) const
 
float generateEnergyLoss (double avCollisions, double anmin, double anmax, const std::vector< float > &collisions) const
 
double generateStep (double avCollisions) const
 
void ionize (double energyTransferred, LocalPoint startHere) const
 
float lnEnergyLoss (float, const std::vector< float > &) const
 
void readCollisionTable ()
 
void writeSummary (int n_steps, double sum_steps, float dedx, float simHiteloss) const
 

Private Attributes

double clusterExtent
 
double deCut
 
double eion
 
double ework
 
double gasDensity
 
const std::string me
 
bool saveGasCollisions
 
std::vector< float > theCollisionTable
 
CSCCrossGaptheCrossGap
 
std::vector< float > theEnergyBins
 
std::vector< float > theGammaBins
 
const ParticleDataTabletheParticleDataTable
 
CLHEP::RandExponential * theRandExponential
 
CLHEP::RandFlat * theRandFlat
 

Detailed Description

Class to handle CSC gas ionization simulation during digitization.

Author
Tim Cox

This class contains the port of the old CMSIM/ORCA Fortran routine mc_clus.F, and the GEANT3 code it relied on. Of course the version here is much improved :)

Definition at line 24 of file CSCGasCollisions.h.

Constructor & Destructor Documentation

CSCGasCollisions::CSCGasCollisions ( )

Definition at line 47 of file CSCGasCollisions.cc.

References clusterExtent, deCut, eion, ework, gasDensity, me, readCollisionTable(), and saveGasCollisions.

47  : me("CSCGasCollisions"),
48  gasDensity( 2.1416e-03 ),
49  deCut( 1.e05 ), eion( 14.95 ), ework( 34.0 ), clusterExtent( 0.001 ),
53  theRandFlat(0),
55  saveGasCollisions ( false )
56 {
57 
58  edm::LogInfo(me) << "Constructing a " << me << ":";
59  edm::LogInfo(me) << "gas density = " << gasDensity << " g/cm3";
60  edm::LogInfo(me) << "max eloss per collision allowed = " << deCut/1000.
61  << " keV (for higher elosses, hits should have been simulated.)";
62  edm::LogInfo(me) << "ionization threshold = " << eion << " eV";
63  edm::LogInfo(me) << "effective work function = " << ework << " eV";
64  edm::LogInfo(me) << "cluster extent = " << clusterExtent*1.e04 << " micrometres";
65  edm::LogInfo(me) << "Save gas collision info to files when using debugV? " << saveGasCollisions;
66 
68 }
static const int N_ENERGY
std::vector< float > theGammaBins
const std::string me
const ParticleDataTable * theParticleDataTable
static const int N_GAMMA
static const int N_ENTRIES
CLHEP::RandExponential * theRandExponential
std::vector< float > theCollisionTable
std::vector< float > theEnergyBins
CSCCrossGap * theCrossGap
CLHEP::RandFlat * theRandFlat
CSCGasCollisions::~CSCGasCollisions ( )
virtual

Definition at line 70 of file CSCGasCollisions.cc.

References me, theCrossGap, theRandExponential, and theRandFlat.

70  {
71  edm::LogInfo(me) << "Destructing a " << me;
72  delete theCrossGap;
73  delete theRandFlat;
74  delete theRandExponential;
75 }
const std::string me
CLHEP::RandExponential * theRandExponential
CSCCrossGap * theCrossGap
CLHEP::RandFlat * theRandFlat

Member Function Documentation

void CSCGasCollisions::fillCollisionsForThisGamma ( float  logGamma,
std::vector< float > &  collisions 
) const
private

Definition at line 510 of file CSCGasCollisions.cc.

References i, LogTrace, me, N_ENERGY, N_GAMMA, theCollisionTable, and theGammaBins.

Referenced by simulate().

512 {
513  std::vector<float>::const_iterator bigger = find_if(theGammaBins.begin(),
514  theGammaBins.end(), bind2nd(greater<float>(), logGamma));
515 
516  if ( bigger == theGammaBins.end() ) {
517  // use highest bin
518  LogTrace(me) << ": using highest gamma bin"
519  << " for logGamma = " << logGamma;
520  for (int i=0; i<N_ENERGY; ++i)
521  collisions[i] = theCollisionTable[i*N_GAMMA];
522  }
523  else {
524  // use bigger and its lower neighbour
525  std::vector<float>::difference_type ihi = bigger - theGammaBins.begin();
526  if ( ihi > 0 ) {
527  double dlg2 = *bigger--; // and decrement after deref
528  //LogTrace(me) << ": using gamma bins "
529  // << ihi-1 << " and " << ihi;
530  double dlg1 = *bigger; // now the preceding element
531  double dlg = (logGamma-dlg1)/(dlg2-dlg1);
532  double omdlg = 1. - dlg;
533  for (int i=0; i<N_ENERGY; ++i)
534  collisions[i] = theCollisionTable[i*N_GAMMA+ihi-1]*omdlg +
535  theCollisionTable[i*N_GAMMA+ihi]*dlg;
536  }
537  else {
538  // bigger has no lower neighbour
539  LogTrace(me) << ": using lowest gamma bin"
540  << " for logGamma = " << logGamma;
541 
542  for (int i=0; i<N_ENERGY; ++i)
543  collisions[i] = theCollisionTable[i*N_GAMMA];
544  }
545  }
546 }
int i
Definition: DBlmapReader.cc:9
static const int N_ENERGY
std::vector< float > theGammaBins
const std::string me
static const int N_GAMMA
#define LogTrace(id)
std::vector< float > theCollisionTable
float CSCGasCollisions::generateEnergyLoss ( double  avCollisions,
double  anmin,
double  anmax,
const std::vector< float > &  collisions 
) const
private

Definition at line 292 of file CSCGasCollisions.cc.

References CSCCrossGap::beta2(), create_public_lumi_plots::exp, CSCCrossGap::gamma(), lnEnergyLoss(), create_public_lumi_plots::log, LogTrace, me, theCrossGap, and theRandFlat.

Referenced by simulate().

294 {
295 // Generate a no. of collisions between collisions[0] and [N_ENERGY-1]
296  float lnColl = log(theRandFlat->fire(anmin, anmax));
297 
298 // Without using CLHEP: approx random between anmin and anmax
299 // double ra = double(rand())/double(RAND_MAX)*avCollisions;
300 // cout << "ra = " << ra << std::endl;
301 // float lnColl = static_cast<float>( log( ra ) );
302 
303  // Find energy loss for that number
304  float lnE = lnEnergyLoss( lnColl, collisions );
305  float eloss = exp(lnE);
306  // Compensate if gamma was actually below 1.1
307  if ( theCrossGap->gamma() < 1.1 ) eloss = eloss * 0.173554/theCrossGap->beta2();
308  LogTrace(me) << "eloss = " << eloss;
309  // Next line only used to fill container of eloss's for later diagnostic dumps
310  // if ( debugV ) theCrossGap->addEloss( eloss );
311  return eloss;
312 }
const std::string me
float lnEnergyLoss(float, const std::vector< float > &) const
double gamma() const
Definition: CSCCrossGap.h:52
#define LogTrace(id)
double beta2() const
Definition: CSCCrossGap.h:51
CSCCrossGap * theCrossGap
CLHEP::RandFlat * theRandFlat
double CSCGasCollisions::generateStep ( double  avCollisions) const
private

Definition at line 277 of file CSCGasCollisions.cc.

References LogTrace, me, relval_parameters_module::step, and theRandExponential.

Referenced by simulate().

278 {
279 // Generate a m.f.p. (1/avCollisions = cm/collision)
280  double step = (theRandExponential->fire())/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;
287  // Next line only used to fill a container of 'step's for later diagnostic dumps
288  //if ( debugV ) theCrossGap->addStep( step );
289  return step;
290 }
const std::string me
#define LogTrace(id)
CLHEP::RandExponential * theRandExponential
void CSCGasCollisions::ionize ( double  energyTransferred,
LocalPoint  startHere 
) const
private

Definition at line 314 of file CSCGasCollisions.cc.

References CSCCrossGap::addCluster(), CSCCrossGap::addElectrons(), CSCCrossGap::addElectronToBack(), clusterExtent, eion, ework, gasDensity, LogTrace, me, CSCCrossGap::noOfClusters(), CSCCrossGap::noOfElectrons(), funct::pow(), theCrossGap, and CSCCrossGap::unitVector().

Referenced by simulate().

315 {
316  while ( energyAvailable > eion ) {
317  LogTrace(me) << " NEW CLUSTER " << theCrossGap->noOfClusters() + 1 <<
318  " AT " << startHere;
319  LocalPoint newCluster( startHere );
320  theCrossGap->addCluster(newCluster);
321 
322  //@@ I consider NOT subtracting eion before calculating range to be a bug.
323  //@@ But this changes tuning of the algorithm so leave it until after the big rush to 7_5_0
324  //@@ energyAvailable -= eion;
325 
326  // Sauli CERN 77-09: delta e range with E in MeV (Sauli references Kobetich & Katz 1968,
327  // but I cannot find this expression in that set of papers.)
328  // Take HALF that range. //@@ Why? Why not...
329  double range = 0.5 * (0.71/gasDensity)*pow( energyAvailable*1.E-6, 1.72);
330  LogTrace(me) << " range = " << range;
331  if ( range < clusterExtent ) {
332 
333  // short-range delta e
334  // How many electrons can we make? Now use *average* energy for ionization (not *minimum*)
335  int nelec = static_cast<int>(energyAvailable/ework);
336  LogTrace(me) << "s-r delta energy in = " << energyAvailable;
337  //energyAvailable -= nelec*(energyAvailable/ework);
338  energyAvailable -= nelec*ework;
339  // If still above eion (minimum, not average) add one more e
340  if ( energyAvailable > eion ) {
341  ++nelec;
342  energyAvailable -= eion;
343  }
344  LogTrace(me) << "s-r delta energy out = " << energyAvailable << ", nelec = " << nelec;
345  theCrossGap->addElectrons( nelec );
346  break;
347 
348  }
349  else {
350  // long-range delta e
351  LogTrace(me) << "l-r delta \n"
352  << "no. of electrons in cluster now = " << theCrossGap->noOfElectrons();
353  theCrossGap->addElectrons( 1 ); // Position is at startHere still
354 
355  bool new_range = false;
356  while ( !new_range && (energyAvailable>ework) ) {
357  energyAvailable -= ework;
358  while ( energyAvailable > eion ) {
359  double range2 = 0.5 * 0.71/gasDensity*pow( 1.E-6*energyAvailable, 1.72);
360  double drange = range - range2;
361  LogTrace(me) << " energy left = " << energyAvailable <<
362  ", range2 = " << range2 << ", drange = " << drange;
363  if ( drange < clusterExtent ) {
364  theCrossGap->addElectronToBack(); // increment last element
365  }
366  else {
367  startHere += drange*theCrossGap->unitVector(); // update delta e start position
368  range = range2; // update range
369  new_range = true; // Test range again
370  LogTrace(me) << "reset range to range2 and iterate";
371  }
372  break; // out of inner while energyAvailable>eion
373 
374  } // inner while energyAvailable>eion
375 
376  } // while !new_range && energyAvailable>ework
377 
378  // energyAvailable now less than ework, but still may be over eion...add an e
379  if ( energyAvailable > eion ) {
380  energyAvailable -= ework; // yes, it may go negative
381  theCrossGap->addElectronToBack(); // add one more e
382  }
383 
384  } // if range
385 
386  } // outer while energyAvailable>eion
387 }
int noOfClusters() const
Definition: CSCCrossGap.h:34
void addElectronToBack()
Definition: CSCCrossGap.h:44
const std::string me
LocalVector unitVector() const
Definition: CSCCrossGap.h:54
int noOfElectrons() const
Definition: CSCCrossGap.h:36
#define LogTrace(id)
CSCCrossGap * theCrossGap
void addCluster(LocalPoint here)
Definition: CSCCrossGap.h:42
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void addElectrons(int nelec=1)
Definition: CSCCrossGap.h:43
float CSCGasCollisions::lnEnergyLoss ( float  lnCollisions,
const std::vector< float > &  collisions 
) const
private

Definition at line 470 of file CSCGasCollisions.cc.

References spr::find(), LogTrace, me, and theEnergyBins.

Referenced by generateEnergyLoss().

471  {
472 
473  float lnE = -1.;
474 
475  // Find collision[] bin in which lnCollisions falls
476  std::vector<float>::const_iterator it = find(collisions.begin(),
477  collisions.end(), lnCollisions );
478 
479  if ( it != collisions.end() ) {
480  // found the value
481  std::vector<float>::difference_type ihi = it - collisions.begin();
482  LogTrace(me) << ": using one energy bin " << ihi << " = "
483  << theEnergyBins[ihi]
484  << " for lnCollisions = " << lnCollisions;
485  lnE = theEnergyBins[ihi];
486  }
487  else {
488  // interpolate the value
489  std::vector<float>::const_iterator loside = find_if(collisions.begin(),
490  collisions.end(), bind2nd(less<float>(), lnCollisions));
491  std::vector<float>::difference_type ilo = loside - collisions.begin();
492  if ( ilo > 0 ) {
493  LogTrace(me) << ": using energy bin "
494  << ilo-1 << " and " << ilo;
495  lnE = theEnergyBins[ilo-1] + (lnCollisions-collisions[ilo-1])*
496  (theEnergyBins[ilo]-theEnergyBins[ilo-1]) /
497  (collisions[ilo]-collisions[ilo-1]);
498  }
499  else {
500  LogTrace(me) << ": using one energy bin 0 = "
501  << theEnergyBins[0]
502  << " for lnCollisions = " << lnCollisions;
503  lnE = theEnergyBins[0]; //@@ WHAT ELSE TO DO?
504  }
505  }
506 
507  return lnE;
508 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
const std::string me
#define LogTrace(id)
std::vector< float > theEnergyBins
void CSCGasCollisions::readCollisionTable ( )
private

Definition at line 77 of file CSCGasCollisions.cc.

References edm::hlt::Exception, python.connectstrParser::f1, groupFilesInBlocks::fin, i, j, LogTrace, me, N_ENERGY, N_ENTRIES, N_GAMMA, FileInPath::name(), getHLTPrescaleColumns::path, theCollisionTable, theEnergyBins, and theGammaBins.

Referenced by CSCGasCollisions().

77  {
78 
79  // I'd prefer to allow comments in the data file which means
80  // complications of reading line-by-line and then item-by-item.
81 
82  // Is float OK? Or do I need double?
83 
84  // Do I need any sort of error trapping?
85 
86  // We use the default CMSSW data file path
87  // We use the default file name 'collisions.dat'
88 
89  // This can be reset in .orcarc by SimpleConfigurable
90  // Muon:Endcap:CollisionsFile
91 
92  string path( getenv( "CMSSW_SEARCH_PATH" ) );
93  // TODO make configurable
94  string colliFile = "SimMuon/CSCDigitizer/data/collisions.dat";
95  FileInPath f1( path, colliFile );
96  if (f1() == 0 ) {
97  string errorMessage = "Input file " + colliFile + "not found";
98  edm::LogError("CSCGasCollisions") << errorMessage << " in path " << path
99  << "\nSet Muon:Endcap:CollisionsFile in .orcarc to the "
100  " location of the file relative to ORCA_DATA_PATH." ;
101  throw cms::Exception( " Endcap Muon gas collisions data file not found.");
102  }
103  else {
104  edm::LogInfo(me) << ": reading " << f1.name();
105  }
106 
107  ifstream & fin = *f1();
108 
109  if (fin == 0) {
110  string errorMessage = "Cannot open input file " + path + colliFile;
111  edm::LogError("CSCGasCollisions") << errorMessage;
112  throw cms::Exception(errorMessage);
113  }
114 
115  fin.clear( ); // Clear eof read status
116  fin.seekg( 0, ios::beg ); // Position at start of file
117 
118  // @@ We had better have the right sizes everywhere or all
119  // hell will break loose. There's no trapping.
120  LogTrace(me) << "Reading gamma bins";
121  int j = 0;
122  for(int i = 0; i<N_GAMMA; ++i) {
123  fin >> theGammaBins[i];
124  LogTrace(me) << ++j << " " << theGammaBins[i];
125  }
126 
127  j = 0;
128  LogTrace(me) << "Reading energy bins \n";
129  for(int i = 0; i<N_ENERGY; ++i) {
130  fin >> theEnergyBins[i];
131  LogTrace(me) << ++j << " " << theEnergyBins[i];
132  }
133 
134  j = 0;
135  LogTrace(me) << "Reading collisions table \n";
136 
137  for(int i = 0; i<N_ENTRIES; ++i) {
138  fin >> theCollisionTable[i];
139  LogTrace(me) << ++j << " " << theCollisionTable[i];
140  }
141 
142  fin.close();
143 }
int i
Definition: DBlmapReader.cc:9
static const int N_ENERGY
std::vector< float > theGammaBins
const std::string me
int j
Definition: DBlmapReader.cc:9
static const int N_GAMMA
#define LogTrace(id)
static const int N_ENTRIES
std::vector< float > theCollisionTable
std::vector< float > theEnergyBins
void CSCGasCollisions::setParticleDataTable ( const ParticleDataTable pdt)

Definition at line 146 of file CSCGasCollisions.cc.

References theParticleDataTable.

Referenced by CSCWireHitSim::setParticleDataTable().

147 {
148  theParticleDataTable = pdt;
149 }
const ParticleDataTable * theParticleDataTable
void CSCGasCollisions::setRandomEngine ( CLHEP::HepRandomEngine &  engine)

Definition at line 152 of file CSCGasCollisions.cc.

References theRandExponential, and theRandFlat.

Referenced by CSCWireHitSim::setRandomEngine().

153 {
154  theRandFlat = new CLHEP::RandFlat(engine);
155  theRandExponential = new CLHEP::RandExponential(engine);
156 }
CLHEP::RandExponential * theRandExponential
CLHEP::RandFlat * theRandFlat
void CSCGasCollisions::simulate ( const PSimHit simHit,
std::vector< LocalPoint > &  clusters,
std::vector< int > &  electrons 
)

Definition at line 159 of file CSCGasCollisions.cc.

References deCut, eion, CSCCrossGap::electrons(), PSimHit::entryPoint(), PSimHit::exitPoint(), create_public_lumi_plots::exp, fillCollisionsForThisGamma(), generateEnergyLoss(), generateStep(), CSCCrossGap::ionClusters(), ionize(), CSCCrossGap::length(), CSCCrossGap::logGamma(), LogTrace, MAX_STEPS, me, N_ENERGY, PSimHit::pabs(), PSimHit::particleType(), relval_parameters_module::step, theCrossGap, theParticleDataTable, and CSCCrossGap::unitVector().

Referenced by CSCWireHitSim::getIonizationClusters().

160  {
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();
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  LogTrace(me) << ": WARNING! simhit entry and exit are very close: \n"
191  << "\n entry = " << simHit.entryPoint()
192  << "\n exit = " << simHit.exitPoint()
193  << "\n particle type = " << iam << ", momentum = " << mom
194  << ", gap length = " << gapSize;
195  return; //@@ Just skip this PSimHit
196  }
197 
198  // Interpolate the table for current gamma value
199  // Extract collisions binned by energy loss values, for this gamma
200  std::vector<float> collisions(N_ENERGY);
201  double loggam = theCrossGap->logGamma();
202  fillCollisionsForThisGamma( static_cast<float>(loggam), collisions );
203 
204  double anmin = exp(collisions[N_ENERGY-1]);
205  double anmax = exp(collisions[0]);
206  double amu = anmax - anmin;
207 
208  LogTrace(me) << "collisions extremes = " << collisions[N_ENERGY-1]
209  << ", " << collisions[0] << "\n"
210  << "anmin = " << anmin << ", anmax = " << anmax << "\n"
211  << "amu = " << amu << "\n";
212 
213  float dedx = 0.; // total energy loss
214  double sum_steps = 0.; // total distance across gap (along simhit direction)
215  int n_steps = 0; // no. of steps/primary collisions
216  int n_try = 0; // no. of tries to generate steps
217  double step = -1.; // Sentinel for start
218 
219  LocalPoint layerLocalPoint( simHit.entryPoint() );
220 
221  // step/primary collision loop
222  while ( sum_steps < gapSize) {
223  ++n_try;
224  if ( n_try > MAX_STEPS ) {
225  LogTrace(me) << ": n_try=" << n_try << " is too large. Skip simhit."
226  << "\n particle type=" << iam << ", momentum= " << mom
227  << "\n gapSize=" << gapSize << ", last step=" << step
228  << ", sum_steps=" << sum_steps << ", n_steps=" << n_steps;
229  break;
230  }
231  step = generateStep( amu );
232  if ( sum_steps + step > gapSize ) break;
233 
234  float eloss = generateEnergyLoss( amu, anmin, anmax, collisions );
235 
236  // Is the eloss too large? (then GEANT should have produced hits!)
237  if ( eloss > deCut ) {
238  LogTrace(me) << "eloss > " << deCut << " = " << eloss;
239  continue; // to generate another collision/step
240  }
241 
242  dedx += eloss; // the energy lost from the ionizing particle
243  sum_steps += step; // the position of the ionizing particle
244  ++n_steps; // the number of primary collisions
245 
246  if (n_steps > MAX_STEPS ) { // Extra-careful trap for bizarreness
247  edm::LogInfo(me) << ": n_steps=" << n_steps << " is too large. Skip simhit.";
248  edm::LogInfo(me) << "particle type=" << iam << ", momentum= " << mom;
249  edm::LogInfo(me) << "gapSize=" << gapSize << ", last step=" << step
250  << ", sum_steps=" << sum_steps;
251  break;
252  }
253  LogTrace(me) << "sum_steps = " << sum_steps << ", dedx = " << dedx;
254 
255  // Generate ionization.
256  // eion is the minimum energy at which ionization can occur in the gas
257  if ( eloss > eion ) {
258  layerLocalPoint += step * theCrossGap->unitVector(); // local point where the collision occurs
259  ionize( eloss, layerLocalPoint );
260  }
261  else {
262  LogTrace(me) << "Energy available = " << eloss <<
263  ", too low for ionization.";
264  }
265 
266  } // step/collision loop
267 
268  //TODO port this
269  //if ( debugV ) writeSummary( n_steps, sum_steps, dedx, simHit.energyLoss() );
270 
271  // Return values in two container arguments
272  positions = theCrossGap->ionClusters();
274  return;
275 }
std::vector< int > electrons() const
Definition: CSCCrossGap.h:35
static const int N_ENERGY
std::vector< LocalPoint > ionClusters() const
Definition: CSCCrossGap.h:33
const std::string me
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
const ParticleDataTable * theParticleDataTable
LocalVector unitVector() const
Definition: CSCCrossGap.h:54
static const int MAX_STEPS
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:63
HepPDT::ParticleData ParticleData
void ionize(double energyTransferred, LocalPoint startHere) const
#define LogTrace(id)
double logGamma(double mass, float momentum)
Definition: CSCCrossGap.cc:23
double generateStep(double avCollisions) const
float generateEnergyLoss(double avCollisions, double anmin, double anmax, const std::vector< float > &collisions) const
int particleType() const
Definition: PSimHit.h:85
float length() const
Definition: CSCCrossGap.h:55
void fillCollisionsForThisGamma(float, std::vector< float > &) const
CSCCrossGap * theCrossGap
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
void CSCGasCollisions::writeSummary ( int  n_steps,
double  sum_steps,
float  dedx,
float  simHiteloss 
) const
private

Definition at line 389 of file CSCGasCollisions.cc.

References filterCSVwithJSON::copy, gather_cfg::cout, eion, CSCCrossGap::electrons(), HI_PhotonSkim_cff::electrons, CSCCrossGap::eLossPerStep(), spr::find(), CSCCrossGap::ionClusters(), CSCCrossGap::noOfSteps(), saveGasCollisions, CSCCrossGap::stepLengths(), relval_steps::steps, and theCrossGap.

390 {
391  std::vector<LocalPoint> ion_clusters = theCrossGap->ionClusters();
392  std::vector<int> electrons = theCrossGap->electrons();
393  std::vector<float> elosses = theCrossGap->eLossPerStep();
394  std::vector<double> steps = theCrossGap->stepLengths();
395 
396  cout << "------------------" << std::endl;
397  cout << "AFTER CROSSING GAP" << std::endl;
398  cout << "No. of steps = " << n_steps << std::endl;
399  cout << "Check: stored steps = " << theCrossGap->noOfSteps() << std::endl;
400 
401  cout << "Lengths of steps: " << std::endl;
402  std::copy( steps.begin(), steps.end(), std::ostream_iterator<float>(cout,"\n"));
403  cout << std::endl;
404 
405  if ( saveGasCollisions ) {
406  ofstream of0("osteplen.dat",ios::app);
407  std::copy( steps.begin(), steps.end(), std::ostream_iterator<float>(of0,"\n"));
408  }
409 
410  cout << "Total sum of steps = " << sum_steps << std::endl;
411  if ( n_steps > 0 ) cout << "Average step length = " <<
412  sum_steps/float(n_steps) << std::endl;
413  cout << std::endl;
414 
415  cout << "Energy loss per collision:" << std::endl;
416  std::copy( elosses.begin(), elosses.end(), std::ostream_iterator<float>(cout,"\n"));
417  cout << std::endl;
418 
419  if ( saveGasCollisions ) {
420  ofstream of1("olperc.dat",ios::app);
421  std::copy( elosses.begin(), elosses.end(), std::ostream_iterator<float>(of1,"\n"));
422  }
423 
424  cout << "Total energy loss across gap = " << dedx << " eV = " <<
425  dedx/1000. << " keV" << std::endl;
426  int n_ic = count_if( elosses.begin(), elosses.end(),
427  bind2nd(greater<float>(), eion) );
428  cout << "No. of primary ionizing collisions across gap = " << n_ic << std::endl;
429  if ( n_steps > 0 ) cout << "Average energy loss/collision = " <<
430  dedx/float(n_steps) << " eV" << std::endl;
431  cout << std::endl;
432 
433  cout << "No. of ion clusters = " << ion_clusters.size() << std::endl;
434  cout << "Positions of clusters:" << std::endl;
435  std::copy( ion_clusters.begin(), ion_clusters.end(),
436  std::ostream_iterator<LocalPoint>(cout,"\n"));
437  cout << std::endl;
438 
439  if ( saveGasCollisions ) {
440  ofstream of2("oclpos.dat",ios::app);
441  std::copy( ion_clusters.begin(), ion_clusters.end(),
442  std::ostream_iterator<LocalPoint>(of2,"\n"));
443  }
444 
445  cout << "No. of electrons per cluster:" << std::endl;
446  std::copy( electrons.begin(), electrons.end(), std::ostream_iterator<int>(cout,"\n"));
447  cout << std::endl;
448 
449  if ( saveGasCollisions ) {
450  ofstream of3("oepercl.dat",ios::app);
451  std::copy( electrons.begin(), electrons.end(), std::ostream_iterator<int>(of3,"\n"));
452  }
453 
454  // Check for zero-e clusters
455  std::vector<int>::const_iterator bigger = find(electrons.begin(),
456  electrons.end(), 0 );
457  if ( bigger != electrons.end() ) {
458  cout << "Error! There is a cluster with 0 electrons." << std::endl;
459  }
460  int n_e = accumulate(electrons.begin(), electrons.end(), 0 );
461  if ( n_steps > 0 ) {
462  cout << "# cm cm keV eV eV eV keV" << std::endl;
463  cout << " " << n_steps << " " << sum_steps << " " << sum_steps/float(n_steps) << " " <<
464  ion_clusters.size() << " " <<
465  dedx/1000. << " " << n_ic << " " << dedx/float(n_steps) << " " << n_e << " " <<
466  dedx/float(n_e) << " " << float(n_e)/float(ion_clusters.size()) << " " << simHiteloss*1.E6 << std::endl;
467  }
468 }
std::vector< int > electrons() const
Definition: CSCCrossGap.h:35
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:7
std::vector< float > eLossPerStep() const
Definition: CSCCrossGap.h:39
std::vector< double > stepLengths() const
Definition: CSCCrossGap.h:37
int noOfSteps() const
Definition: CSCCrossGap.h:38
CSCCrossGap * theCrossGap
tuple cout
Definition: gather_cfg.py:121

Member Data Documentation

double CSCGasCollisions::clusterExtent
private

Definition at line 64 of file CSCGasCollisions.h.

Referenced by CSCGasCollisions(), and ionize().

double CSCGasCollisions::deCut
private

Definition at line 61 of file CSCGasCollisions.h.

Referenced by CSCGasCollisions(), and simulate().

double CSCGasCollisions::eion
private

Definition at line 62 of file CSCGasCollisions.h.

Referenced by CSCGasCollisions(), ionize(), simulate(), and writeSummary().

double CSCGasCollisions::ework
private

Definition at line 63 of file CSCGasCollisions.h.

Referenced by CSCGasCollisions(), and ionize().

double CSCGasCollisions::gasDensity
private

Definition at line 55 of file CSCGasCollisions.h.

Referenced by CSCGasCollisions(), and ionize().

const int CSCGasCollisions::MAX_STEPS = 400
static

Definition at line 40 of file CSCGasCollisions.h.

Referenced by simulate().

const std::string CSCGasCollisions::me
private
const int CSCGasCollisions::N_ENERGY = 63
static

Definition at line 38 of file CSCGasCollisions.h.

Referenced by fillCollisionsForThisGamma(), readCollisionTable(), and simulate().

const int CSCGasCollisions::N_ENTRIES = N_GAMMA*N_ENERGY
static

Definition at line 39 of file CSCGasCollisions.h.

Referenced by readCollisionTable().

const int CSCGasCollisions::N_GAMMA = 21
static

Definition at line 37 of file CSCGasCollisions.h.

Referenced by fillCollisionsForThisGamma(), and readCollisionTable().

bool CSCGasCollisions::saveGasCollisions
private

Definition at line 74 of file CSCGasCollisions.h.

Referenced by CSCGasCollisions(), and writeSummary().

std::vector<float> CSCGasCollisions::theCollisionTable
private

Definition at line 68 of file CSCGasCollisions.h.

Referenced by fillCollisionsForThisGamma(), and readCollisionTable().

CSCCrossGap* CSCGasCollisions::theCrossGap
private
std::vector<float> CSCGasCollisions::theEnergyBins
private

Definition at line 67 of file CSCGasCollisions.h.

Referenced by lnEnergyLoss(), and readCollisionTable().

std::vector<float> CSCGasCollisions::theGammaBins
private

Definition at line 66 of file CSCGasCollisions.h.

Referenced by fillCollisionsForThisGamma(), and readCollisionTable().

const ParticleDataTable* CSCGasCollisions::theParticleDataTable
private

Definition at line 71 of file CSCGasCollisions.h.

Referenced by setParticleDataTable(), and simulate().

CLHEP::RandExponential* CSCGasCollisions::theRandExponential
private

Definition at line 73 of file CSCGasCollisions.h.

Referenced by generateStep(), setRandomEngine(), and ~CSCGasCollisions().

CLHEP::RandFlat* CSCGasCollisions::theRandFlat
private

Definition at line 72 of file CSCGasCollisions.h.

Referenced by generateEnergyLoss(), setRandomEngine(), and ~CSCGasCollisions().