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
EnergyLossUpdator Class Reference

#include <EnergyLossUpdator.h>

Inheritance diagram for EnergyLossUpdator:
MaterialEffectsUpdator

Public Member Functions

virtual EnergyLossUpdatorclone () const
 
 EnergyLossUpdator (double mass)
 
- Public Member Functions inherited from MaterialEffectsUpdator
virtual const
AlgebraicSymMatrix55
deltaLocalError (const TrajectoryStateOnSurface &TSoS, const PropagationDirection propDir) const
 
virtual double deltaP (const TrajectoryStateOnSurface &TSoS, const PropagationDirection propDir) const
 
double mass () const
 
 MaterialEffectsUpdator (double mass)
 
virtual TrajectoryStateOnSurface updateState (const TrajectoryStateOnSurface &TSoS, const PropagationDirection propDir) const
 
virtual bool updateStateInPlace (TrajectoryStateOnSurface &TSoS, const PropagationDirection propDir) const
 
virtual ~MaterialEffectsUpdator ()
 

Private Member Functions

virtual void compute (const TrajectoryStateOnSurface &, const PropagationDirection) const
 
void computeBetheBloch (const LocalVector &, const MediumProperties &) const
 
void computeElectrons (const LocalVector &, const MediumProperties &, const PropagationDirection) const
 

Additional Inherited Members

- Protected Attributes inherited from MaterialEffectsUpdator
AlgebraicSymMatrix55 theDeltaCov
 
double theDeltaP
 
- Static Protected Attributes inherited from MaterialEffectsUpdator
static AlgebraicSymMatrix55 theNullMatrix
 

Detailed Description

Energy loss according to Bethe-Bloch + special treatment for electrons. Adds effects from energy loss according to Bethe-Bloch formula without density effect. Assumes silicon as material. For electrons energy loss due to radiation added according to formulae by Bethe & Heitler. Ported from ORCA.

Definition at line 20 of file EnergyLossUpdator.h.

Constructor & Destructor Documentation

EnergyLossUpdator::EnergyLossUpdator ( double  mass)
inline

Definition at line 28 of file EnergyLossUpdator.h.

Referenced by clone().

Member Function Documentation

virtual EnergyLossUpdator* EnergyLossUpdator::clone ( void  ) const
inlinevirtual

Implements MaterialEffectsUpdator.

Definition at line 23 of file EnergyLossUpdator.h.

References EnergyLossUpdator().

23  {
24  return new EnergyLossUpdator(*this);
25  }
EnergyLossUpdator(double mass)
void EnergyLossUpdator::compute ( const TrajectoryStateOnSurface TSoS,
const PropagationDirection  propDir 
) const
privatevirtual

Implements MaterialEffectsUpdator.

Definition at line 9 of file EnergyLossUpdator.cc.

References alongMomentum, computeBetheBloch(), computeElectrons(), TrajectoryStateOnSurface::localMomentum(), MaterialEffectsUpdator::mass(), Surface::mediumProperties(), TrajectoryStateOnSurface::surface(), MaterialEffectsUpdator::theDeltaCov, and MaterialEffectsUpdator::theDeltaP.

11 {
12  //
13  // Get surface
14  //
15  const Surface& surface = TSoS.surface();
16  //
17  // Initialise dP and the update to the covariance matrix
18  //
19  theDeltaP = 0.;
20  theDeltaCov(0,0) = 0.;
21  //
22  // Now get information on medium
23  //
24  if (surface.mediumProperties()) {
25  //
26  // Bethe-Bloch
27  //
28  if ( mass()>0.001 )
30  //
31  // Special treatment for electrons (currently rather crude
32  // distinction using mass)
33  //
34  else
36  propDir);
37  if (propDir != alongMomentum) theDeltaP *= -1.;
38  }
39 }
void computeElectrons(const LocalVector &, const MediumProperties &, const PropagationDirection) const
const MediumProperties * mediumProperties() const
Definition: Surface.h:93
LocalVector localMomentum() const
const Surface & surface() const
void computeBetheBloch(const LocalVector &, const MediumProperties &) const
AlgebraicSymMatrix55 theDeltaCov
void EnergyLossUpdator::computeBetheBloch ( const LocalVector localP,
const MediumProperties materialConstants 
) const
private

Definition at line 44 of file EnergyLossUpdator.cc.

