00001
00002 #include "FastSimulation/MaterialEffects/interface/BremsstrahlungSimulator.h"
00003 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00004
00005 #include <cmath>
00006
00007 BremsstrahlungSimulator::BremsstrahlungSimulator(double photonEnergyCut,
00008 double photonFractECut,
00009 const RandomEngine* engine) :
00010 MaterialEffectsSimulator(engine)
00011 {
00012
00013
00014 photonEnergy = photonEnergyCut;
00015 photonFractE = photonFractECut;
00016
00017 }
00018
00019
00020 void
00021 BremsstrahlungSimulator::compute(ParticlePropagator &Particle)
00022 {
00023
00024
00025
00026
00027
00028
00029 if ( radLengths > 4. ) Particle.SetXYZT(0.,0.,0.,0.);
00030
00031
00032 xmin = std::max(photonEnergy/Particle.e(),photonFractE);
00033 if ( xmin >=1. || xmin <=0. ) return;
00034
00035 double bremProba = radLengths * ( 4./3. * std::log(1./xmin)
00036 - 4./3. * (1.-xmin)
00037 + 1./2. * (1.-xmin*xmin) );
00038
00039
00040
00041 unsigned int nPhotons = poisson(bremProba);
00042 _theUpdatedState.reserve(nPhotons);
00043
00044 if ( !nPhotons ) return;
00045
00046
00047 double chi = Particle.theta();
00048 double psi = Particle.phi();
00049 RawParticle::RotationZ rotZ(psi);
00050 RawParticle::RotationY rotY(chi);
00051
00052
00053 for ( unsigned int i=0; i<nPhotons; ++i ) {
00054
00055
00056 if ( Particle.e() < photonEnergy ) break;
00057
00058
00059 RawParticle thePhoton(22,brem(Particle));
00060 thePhoton.rotate(rotY);
00061 thePhoton.rotate(rotZ);
00062 _theUpdatedState.push_back(thePhoton);
00063
00064
00065 Particle -= thePhoton.momentum();
00066
00067 }
00068 }
00069
00070 XYZTLorentzVector
00071 BremsstrahlungSimulator::brem(ParticlePropagator& pp) const {
00072
00073
00074
00075
00076 double emass = 0.0005109990615;
00077 double xp=0;
00078 double weight = 0.;
00079
00080 do {
00081 xp = xmin * std::exp ( -std::log(xmin) * random->flatShoot() );
00082 weight = 1. - xp + 3./4.*xp*xp;
00083 } while ( weight < random->flatShoot() );
00084
00085
00086
00087
00088
00089
00090 const double phi = random->flatShoot()*2*M_PI;
00091
00092 const double theta = gbteth(pp.e(),emass,xp)*emass/pp.e();
00093
00094
00095 double stheta = std::sin(theta);
00096 double ctheta = std::cos(theta);
00097 double sphi = std::sin(phi);
00098 double cphi = std::cos(phi);
00099
00100 return xp * pp.e() * XYZTLorentzVector(stheta*cphi,stheta*sphi,ctheta,1.);
00101
00102 }
00103
00104 double
00105 BremsstrahlungSimulator::gbteth(const double ener,
00106 const double partm,
00107 const double efrac) const {
00108 const double alfa = 0.625;
00109
00110 const double d = 0.13*(0.8+1.3/theZ())*(100.0+(1.0/ener))*(1.0+efrac);
00111 const double w1 = 9.0/(9.0+d);
00112 const double umax = ener*M_PI/partm;
00113 double u;
00114
00115 do {
00116 double beta = (random->flatShoot()<=w1) ? alfa : 3.0*alfa;
00117 u = -std::log(random->flatShoot()*random->flatShoot())/beta;
00118 } while (u>=umax);
00119
00120 return u;
00121 }
00122
00123
00124 unsigned int
00125 BremsstrahlungSimulator::poisson(double ymu) {
00126
00127 unsigned int n = 0;
00128 double prob = std::exp(-ymu);
00129 double proba = prob;
00130 double x = random->flatShoot();
00131
00132 while ( proba <= x ) {
00133 prob *= ymu / double(++n);
00134 proba += prob;
00135 }
00136
00137 return n;
00138
00139 }
00140