CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Protected Attributes | Private Member Functions
MaterialEffectsSimulator Class Referenceabstract

#include <MaterialEffectsSimulator.h>

Inheritance diagram for MaterialEffectsSimulator:
BremsstrahlungSimulator EnergyLossSimulator MultipleScatteringSimulator MuonBremsstrahlungSimulator NuclearInteractionFTFSimulator NuclearInteractionSimulator PairProductionSimulator

Public Types

typedef std::vector
< RawParticle >
::const_iterator 
RHEP_const_iter
 

Public Member Functions

RHEP_const_iter beginDaughters () const
 Returns const iterator to the beginning of the daughters list. More...
 
int closestDaughterId ()
 The id of the closest charged daughter (filled for nuclear interactions only) More...
 
double eMass () const
 Electron mass in GeV/c2. More...
 
RHEP_const_iter endDaughters () const
 Returns const iterator to the end of the daughters list. More...
 
double excitE () const
 Mean excitation energy (in GeV) More...
 
 MaterialEffectsSimulator (double A=28.0855, double Z=14.0000, double density=2.329, double radLen=9.360)
 
unsigned nDaughters () const
 Returns the number of daughters. More...
 
XYZVector orthogonal (const XYZVector &) const
 A vector orthogonal to another one (because it's not in XYZTLorentzVector) More...
 
double radLenIncm () const
 One radiation length in cm. More...
 
double rho () const
 Density in g/cm3. More...
 
virtual void save ()
 Used by NuclearInteractionSimulator to save last sampled event. More...
 
void setNormalVector (const GlobalVector &normal)
 Sets the vector normal to the surface traversed. More...
 
double theA () const
 A. More...
 
double theZ () const
 Z. More...
 
void updateState (ParticlePropagator &myTrack, double radlen, RandomEngineAndDistribution const *)
 Compute the material effect (calls the sub class) More...
 
virtual ~MaterialEffectsSimulator ()
 

Protected Attributes

std::vector< RawParticle_theUpdatedState
 
double A
 
double density
 
double radLen
 
double radLengths
 
int theClosestChargedDaughterId
 
GlobalVector theNormalVector
 
double Z
 

Private Member Functions

virtual void compute (ParticlePropagator &Particle, RandomEngineAndDistribution const *)=0
 Overloaded in all material effects updtators. More...
 
double radiationLength () const
 Returns the fraction of radiation lengths traversed. More...
 

Detailed Description

This is the generic class for Material Effects in the tracker material, from which FamosPairProductionSimulator, FamosBremsstrahlungSimulator, FamosEnergyLossSimulator and FamosMultipleScatteringSimulator inherit. It determines the fraction of radiation lengths traversed by the particle in this tracker layer, defines the tracker layer characteristics (hmmm, 100% Silicon, very crude, but what can we do?) and returns a list of new RawParticles if needed.

Author
Stephan Wynhoff, Florian Beaudette, Patrick Janot $Date: Last update 8-Jan-2004

Definition at line 25 of file MaterialEffectsSimulator.h.

Member Typedef Documentation

typedef std::vector<RawParticle>::const_iterator MaterialEffectsSimulator::RHEP_const_iter

Definition at line 29 of file MaterialEffectsSimulator.h.

Constructor & Destructor Documentation

MaterialEffectsSimulator::MaterialEffectsSimulator ( double  A = 28.0855,
double  Z = 14.0000,
double  density = 2.329,
double  radLen = 9.360 
)
MaterialEffectsSimulator::~MaterialEffectsSimulator ( )
virtual

Definition at line 16 of file MaterialEffectsSimulator.cc.

References _theUpdatedState.

16  {
17  // Don't delete the objects contained in the list
18  _theUpdatedState.clear();
19 }
std::vector< RawParticle > _theUpdatedState

Member Function Documentation

RHEP_const_iter MaterialEffectsSimulator::beginDaughters ( ) const
inline

Returns const iterator to the beginning of the daughters list.

Definition at line 59 of file MaterialEffectsSimulator.h.

References _theUpdatedState.

Referenced by MaterialEffects::interact().

59 {return _theUpdatedState.begin();}
std::vector< RawParticle > _theUpdatedState
int MaterialEffectsSimulator::closestDaughterId ( )
inline

The id of the closest charged daughter (filled for nuclear interactions only)

Definition at line 74 of file MaterialEffectsSimulator.h.

References theClosestChargedDaughterId.

Referenced by MaterialEffects::interact().

virtual void MaterialEffectsSimulator::compute ( ParticlePropagator Particle,
RandomEngineAndDistribution const *   
)
privatepure virtual
double MaterialEffectsSimulator::eMass ( ) const
inline

Electron mass in GeV/c2.

Definition at line 53 of file MaterialEffectsSimulator.h.

Referenced by PairProductionSimulator::compute(), and EnergyLossSimulator::compute().

53 { return 0.000510998902; }
RHEP_const_iter MaterialEffectsSimulator::endDaughters ( ) const
inline

Returns const iterator to the end of the daughters list.

Definition at line 62 of file MaterialEffectsSimulator.h.

References _theUpdatedState.

Referenced by MaterialEffects::interact().

62 {return _theUpdatedState.end();}
std::vector< RawParticle > _theUpdatedState
double MaterialEffectsSimulator::excitE ( ) const
inline

Mean excitation energy (in GeV)

Definition at line 51 of file MaterialEffectsSimulator.h.

References theZ().

Referenced by EnergyLossSimulator::compute().

51 { return 12.5E-9*theZ(); }
unsigned MaterialEffectsSimulator::nDaughters ( ) const
inline

Returns the number of daughters.

Definition at line 65 of file MaterialEffectsSimulator.h.

References _theUpdatedState.

Referenced by MaterialEffects::interact().

65 {return _theUpdatedState.size();}
std::vector< RawParticle > _theUpdatedState
XYZVector MaterialEffectsSimulator::orthogonal ( const XYZVector aVector) const

A vector orthogonal to another one (because it's not in XYZTLorentzVector)

Definition at line 35 of file MaterialEffectsSimulator.cc.

References x, y, and z.

Referenced by MultipleScatteringSimulator::compute(), and NuclearInteractionSimulator::compute().

35  {
36 
37  double x = fabs(aVector.X());
38  double y = fabs(aVector.Y());
39  double z = fabs(aVector.Z());
40 
41  if ( x < y )
42  return x < z ?
43  XYZVector(0.,aVector.Z(),-aVector.Y()) :
44  XYZVector(aVector.Y(),-aVector.X(),0.);
45  else
46  return y < z ?
47  XYZVector(-aVector.Z(),0.,aVector.X()) :
48  XYZVector(aVector.Y(),-aVector.X(),0.);
49 
50 }
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
double MaterialEffectsSimulator::radiationLength ( ) const
inlineprivate

Returns the fraction of radiation lengths traversed.

Definition at line 85 of file MaterialEffectsSimulator.h.

References radLengths.

double MaterialEffectsSimulator::radLenIncm ( ) const
inline

One radiation length in cm.

Definition at line 49 of file MaterialEffectsSimulator.h.

References radLen.

Referenced by MultipleScatteringSimulator::compute(), and EnergyLossSimulator::compute().

double MaterialEffectsSimulator::rho ( ) const
inline

Density in g/cm3.

Definition at line 47 of file MaterialEffectsSimulator.h.

References density.

Referenced by Lepton.Lepton::absIsoFromEA(), Muon.Muon::absIsoWithFSR(), and EnergyLossSimulator::compute().

virtual void MaterialEffectsSimulator::save ( )
inlinevirtual
void MaterialEffectsSimulator::setNormalVector ( const GlobalVector normal)
inline

Sets the vector normal to the surface traversed.

Definition at line 68 of file MaterialEffectsSimulator.h.

References theNormalVector.

Referenced by MuonSimHitProducer::applyMaterialEffects(), and MaterialEffects::interact().

68 { theNormalVector = normal; }
double MaterialEffectsSimulator::theA ( ) const
inline

A.

Functions to return atomic properties of the material Here the tracker material is assumed to be 100% Silicon

Definition at line 43 of file MaterialEffectsSimulator.h.

References A.

Referenced by EnergyLossSimulator::compute(), and NuclearInteractionSimulator::compute().

43 { return A; }
double MaterialEffectsSimulator::theZ ( ) const
inline
void MaterialEffectsSimulator::updateState ( ParticlePropagator myTrack,
double  radlen,
RandomEngineAndDistribution const *  random 
)

Compute the material effect (calls the sub class)

Definition at line 21 of file MaterialEffectsSimulator.cc.

References _theUpdatedState, compute(), radLengths, and theClosestChargedDaughterId.

Referenced by MuonSimHitProducer::applyMaterialEffects(), MaterialEffects::interact(), and CalorimetryManager::MuonMipSimulation().

24 {
25 
26  _theUpdatedState.clear();
28 
29  radLengths = radlen;
30  if ( radLengths > 0. ) compute(Particle, random);
31 
32 }
virtual void compute(ParticlePropagator &Particle, RandomEngineAndDistribution const *)=0
Overloaded in all material effects updtators.
TRandom random
Definition: MVATrainer.cc:138
std::vector< RawParticle > _theUpdatedState

Member Data Documentation

std::vector<RawParticle> MaterialEffectsSimulator::_theUpdatedState
protected
double MaterialEffectsSimulator::A
protected

Definition at line 94 of file MaterialEffectsSimulator.h.

Referenced by MuonBremsstrahlungSimulator::compute(), and theA().

double MaterialEffectsSimulator::density
protected

Definition at line 96 of file MaterialEffectsSimulator.h.

Referenced by MuonBremsstrahlungSimulator::compute(), and rho().

double MaterialEffectsSimulator::radLen
protected

Definition at line 97 of file MaterialEffectsSimulator.h.

Referenced by MuonBremsstrahlungSimulator::compute(), and radLenIncm().

double MaterialEffectsSimulator::radLengths
protected
int MaterialEffectsSimulator::theClosestChargedDaughterId
protected
GlobalVector MaterialEffectsSimulator::theNormalVector
protected
double MaterialEffectsSimulator::Z
protected