CMS 3D CMS Logo

MaterialEffectsSimulator.h
Go to the documentation of this file.
1 #ifndef MATERIALEFFECTSSIMULATOR_H
2 #define MATERIALEFFECTSSIMULATOR_H
3 
5 
7 
8 #include <vector>
9 
11 
26 public:
27  typedef std::vector<RawParticle>::const_iterator RHEP_const_iter;
28 
29  // Constructor : default values are for Silicon
30  MaterialEffectsSimulator(double A = 28.0855, double Z = 14.0000, double density = 2.329, double radLen = 9.360);
31 
32  virtual ~MaterialEffectsSimulator();
33 
36 
38  inline double theA() const { return A; }
40  inline double theZ() const { return Z; }
42  inline double rho() const { return density; }
44  inline double radLenIncm() const { return radLen; }
46  inline double excitE() const { return 12.5E-9 * theZ(); }
48  inline double eMass() const { return 0.000510998902; }
49 
51  void updateState(ParticlePropagator& myTrack, double radlen, RandomEngineAndDistribution const*);
52 
54  inline RHEP_const_iter beginDaughters() const { return _theUpdatedState.begin(); }
55 
57  inline RHEP_const_iter endDaughters() const { return _theUpdatedState.end(); }
58 
60  inline unsigned nDaughters() const { return _theUpdatedState.size(); }
61 
63  inline void setNormalVector(const GlobalVector& normal) { theNormalVector = normal; }
64 
66  XYZVector orthogonal(const XYZVector&) const;
67 
70 
72  virtual void save(){};
73 
74 private:
77 
79  inline double radiationLength() const { return radLengths; }
80 
81 protected:
82  std::vector<RawParticle> _theUpdatedState;
83 
84  double radLengths;
85 
86  // Material properties
87  double A;
88  double Z;
89  double density;
90  double radLen;
91 
93 
95 };
96 
97 #endif
double radLenIncm() const
One radiation length in cm.
virtual void compute(ParticlePropagator &Particle, RandomEngineAndDistribution const *)=0
Overloaded in all material effects updtators.
RHEP_const_iter endDaughters() const
Returns const iterator to the end of the daughters list.
double radiationLength() const
Returns the fraction of radiation lengths traversed.
int closestDaughterId()
The id of the closest charged daughter (filled for nuclear interactions only)
void setNormalVector(const GlobalVector &normal)
Sets the vector normal to the surface traversed.
double rho() const
Density in g/cm3.
void updateState(ParticlePropagator &myTrack, double radlen, RandomEngineAndDistribution const *)
Compute the material effect (calls the sub class)
std::vector< RawParticle >::const_iterator RHEP_const_iter
double eMass() const
Electron mass in GeV/c2.
MaterialEffectsSimulator(double A=28.0855, double Z=14.0000, double density=2.329, double radLen=9.360)
XYZVector orthogonal(const XYZVector &) const
A vector orthogonal to another one (because it&#39;s not in XYZTLorentzVector)
double excitE() const
Mean excitation energy (in GeV)
RHEP_const_iter beginDaughters() const
Returns const iterator to the beginning of the daughters list.
std::vector< RawParticle > _theUpdatedState
math::XYZVector XYZVector
Definition: RawParticle.h:26
unsigned nDaughters() const
Returns the number of daughters.
virtual void save()
Used by NuclearInteractionSimulator to save last sampled event.