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(
nullptr ),
61 saveGasCollisions_ (
false ), dumpGasCollisions_(
false )
69 <<
" keV (for higher elosses, hits should have been simulated.)";
100 string colliFile =
"SimMuon/CSCDigitizer/data/collisions.dat";
104 ifstream
fin{
f1.fullPath() };
107 string errorMessage =
"Cannot open input file " +
f1.fullPath();
113 fin.seekg( 0, ios::beg );
132 LogTrace(
me) <<
"Reading collisions table \n";
136 LogTrace(
me) << ++j <<
" " << theCollisionTable[
i];
150 std::vector<LocalPoint>& positions, std::vector<int>&
electrons, CLHEP::HepRandomEngine* engine ) {
152 const float epsilonL = 0.01;
158 double mom = simhit.
pabs();
163 double mass = 0.105658;
164 if(particle ==
nullptr)
171 mass = particle->mass();
179 if ( gapSize <= epsilonL ) {
181 <<
"[CSCGasCollisions] WARNING! simhit entry and exit are too close - skipping simhit:" 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)";
192 std::vector<float> collisions(
N_ENERGY);
197 double anmax =
exp(collisions[0]);
198 double amu = anmax - anmin;
201 <<
", " << collisions[0] <<
"\n" 202 <<
"anmin = " << anmin <<
", anmax = " << anmax <<
"\n" 203 <<
"amu = " << amu <<
"\n";
206 double sum_steps = 0.;
214 while ( sum_steps < gapSize) {
218 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! n_try = " << n_try
219 <<
" is > MAX_STEPS = " << maxst <<
" - skipping simhit:" 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;
229 if ( sum_steps + step > gapSize )
break;
234 if ( eloss >
deCut ) {
235 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! eloss = " << eloss
236 <<
" eV is too large (> " <<
deCut <<
" eV) - trying another collision";
245 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! " << n_steps
246 <<
" is too many steps - skipping simhit:" 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";
254 LogTrace(
me) <<
"sum_steps = " << sum_steps <<
" cm , dedx = " << dedx <<
" eV";
258 if ( eloss >
eion ) {
260 ionize( eloss, layerLocalPoint );
263 LogTrace(
me) <<
"Energy available = " << eloss <<
" eV is too low for ionization (< eion = " <<
eion <<
" eV)";
279 double step = (CLHEP::RandExponential::shoot(engine))/avCollisions;
292 const std::vector<float>& collisions, CLHEP::HepRandomEngine* engine )
const 295 float lnColl =
log(CLHEP::RandFlat::shoot(engine, anmin, anmax));
304 float eloss =
exp(lnE);
315 while ( energyAvailable >
eion ) {
332 double range = 0.5 * (0.71/
gasDensity)*
pow( energyAvailable*1.E-6, 1.72);
333 LogTrace(
me) <<
" range = " << range <<
" cm";
338 int nelec =
static_cast<int>(energyAvailable/
ework);
339 LogTrace(
me) <<
"short-range delta energy in = " << energyAvailable <<
" eV";
341 energyAvailable -= nelec*
ework;
343 if ( energyAvailable >
eion ) {
345 energyAvailable -=
eion;
347 LogTrace(
me) <<
"short-range delta energy out = " << energyAvailable <<
" eV, nelec = " << nelec;
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";
373 LogTrace(
me) <<
"reset range to range2 = " << range
374 <<
" from startHere = " << startHere <<
" and iterate";
383 if ( energyAvailable >
eion ) {
384 energyAvailable -=
ework;
404 <<
"\nAFTER CROSSING GAP";
416 size_t nsteps = n_steps;
417 size_t mstep = steps.size() - 1;
418 if ( (nsteps !=
MAX_STEPS) && (nsteps != mstep) ) {
419 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! no. of steps = " << nsteps
420 <<
" .ne. steps.size()-1 = " << mstep;
422 size_t meloss = elosses.size();
424 if ( nsteps != meloss ){
425 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! no. of steps = " << nsteps
426 <<
" .ne. no. of elosses = " << meloss;
429 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] # / length of step / energy loss per collision:";
430 for (
size_t i = 0;
i != nsteps; ++
i ) {
435 size_t mclus = ion_clusters.size();
436 size_t melec = electrons.size();
438 if ( mclus != melec ){
439 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! size of cluster vector = " << mclus
440 <<
" .ne. size of electrons vector = " << melec;
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];
449 int n_ic = count_if( elosses.begin(), elosses.end(), [&](
auto c){
return c > this->
eion;} );
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";
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.";
467 int n_e = accumulate(electrons.begin(), electrons.end(), 0 );
472 <<
"\n particle type = " << simhit.
particleType() <<
" : momentum = " << simhit.
pabs()
473 <<
" GeV/c : energy loss = " << simhit.
energyLoss()*1.E06 <<
" keV";
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";
485 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisons] ERROR? no collision steps!";
506 const std::vector<float>& collisions )
const {
511 std::vector<float>::const_iterator it =
find(collisions.begin(),
512 collisions.end(), lnCollisions );
514 if ( it != collisions.end() ) {
516 std::vector<float>::difference_type ihi = it - collisions.begin();
517 LogTrace(
me) <<
": using one energy bin " << ihi <<
" = " 519 <<
" for lnCollisions = " << lnCollisions;
524 std::vector<float>::const_iterator loside = find_if(collisions.begin(),
525 collisions.end(), [&lnCollisions](
auto c){
return c < lnCollisions;});
526 std::vector<float>::difference_type ilo = loside - collisions.begin();
529 << ilo-1 <<
" and " << ilo;
530 lnE =
theEnergyBins[ilo-1] + (lnCollisions-collisions[ilo-1])*
532 (collisions[ilo]-collisions[ilo-1]);
537 <<
" for lnCollisions = " << lnCollisions;
546 std::vector<float>& collisions )
const 548 std::vector<float>::const_iterator bigger = find_if(
theGammaBins.begin(),
554 <<
" for logGamma = " << logGamma;
560 std::vector<float>::difference_type ihi = bigger -
theGammaBins.begin();
562 double dlg2 = *bigger--;
565 double dlg1 = *bigger;
566 double dlg = (logGamma-dlg1)/(dlg2-dlg1);
567 double omdlg = 1. - dlg;
575 <<
" 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
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
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)