Go to the documentation of this file.00001
00002
00003
00004
00005
00007 #include "FastSimulation/MaterialEffects/interface/MuonBremsstrahlungSimulator.h"
00008 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 #include "FastSimulation/MaterialEffects/interface/PetrukhinModel.h"
00011
00012 #include <cmath>
00013 #include <string>
00014 #include <iostream>
00015
00016 #include "TF1.h"
00017
00018
00019
00020
00022 MuonBremsstrahlungSimulator::MuonBremsstrahlungSimulator(const RandomEngine* engine,
00023 double A,double Z, double density,double radLen,
00024 double photonEnergyCut,double photonFractECut) :
00025 MaterialEffectsSimulator(engine,A,Z,density,radLen)
00026 {
00027
00028 photonEnergy = photonEnergyCut;
00029 photonFractE = photonFractECut;
00030 d = 0.;
00031 LogDebug("MuonBremsstrahlungSimulator")<< "Starting the MuonBremsstrahlungSimulator"<< std::endl;
00032 }
00033
00034
00035
00037 void
00038 MuonBremsstrahlungSimulator::compute(ParticlePropagator &Particle)
00039 {
00040
00041 double NA = 6.022e+23;
00042
00043 if ( radLengths > 4. ) {
00044 Particle.SetXYZT(0.,0.,0.,0.);
00045 deltaPMuon.SetXYZT(0.,0.,0.,0.);
00046 brem_photon.SetXYZT(0.,0.,0.,0.);
00047 }
00048
00049
00050
00051 double EMuon = Particle.e();
00052 if (EMuon<photonEnergy) return;
00053 xmin = std::max(photonEnergy/EMuon,photonFractE);
00054
00055
00056 if ( xmin >=1. || xmin <=0. ) return;
00057
00058 xmax = 1.;
00059 npar = 3 ;
00060
00061
00062 f1 = new TF1("f1",PetrukhinFunc,xmin,xmax,npar);
00063
00064 f1->SetParameters(EMuon,A,Z);
00066
00067
00068
00069 d = radLengths * radLen ;
00070
00071 bremProba = density * d *(NA/A)* (f1->Integral(0.,1.));
00072
00073
00074
00075 unsigned int nPhotons = random->poissonShoot(bremProba);
00076 _theUpdatedState.reserve(nPhotons);
00077
00078
00079 if ( !nPhotons ) return;
00080
00081
00082 double chi = Particle.theta();
00083 double psi = Particle.phi();
00084 RawParticle::RotationZ rotZ(psi);
00085 RawParticle::RotationY rotY(chi);
00086
00087
00088
00089
00090 for ( unsigned int i=0; i<nPhotons; ++i ) {
00091
00092
00093 if ( Particle.e() < photonEnergy ) break;
00094 LogDebug("MuonBremsstrahlungSimulator")<< "MuonBremsstrahlungSimulator parameters:"<< std::endl;
00095 LogDebug("MuonBremsstrahlungSimulator")<< "xmin-> " << xmin << std::endl;
00096 LogDebug("MuonBremsstrahlungSimulator")<< "Atomic Weight-> " << A << std::endl;
00097 LogDebug("MuonBremsstrahlungSimulator")<< "Density-> " << density << std::endl;
00098 LogDebug("MuonBremsstrahlungSimulator")<< "Distance-> " << d << std::endl;
00099 LogDebug("MuonBremsstrahlungSimulator")<< "bremProba->" << bremProba << std::endl;
00100 LogDebug("MuonBremsstrahlungSimulator")<< "nPhotons->" << nPhotons << std::endl;
00101 LogDebug("MuonBremsstrahlungSimulator")<< " Muon_Energy-> " << EMuon << std::endl;
00102 LogDebug("MuonBremsstrahlungSimulator")<< "X0-> "<< radLen << std::endl;
00103 LogDebug("MuonBremsstrahlungSimulator")<< " radLengths-> " << radLengths << std::endl;
00104
00105
00106 RawParticle thePhoton(22,brem(Particle));
00107 if (thePhoton.E()>0.){
00108
00109 thePhoton.rotate(rotY);
00110 thePhoton.rotate(rotZ);
00111
00112 _theUpdatedState.push_back(thePhoton);
00113
00114
00115 deltaPMuon = Particle -= thePhoton.momentum();
00116
00117 brem_photon.SetXYZT(thePhoton.Px(),thePhoton.Py(),thePhoton.Pz(),thePhoton.E());
00118
00119 LogDebug("MuonBremsstrahlungSimulator")<< " Muon Bremsstrahlung: photon_energy-> " << thePhoton.E() << std::endl;
00120 LogDebug("MuonBremsstrahlungSimulator")<< "photon_px->" << thePhoton.Px() << std::endl;
00121 LogDebug("MuonBremsstrahlungSimulator")<< "photon_py->" << thePhoton.Py() << std::endl;
00122 LogDebug("MuonBremsstrahlungSimulator")<< "photon_pz->" << thePhoton.Pz() << std::endl;
00123
00124 }
00125
00126 }
00127
00128 }
00130 XYZTLorentzVector
00131 MuonBremsstrahlungSimulator::brem(ParticlePropagator& pp) const {
00132
00133
00134
00135 double mumass = 0.105658367;
00136 double xp = f1->GetRandom();
00137 LogDebug("MuonBremsstrahlungSimulator")<< "MuonBremsstrahlungSimulator: xp->" << xp << std::endl;
00138 std::cout << "MuonBremsstrahlungSimulator: xp->" << xp << std::endl;
00139
00140
00141
00142
00143
00144
00145 const double phi = random->flatShoot()*2*M_PI;
00146
00147 const double theta = gbteth(pp.e(),mumass,xp)*mumass/pp.e();
00148
00149
00150 double stheta = std::sin(theta);
00151 double ctheta = std::cos(theta);
00152 double sphi = std::sin(phi);
00153 double cphi = std::cos(phi);
00154
00155 return xp * pp.e() * XYZTLorentzVector(stheta*cphi,stheta*sphi,ctheta,1.);
00156
00157 }
00159 double
00160 MuonBremsstrahlungSimulator::gbteth(const double ener,
00161 const double partm,
00162 const double efrac) const {
00163 const double alfa = 0.625;
00164
00165 const double d = 0.13*(0.8+1.3/theZ())*(100.0+(1.0/ener))*(1.0+efrac);
00166 const double w1 = 9.0/(9.0+d);
00167 const double umax = ener*M_PI/partm;
00168 double u;
00169
00170 do {
00171 double beta = (random->flatShoot()<=w1) ? alfa : 3.0*alfa;
00172 u = -std::log(random->flatShoot()*random->flatShoot())/beta;
00173 } while (u>=umax);
00174
00175 return u;
00176 }
00177
00178