22 #include "CLHEP/Random/RandExponential.h"
23 #include "CLHEP/Random/RandFlat.h"
55 :
me(
"CSCGasCollisions"),
56 gasDensity( 2.1416
e-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 )
69 <<
" keV (for higher elosses, hits should have been simulated.)";
99 string path( getenv(
"CMSSW_SEARCH_PATH" ) );
101 string colliFile =
"SimMuon/CSCDigitizer/data/collisions.dat";
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.");
114 ifstream &
fin = *
f1();
117 string errorMessage =
"Cannot open input file " + path + colliFile;
123 fin.seekg( 0, ios::beg );
142 LogTrace(
me) <<
"Reading collisions table \n";
146 LogTrace(
me) << ++j <<
" " << theCollisionTable[
i];
160 std::vector<LocalPoint>& positions, std::vector<int>&
electrons, CLHEP::HepRandomEngine* engine ) {
162 const float epsilonL = 0.01;
168 double mom = simhit.
pabs();
173 double mass = 0.105658;
181 mass = particle->mass();
189 if ( gapSize <= epsilonL ) {
191 <<
"[CSCGasCollisions] WARNING! simhit entry and exit are too close - skipping simhit:"
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)";
202 std::vector<float> collisions(
N_ENERGY);
207 double anmax =
exp(collisions[0]);
208 double amu = anmax - anmin;
211 <<
", " << collisions[0] <<
"\n"
212 <<
"anmin = " << anmin <<
", anmax = " << anmax <<
"\n"
213 <<
"amu = " << amu <<
"\n";
216 double sum_steps = 0.;
224 while ( sum_steps < gapSize) {
228 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! n_try = " << n_try
229 <<
" is > MAX_STEPS = " << maxst <<
" - skipping simhit:"
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;
239 if ( sum_steps + step > gapSize )
break;
244 if ( eloss >
deCut ) {
245 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! eloss = " << eloss
246 <<
" eV is too large (> " <<
deCut <<
" eV) - trying another collision";
255 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! " << n_steps
256 <<
" is too many steps - skipping simhit:"
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";
264 LogTrace(
me) <<
"sum_steps = " << sum_steps <<
" cm , dedx = " << dedx <<
" eV";
268 if ( eloss >
eion ) {
270 ionize( eloss, layerLocalPoint );
273 LogTrace(
me) <<
"Energy available = " << eloss <<
" eV is too low for ionization (< eion = " <<
eion <<
" eV)";
289 double step = (CLHEP::RandExponential::shoot(engine))/avCollisions;
302 const std::vector<float>& collisions, CLHEP::HepRandomEngine* engine )
const
305 float lnColl =
log(CLHEP::RandFlat::shoot(engine, anmin, anmax));
314 float eloss =
exp(lnE);
325 while ( energyAvailable >
eion ) {
342 double range = 0.5 * (0.71/
gasDensity)*
pow( energyAvailable*1.E-6, 1.72);
343 LogTrace(
me) <<
" range = " << range <<
" cm";
348 int nelec =
static_cast<int>(energyAvailable/
ework);
349 LogTrace(
me) <<
"short-range delta energy in = " << energyAvailable <<
" eV";
351 energyAvailable -= nelec*
ework;
353 if ( energyAvailable >
eion ) {
355 energyAvailable -=
eion;
357 LogTrace(
me) <<
"short-range delta energy out = " << energyAvailable <<
" eV, nelec = " << nelec;
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";
383 LogTrace(
me) <<
"reset range to range2 = " << range
384 <<
" from startHere = " << startHere <<
" and iterate";
393 if ( energyAvailable >
eion ) {
394 energyAvailable -=
ework;
414 <<
"\nAFTER CROSSING GAP";
426 size_t nsteps = n_steps;
427 size_t mstep = steps.size() - 1;
428 if ( (nsteps !=
MAX_STEPS) && (nsteps != mstep) ) {
429 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! no. of steps = " << nsteps
430 <<
" .ne. steps.size()-1 = " << mstep;
432 size_t meloss = elosses.size();
434 if ( nsteps != meloss ){
435 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! no. of steps = " << nsteps
436 <<
" .ne. no. of elosses = " << meloss;
439 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] # / length of step / energy loss per collision:";
440 for (
size_t i = 0;
i != nsteps; ++
i ) {
445 size_t mclus = ion_clusters.size();
446 size_t melec = electrons.size();
448 if ( mclus != melec ){
449 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! size of cluster vector = " << mclus
450 <<
" .ne. size of electrons vector = " << melec;
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];
459 int n_ic = count_if( elosses.begin(), elosses.end(), bind2nd(greater<float>(),
eion) );
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";
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.";
477 int n_e = accumulate(electrons.begin(), electrons.end(), 0 );
482 <<
"\n particle type = " << simhit.
particleType() <<
" : momentum = " << simhit.
pabs()
483 <<
" GeV/c : energy loss = " << simhit.
energyLoss()*1.E06 <<
" keV";
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";
495 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisons] ERROR? no collision steps!";
516 const std::vector<float>& collisions )
const {
521 std::vector<float>::const_iterator it =
find(collisions.begin(),
522 collisions.end(), lnCollisions );
524 if ( it != collisions.end() ) {
526 std::vector<float>::difference_type ihi = it - collisions.begin();
527 LogTrace(
me) <<
": using one energy bin " << ihi <<
" = "
529 <<
" for lnCollisions = " << lnCollisions;
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();
539 << ilo-1 <<
" and " << ilo;
540 lnE =
theEnergyBins[ilo-1] + (lnCollisions-collisions[ilo-1])*
542 (collisions[ilo]-collisions[ilo-1]);
547 <<
" for lnCollisions = " << lnCollisions;
556 std::vector<float>& collisions )
const
559 theGammaBins.end(), bind2nd(greater<float>(), logGamma));
564 <<
" for logGamma = " << logGamma;
570 std::vector<float>::difference_type ihi = bigger -
theGammaBins.begin();
572 double dlg2 = *bigger--;
576 double dlg = (logGamma-dlg1)/(dlg2-dlg1);
577 double omdlg = 1. - dlg;
585 <<
" for logGamma = " << logGamma;
T getUntrackedParameter(std::string const &, T const &) const
std::vector< int > electrons() const
HepPDT::ParticleDataTable ParticleDataTable
bool saveGasCollisions(void) const
void addEloss(float eloss)
static const int N_ENERGY
std::vector< LocalPoint > ionClusters() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
CSCGasCollisions(const edm::ParameterSet &pset)
std::vector< float > eLossPerStep() const
std::vector< double > stepLengths() const
std::vector< float > theGammaBins
float lnEnergyLoss(float, const std::vector< float > &) const
Local3DPoint exitPoint() const
Exit point in the local Det frame.
bool dumpGasCollisions(void) const
const ParticleDataTable * theParticleDataTable
LocalVector unitVector() const
int noOfElectrons() const
static const int MAX_STEPS
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
HepPDT::ParticleData ParticleData
void setParticleDataTable(const ParticleDataTable *pdt)
void ionize(double energyTransferred, LocalPoint startHere) const
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)
void readCollisionTable()
std::vector< float > theCollisionTable
double generateStep(double avCollisions, CLHEP::HepRandomEngine *) const
std::vector< float > theEnergyBins
const String & name() const
return full name
float energyLoss() const
The energy deposit in the PSimHit, in ???.
float generateEnergyLoss(double avCollisions, double anmin, double anmax, const std::vector< float > &collisions, CLHEP::HepRandomEngine *) const
void addStep(double step)
void fillCollisionsForThisGamma(float, std::vector< float > &) const
CSCCrossGap * theCrossGap
volatile std::atomic< bool > shutdown_flag false
void simulate(const PSimHit &, std::vector< LocalPoint > &clusters, std::vector< int > &electrons, CLHEP::HepRandomEngine *)
void addCluster(LocalPoint here)
virtual ~CSCGasCollisions()
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Power< A, B >::type pow(const A &a, const B &b)
void addElectrons(int nelec=1)