CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 {
27  public:
28 
29  typedef std::vector<RawParticle>::const_iterator RHEP_const_iter;
30 
31  // Constructor : default values are for Silicon
32  MaterialEffectsSimulator(double A = 28.0855,
33  double Z = 14.0000,
34  double density = 2.329,
35  double radLen = 9.360);
36 
37  virtual ~MaterialEffectsSimulator();
38 
41 
43  inline double theA() const { return A; }
45  inline double theZ() const { return Z; }
47  inline double rho() const { return density; }
49  inline double radLenIncm() const { return radLen; }
51  inline double excitE() const { return 12.5E-9*theZ(); }
53  inline double eMass() const { return 0.000510998902; }
54 
55 
57  void updateState(ParticlePropagator& myTrack, double radlen, RandomEngineAndDistribution const*);
58 
60  inline RHEP_const_iter beginDaughters() const {return _theUpdatedState.begin();}
61 
63  inline RHEP_const_iter endDaughters() const {return _theUpdatedState.end();}
64 
66  inline unsigned nDaughters() const {return _theUpdatedState.size();}
67 
69  inline void setNormalVector(const GlobalVector& normal) { theNormalVector = normal; }
70 
72  XYZVector orthogonal(const XYZVector&) const;
73 
76 
77  private:
78 
80  virtual void compute(ParticlePropagator& Particle, RandomEngineAndDistribution const*) = 0;
81 
83  inline double radiationLength() const {return radLengths;}
84 
85 
86  protected:
87 
88  std::vector<RawParticle> _theUpdatedState;
89 
90  double radLengths;
91 
92  // Material properties
93  double A;
94  double Z;
95  double density;
96  double radLen;
97 
99 
101 
102 };
103 
104 #endif
virtual void compute(ParticlePropagator &Particle, RandomEngineAndDistribution const *)=0
Overloaded in all material effects updtators.
RHEP_const_iter beginDaughters() const
Returns const iterator to the beginning of the daughters list.
unsigned nDaughters() const
Returns the number of daughters.
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.
XYZVector orthogonal(const XYZVector &) const
A vector orthogonal to another one (because it&#39;s not in XYZTLorentzVector)
double rho() const
Density in g/cm3.
void updateState(ParticlePropagator &myTrack, double radlen, RandomEngineAndDistribution const *)
Compute the material effect (calls the sub class)
double radiationLength() const
Returns the fraction of radiation lengths traversed.
math::XYZVector XYZVector
std::vector< RawParticle >::const_iterator RHEP_const_iter
MaterialEffectsSimulator(double A=28.0855, double Z=14.0000, double density=2.329, double radLen=9.360)
RHEP_const_iter endDaughters() const
Returns const iterator to the end of the daughters list.
double eMass() const
Electron mass in GeV/c2.
double excitE() const
Mean excitation energy (in GeV)
std::vector< RawParticle > _theUpdatedState
double radLenIncm() const
One radiation length in cm.