References beta, alignCSCRings::e, create_public_lumi_plots::log, m, PV3DBase< T, PVType, FrameType >::mag(), MaterialEffectsUpdator::mass(), AlCaHLTBitMon_ParallelJobs::p, mathSSE::sqrt(), MaterialEffectsUpdator::theDeltaCov, MaterialEffectsUpdator::theDeltaP, MediumProperties::xi(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by compute().

45  {
46  //
47  // calculate absolute momentum and correction to path length from angle
48  // of incidence
49  //
50  double p = localP.mag();
51  double xf = fabs(p/localP.z());
52  // constants
53  const double m = mass(); // use mass hypothesis from constructor
54 
55  const double emass = 0.511e-3;
56  const double poti = 16.e-9 * 10.75; // = 16 eV * Z**0.9, for Si Z=14
57  const double eplasma = 28.816e-9 * sqrt(2.33*0.498); // 28.816 eV * sqrt(rho*(Z/A)) for Si
58  const double delta0 = 2*log(eplasma/poti) - 1.;
59 
60  // calculate general physics things
61  double e = sqrt(p*p + m*m);
62  double beta = p/e;
63  double gamma = e/m;
64  double eta2 = beta*gamma; eta2 *= eta2;
65  // double lnEta2 = log(eta2);
66  double ratio = emass/m;
67  double emax = 2.*emass*eta2/(1. + 2.*ratio*gamma + ratio*ratio);
68  // double delta = delta0 + lnEta2;
69 
70  // calculate the mean and sigma of energy loss
71  // xi = d[g/cm2] * 0.307075MeV/(g/cm2) * Z/A * 1/2
72  double xi = materialConstants.xi()*xf; xi /= (beta*beta);
73 
74 // double dEdx = xi*(log(2.*emass*eta2*emax/(poti*poti)) - 2.*(beta*beta));
75  //double dEdx = xi*(log(2.*emass*emax/(poti*poti))+lnEta2 - 2.*(beta*beta) - delta);
76 
77  double dEdx = xi*(log(2.*emass*emax/(poti*poti)) - 2.*(beta*beta) - delta0);
78  double dEdx2 = xi*emax*(1.-0.5*(beta*beta));
79  double dP = dEdx/beta;
80  double sigp2 = dEdx2*e*e/(p*p*p*p*p*p);
81  theDeltaP += -dP;
82  theDeltaCov(0,0) += sigp2;
83 }
const double beta
T mag() const
Definition: PV3DBase.h:66
T sqrt(T t)
Definition: SSEVec.h:46
T z() const
Definition: PV3DBase.h:63
float xi() const
AlgebraicSymMatrix55 theDeltaCov
void EnergyLossUpdator::computeElectrons ( const LocalVector localP,
const MediumProperties materialConstants,
const PropagationDirection  propDir 
) const
private

Definition at line 88 of file EnergyLossUpdator.cc.

References create_public_lumi_plots::exp, f, create_public_lumi_plots::log, PV3DBase< T, PVType, FrameType >::mag(), oppositeToMomentum, AlCaHLTBitMon_ParallelJobs::p, MediumProperties::radLen(), MaterialEffectsUpdator::theDeltaCov, MaterialEffectsUpdator::theDeltaP, detailsBasic3DVector::z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by compute().

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 = fabs(p/localP.z())*materialConstants.radLen();
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  z*z;
104  // exp(-2*normalisedPath);
105 
106  if ( propDir==oppositeToMomentum ) {
107  //
108  // for backward propagation: delta(1/p) is linear in z=p_outside/p_inside
109  // convert to obtain equivalent delta(p). Sign of deltaP is corrected
110  // in method compute -> deltaP<0 at this place!!!
111  //
112  theDeltaP += -p*(1/z-1);
113  theDeltaCov(0,0) += varz/(p*p);
114  }
115  else {
116  //
117  // for forward propagation: calculate in p (linear in 1/z=p_inside/p_outside),
118  // then convert sig(p) to sig(1/p).
119  //
120  theDeltaP += p*(z-1);
121  // double f = 1/p/z/z;
122  // patch to ensure consistency between for- and backward propagation
123  double f = 1./(p*z);
124  theDeltaCov(0,0) += f*f*varz;
125  }
126 }
float radLen() const
double double double z
T mag() const
Definition: PV3DBase.h:66
T z() const
Definition: PV3DBase.h:63
double f[11][100]
AlgebraicSymMatrix55 theDeltaCov