CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions
VolumeEnergyLossEstimator Class Reference

#include <VolumeEnergyLossEstimator.h>

Inheritance diagram for VolumeEnergyLossEstimator:
VolumeMaterialEffectsEstimator

Public Member Functions

virtual VolumeEnergyLossEstimatorclone () const
 
virtual
VolumeMaterialEffectsEstimate 
estimate (const TrajectoryStateOnSurface refTSOS, double pathLength, const VolumeMediumProperties &medium) const
 Creates an estimate. More...
 
 VolumeEnergyLossEstimator (float mass)
 Constructor with explicit mass hypothesis. More...
 
virtual ~VolumeEnergyLossEstimator ()
 
- Public Member Functions inherited from VolumeMaterialEffectsEstimator
virtual float mass () const
 Particle mass assigned at construction. More...
 
 VolumeMaterialEffectsEstimator (float mass)
 Constructor with explicit mass hypothesis. More...
 
virtual ~VolumeMaterialEffectsEstimator ()
 

Private Member Functions

void computeBetheBloch (const LocalVector &localP, double pathLength, const VolumeMediumProperties &medium, double &deltaP, double &deltaCov00) const
 
void computeElectrons (const LocalVector &localP, double pathLength, const VolumeMediumProperties &medium, double &deltaP, double &deltaCov00) const
 

Detailed Description

Estimation of energy loss for a finite step size in a volume.

Definition at line 16 of file VolumeEnergyLossEstimator.h.

Constructor & Destructor Documentation

VolumeEnergyLossEstimator::VolumeEnergyLossEstimator ( float  mass)
inline

Constructor with explicit mass hypothesis.

Definition at line 20 of file VolumeEnergyLossEstimator.h.

Referenced by clone().

20  :
virtual float mass() const
Particle mass assigned at construction.
VolumeMaterialEffectsEstimator(float mass)
Constructor with explicit mass hypothesis.
virtual VolumeEnergyLossEstimator::~VolumeEnergyLossEstimator ( )
inlinevirtual

Definition at line 23 of file VolumeEnergyLossEstimator.h.

23 {}

Member Function Documentation

VolumeEnergyLossEstimator * VolumeEnergyLossEstimator::clone ( void  ) const
virtual

Implements VolumeMaterialEffectsEstimator.

Definition at line 40 of file VolumeEnergyLossEstimator.cc.

References VolumeEnergyLossEstimator().

41 {
42  return new VolumeEnergyLossEstimator(*this);
43 }
VolumeEnergyLossEstimator(float mass)
Constructor with explicit mass hypothesis.
void VolumeEnergyLossEstimator::computeBetheBloch ( const LocalVector localP,
double  pathLength,
const VolumeMediumProperties medium,
double &  deltaP,
double &  deltaCov00 
) const
private

Definition at line 49 of file VolumeEnergyLossEstimator.cc.

References beta, delta, alignCSCRings::e, create_public_lumi_plots::log, contentValuesFiles::m, PV3DBase< T, PVType, FrameType >::mag(), VolumeMaterialEffectsEstimator::mass(), AlCaHLTBitMon_ParallelJobs::p, mathSSE::sqrt(), and VolumeMediumProperties::xi().

Referenced by estimate().

51  {
52  //
53  // calculate absolute momentum and correction to path length from angle
54  // of incidence
55  //
56  double p = localP.mag();
57  // constants
58  const double m = mass(); // use mass hypothesis from constructor
59  const double emass = 0.511e-3;
60  // FIXME: replace constants for Si
61  const double poti = 16.e-9 * 10.75; // = 16 eV * Z**0.9, for Si Z=14
62  const double eplasma = 28.816e-9 * sqrt(2.33*0.498); // 28.816 eV * sqrt(rho*(Z/A)) for Si
63  const double delta0 = 2*log(eplasma/poti) - 1.;
64  // calculate general physics things
65  double e = sqrt(p*p + m*m);
66  double beta = p/e;
67  double gamma = e/m;
68  double eta2 = beta*gamma; eta2 *= eta2;
69  double lnEta2 = log(eta2);
70  double ratio = emass/m;
71  double emax = 2.*emass*eta2/(1. + 2.*ratio*gamma + ratio*ratio);
72  double delta = delta0 + lnEta2;
73  // calculate the mean and sigma of energy loss
74  // xi = d[g/cm2] * 0.307075MeV/(g/cm2) * Z/A * 1/2
75  double xi = pathLength * medium.xi(); xi /= (beta*beta);
76 // double dEdx = xi*(log(2.*emass*eta2*emax/(poti*poti)) - 2.*(beta*beta));
77  double dEdx = xi*(log(2.*emass*emax/(poti*poti))+lnEta2 - 2.*(beta*beta) - delta);
78  double dEdx2 = xi*emax*(1.-beta*beta/2.);
79  double dP = dEdx/beta;
80  double sigp2 = dEdx2*e*e/(p*p*p*p*p*p);
81  deltaP = -dP;
82  deltaCov00 = sigp2;
83 }
const double beta
dbl * delta
Definition: mlp_gen.cc:36
virtual float mass() const
Particle mass assigned at construction.
T mag() const
Definition: PV3DBase.h:67
T sqrt(T t)
Definition: SSEVec.h:48
void VolumeEnergyLossEstimator::computeElectrons ( const LocalVector localP,
double  pathLength,
const VolumeMediumProperties medium,
double &  deltaP,
double &  deltaCov00 
) const
private

