CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCGasCollisions.h
Go to the documentation of this file.
1 #ifndef MU_END_CSC_GAS_COLLISIONS_H
2 #define MU_END_CSC_GAS_COLLISIONS_H
3 
18 #include "CLHEP/Random/RandFlat.h"
19 #include "CLHEP/Random/RandExponential.h"
20 
21 #include <vector>
22 #include <string>
23 
25 public:
26 
28  virtual ~CSCGasCollisions();
29 
30  void setParticleDataTable(const ParticleDataTable * pdt);
31 
32  void setRandomEngine(CLHEP::HepRandomEngine & engine);
33 
34  void simulate(const PSimHit&,
35  std::vector<LocalPoint>& clusters, std::vector<int>& electrons );
36 
37  static const int N_GAMMA = 21;
38  static const int N_ENERGY = 63;
39  static const int N_ENTRIES = N_GAMMA*N_ENERGY;
40  static const int MAX_STEPS = 400;
41 
42 private:
43  void readCollisionTable();
44  void fillCollisionsForThisGamma( float, std::vector<float>& ) const;
45  float lnEnergyLoss( float, const std::vector<float>& ) const;
46  double generateStep( double avCollisions ) const;
47  float generateEnergyLoss( double avCollisions,
48  double anmin, double anmax, const std::vector<float>& collisions ) const;
49 
50  void ionize( double energyTransferred, LocalPoint startHere) const;
51 
52  void writeSummary( int n_steps, double sum_steps, float dedx, float simHiteloss ) const;
53 
54  const std::string me; // class name
55  double gasDensity; // Density of CSC gas mix
56  // The question of what is reasonable for deCut is complex. But it seems clear that
57  // this simulation is not credible if any delta electrons generated here
58  // have ranges more than a few mm, or equivalently, energies above a few keV.
59  // deCut = 1000 = 1 keV
60  // deCut = 10000 = 10 keV
61  double deCut; // Delta electron cutoff in eV (Match GEANT!)
62  double eion; // ionization threshold (eV) (min. E for ionizatiom)
63  double ework; // effective work function (av. energy to create one ion pair)
64  double clusterExtent; // Precision of localization of ion clus. Typically 10 microns.
65 
66  std::vector<float> theGammaBins;
67  std::vector<float> theEnergyBins;
68  std::vector<float> theCollisionTable;
69 
70  CSCCrossGap* theCrossGap; // Owned by CSCGasCollisions
72  CLHEP::RandFlat * theRandFlat;
73  CLHEP::RandExponential * theRandExponential;
74  bool saveGasCollisions; // Simple Configurable to flag saving info w. debugV
75 };
76 
77 #endif
void simulate(const PSimHit &, std::vector< LocalPoint > &clusters, std::vector< int > &electrons)
HepPDT::ParticleDataTable ParticleDataTable
static const int N_ENERGY
void writeSummary(int n_steps, double sum_steps, float dedx, float simHiteloss) const
std::vector< float > theGammaBins
const std::string me
float lnEnergyLoss(float, const std::vector< float > &) const
const ParticleDataTable * theParticleDataTable
static const int N_GAMMA
static const int MAX_STEPS
void setParticleDataTable(const ParticleDataTable *pdt)
void ionize(double energyTransferred, LocalPoint startHere) const
static const int N_ENTRIES
CLHEP::RandExponential * theRandExponential
std::vector< float > theCollisionTable
double generateStep(double avCollisions) const
std::vector< float > theEnergyBins
float generateEnergyLoss(double avCollisions, double anmin, double anmax, const std::vector< float > &collisions) const
void fillCollisionsForThisGamma(float, std::vector< float > &) const
CSCCrossGap * theCrossGap
virtual ~CSCGasCollisions()
CLHEP::RandFlat * theRandFlat
void setRandomEngine(CLHEP::HepRandomEngine &engine)