26 #include "CLHEP/Random/RandExponential.h"
27 #include "CLHEP/Random/RandFlat.h"
59 :
me(
"CSCGasCollisions"),
60 gasDensity(2.1416
e-03),
65 theGammaBins(N_GAMMA, 0.),
66 theEnergyBins(N_ENERGY, 0.),
67 theCollisionTable(N_ENTRIES, 0.),
69 theParticleDataTable(nullptr),
70 saveGasCollisions_(
false),
71 dumpGasCollisions_(
false) {
77 <<
" keV (for higher elosses, hits should have been simulated.)";
107 string colliFile =
"SimMuon/CSCDigitizer/data/collisions.dat";
111 ifstream
fin{
f1.fullPath()};
114 string errorMessage =
"Cannot open input file " +
f1.fullPath();
120 fin.seekg(0, ios::beg);
136 LogTrace(
me) <<
"Reading collisions table \n";
149 std::vector<LocalPoint> &positions,
151 CLHEP::HepRandomEngine *engine) {
152 const float epsilonL = 0.01;
160 double mom = simhit.
pabs();
165 double mass = 0.105658;
166 if (particle ==
nullptr) {
169 mass = particle->mass();
171 <<
" in the PDT; mass = " <<
mass;
179 if (gapSize <= epsilonL) {
180 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! simhit entry and exit are too close - "
184 <<
" : 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;
200 LogTrace(
me) <<
"collisions extremes = " << collisions[
N_ENERGY - 1] <<
", " << collisions[0] <<
"\n"
201 <<
"anmin = " << anmin <<
", anmax = " << anmax <<
"\n"
202 <<
"amu = " << amu <<
"\n";
205 double sum_steps = 0.;
213 while (sum_steps < gapSize) {
217 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! n_try = " << n_try
218 <<
" is > MAX_STEPS = " << maxst <<
" - skipping simhit:"
221 <<
" : momentum = " << simhit.
pabs()
222 <<
" GeV/c : energy loss = " << simhit.
energyLoss() * 1.E06 <<
" keV"
223 <<
"\n gapSize = " << gapSize <<
" cm, last step = " << step
224 <<
" cm, sum_steps = " << sum_steps <<
" cm, n_steps = " << n_steps;
228 if (sum_steps + step > gapSize)
235 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! eloss = " << eloss <<
" eV is too large (> "
236 <<
deCut <<
" eV) - trying another collision";
245 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! " << n_steps
246 <<
" is too many steps - skipping simhit:"
249 <<
" : momentum = " << simhit.
pabs()
250 <<
" GeV/c : energy loss = " << simhit.
energyLoss() * 1.E06 <<
" keV"
251 <<
"\ngapSize=" << gapSize <<
" cm, last step=" << step
252 <<
" cm, sum_steps=" << sum_steps <<
" cm";
255 LogTrace(
me) <<
"sum_steps = " << sum_steps <<
" cm , dedx = " << dedx <<
" eV";
261 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;
296 const std::vector<float> &collisions,
297 CLHEP::HepRandomEngine *engine)
const {
299 float lnColl =
log(CLHEP::RandFlat::shoot(engine, anmin, anmax));
308 float eloss =
exp(lnE);
320 while (energyAvailable >
eion) {
339 LogTrace(
me) <<
" range = " << range <<
" cm";
344 int nelec =
static_cast<int>(energyAvailable /
ework);
345 LogTrace(
me) <<
"short-range delta energy in = " << energyAvailable <<
" eV";
347 energyAvailable -= nelec *
ework;
349 if (energyAvailable >
eion) {
351 energyAvailable -=
eion;
353 LogTrace(
me) <<
"short-range delta energy out = " << energyAvailable <<
" eV, nelec = " << nelec;
363 bool new_range =
false;
364 while (!new_range && (energyAvailable >
ework)) {
365 energyAvailable -=
ework;
366 while (energyAvailable >
eion) {
367 double range2 = 0.5 * 0.71 /
gasDensity *
pow(1.E-6 * energyAvailable, 1.72);
368 double drange = range - range2;
369 LogTrace(
me) <<
" energy left = " << energyAvailable <<
" eV, range2 = " << range2
370 <<
" cm, drange = " << drange <<
" cm";
377 LogTrace(
me) <<
"reset range to range2 = " << range <<
" from startHere = " << startHere <<
" and iterate";
387 if (energyAvailable >
eion) {
388 energyAvailable -=
ework;
406 <<
"\nAFTER CROSSING GAP";
425 size_t nsteps = n_steps;
426 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;
438 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] # / length of step / energy loss per collision:";
439 for (
size_t i = 0;
i != nsteps; ++
i) {
440 edm::LogVerbatim(
"CSCDigitizer") <<
i + 1 <<
" / S: " << steps[
i] <<
" / E: " << elosses[
i];
444 size_t mclus = ion_clusters.size();
445 size_t melec = electrons.size();
447 if (mclus != melec) {
448 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] WARNING! size of cluster vector = " << mclus
449 <<
" .ne. size of electrons vector = " << melec;
452 "cluster / electrons per cluster: ";
453 for (
size_t i = 0;
i != mclus; ++
i) {
454 edm::LogVerbatim(
"CSCDigitizer") <<
i + 1 <<
" / I: " << ion_clusters[
i] <<
" / E: " << electrons[
i];
458 int n_ic = count_if(elosses.begin(), elosses.end(), [&](
auto c) {
return c > this->
eion; });
460 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] total no. of collision steps = " << n_steps;
461 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] total sum of steps = " << sum_steps <<
" cm";
463 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] average step length = " << sum_steps / float(nsteps)
465 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] total energy loss across gap = " << dedx
466 <<
" eV = " << dedx / 1000. <<
" keV";
467 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] no. of primary ionizing collisions across gap = "
468 "no. with eloss > eion = "
469 <<
eion <<
" eV = " << n_ic;
471 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] average energy loss/collision = " << dedx / float(nsteps)
474 std::vector<int>::const_iterator bigger =
find(electrons.begin(), electrons.end(), 0);
475 if (bigger != electrons.end()) {
476 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] TROUBLE! There is a cluster with 0 electrons.";
479 int n_e = accumulate(electrons.begin(), electrons.end(), 0);
484 <<
" : momentum = " << simhit.
pabs()
485 <<
" GeV/c : energy loss = " << simhit.
energyLoss() * 1.E06 <<
" keV";
488 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisions] SUMMARY: ionization"
489 <<
" : steps= " << nsteps <<
" : sum(steps)= " << sum_steps
490 <<
" cm : <step>= " << sum_steps / float(nsteps) <<
" cm"
491 <<
" : ionizing= " << n_ic <<
" : ionclus= " << mclus <<
" : total e= " << n_e
492 <<
" : <dedx/step>= " << dedx / float(nsteps)
493 <<
" eV : <e/ionclus>= " << float(n_e) / float(mclus)
494 <<
" : dedx= " << dedx / 1000. <<
" keV";
496 edm::LogVerbatim(
"CSCDigitizer") <<
"[CSCGasCollisons] ERROR? no collision steps!";
523 std::vector<float>::const_iterator it =
find(collisions.begin(), collisions.end(), lnCollisions);
525 if (it != collisions.end()) {
527 std::vector<float>::difference_type ihi = it - collisions.begin();
529 <<
" for lnCollisions = " << lnCollisions;
533 std::vector<float>::const_iterator loside =
534 find_if(collisions.begin(), collisions.end(), [&lnCollisions](
auto c) {
return c < lnCollisions; });
535 std::vector<float>::difference_type ilo = loside - collisions.begin();
537 LogTrace(
me) <<
": using energy bin " << ilo - 1 <<
" and " << ilo;
538 lnE =
theEnergyBins[ilo - 1] + (lnCollisions - collisions[ilo - 1]) *
540 (collisions[ilo] - collisions[ilo - 1]);
551 std::vector<float>::const_iterator bigger =
557 <<
" for logGamma = " << logGamma;
562 std::vector<float>::difference_type ihi = bigger -
theGammaBins.begin();
564 double dlg2 = *bigger--;
567 double dlg1 = *bigger;
568 double dlg = (logGamma - dlg1) / (dlg2 - dlg1);
569 double omdlg = 1. - dlg;
575 <<
" for logGamma = " << logGamma;
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
static std::vector< std::string > checklist log
const edm::EventSetup & c
std::vector< int > electrons() const
HepPDT::ParticleDataTable ParticleDataTable
bool saveGasCollisions(void) const
void addEloss(float eloss)
static const int N_ENERGY
Exp< T >::type exp(const T &t)
std::vector< LocalPoint > ionClusters() const
Log< level::Error, false > LogError
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
const uint16_t range(const Frame &aFrame)
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
Log< level::Info, false > LogInfo
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)