Definition at line 88 of file VolumeEnergyLossEstimator.cc.

References create_public_lumi_plots::exp, f, create_public_lumi_plots::log, PV3DBase< T, PVType, FrameType >::mag(), AlCaHLTBitMon_ParallelJobs::p, VolumeMediumProperties::x0(), and detailsBasic3DVector::z.

Referenced by estimate().

90  {
91  //
92  // calculate absolute momentum and correction to path length from angle
93  // of incidence
94  //
95  double p = localP.mag();
96  double normalisedPath = pathLength / medium.x0();
97  //
98  // Energy loss and variance according to Bethe and Heitler, see also
99  // Comp. Phys. Comm. 79 (1994) 157.
100  //
101  double z = exp(-normalisedPath);
102  double varz = (exp(-normalisedPath*log(3.)/log(2.))-
103  exp(-2*normalisedPath));
104  // FIXME: need to know propagation direction at this point
105 // if ( propDir==oppositeToMomentum ) {
106 // //
107 // // for backward propagation: delta(1/p) is linear in z=p_outside/p_inside
108 // // convert to obtain equivalent delta(p). Sign of deltaP is corrected
109 // // in method compute -> deltaP<0 at this place!!!
110 // //
111 // deltaP = -p*(1/z-1);
112 // deltaCov00 = varz/p/p;
113 // }
114 // else {
115  //
116  // for forward propagation: calculate in p (linear in 1/z=p_inside/p_outside),
117  // then convert sig(p) to sig(1/p).
118  //
119  deltaP = p*(z-1);
120  // double f = 1/p/z/z;
121  // patch to ensure consistency between for- and backward propagation
122  double f = 1./p/z;
123  deltaCov00 += f*f*varz;
124 }
float float float z
T mag() const
Definition: PV3DBase.h:67
double f[11][100]
VolumeMaterialEffectsEstimate VolumeEnergyLossEstimator::estimate ( const TrajectoryStateOnSurface  refTSOS,
double  pathLength,
const VolumeMediumProperties medium 
) const
virtual

Creates an estimate.

Implements VolumeMaterialEffectsEstimator.

Definition at line 10 of file VolumeEnergyLossEstimator.cc.

References computeBetheBloch(), computeElectrons(), TrajectoryStateOnSurface::localMomentum(), and VolumeMaterialEffectsEstimator::mass().

13 {
14  //
15  // Initialise dP and the update to the covariance matrix
16  //
17  double deltaP = 0.;
18  double deltaCov00 = 0.;
19  //
20  // Bethe-Bloch
21  //
22  if ( mass()>0.001 )
23  computeBetheBloch(refTSOS.localMomentum(),pathLength,medium,
24  deltaP,deltaCov00);
25  //
26  // Special treatment for electrons (currently rather crude
27  // distinction using mass)
28  //
29  else
30  computeElectrons(refTSOS.localMomentum(),pathLength,medium,
31  deltaP,deltaCov00);
32 
33  AlgebraicSymMatrix55 deltaCov;
34  deltaCov(0,0) = deltaCov00;
35  return VolumeMaterialEffectsEstimate(deltaP,deltaCov);
36 }
virtual float mass() const
Particle mass assigned at construction.
void computeBetheBloch(const LocalVector &localP, double pathLength, const VolumeMediumProperties &medium, double &deltaP, double &deltaCov00) const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
LocalVector localMomentum() const
void computeElectrons(const LocalVector &localP, double pathLength, const VolumeMediumProperties &medium, double &deltaP, double &deltaCov00) const