CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonBremsstrahlungSimulator.cc
Go to the documentation of this file.
1 //MuonBremsstrahlungSimulator
3 // Description: Implementation of Muon bremsstrahlung using the Petrukhin model
4 //Authors :Sandro Fonseca de Souza and Andre Sznajder (UERJ/Brazil)
5 // Date: 23-Nov-2010
11 
12 #include <cmath>
13 #include <string>
14 #include <iostream>
15 
16 #include "TF1.h"
17 
18 
19 
20 
23  double A,double Z, double density,double radLen,
24  double photonEnergyCut,double photonFractECut) :
25  MaterialEffectsSimulator(engine,A,Z,density,radLen)
26 {
27  // Set the minimal photon energy for a Brem from mu+/-
28  photonEnergy = photonEnergyCut;
29  photonFractE = photonFractECut;
30  d = 0.; //distance
31  LogDebug("MuonBremsstrahlungSimulator")<< "Starting the MuonBremsstrahlungSimulator"<< std::endl;
32 }
33 
34 
35 
37 void
39 {
40 
41  double NA = 6.022e+23; //Avogadro's number
42 
43  if ( radLengths > 4. ) {
44  Particle.SetXYZT(0.,0.,0.,0.);
45  deltaPMuon.SetXYZT(0.,0.,0.,0.);
46  brem_photon.SetXYZT(0.,0.,0.,0.);
47  }
48 
49 // Hard brem probability with a photon Energy above photonEnergy.
50 
51  double EMuon = Particle.e();//Muon Energy
52  if (EMuon<photonEnergy) return;
53  xmin = std::max(photonEnergy/EMuon,photonFractE);//fraction of muon's energy transferred to the photon
54 
55  // xmax = photonEnergy/Particle.e();
56  if ( xmin >=1. || xmin <=0. ) return;
57 
58  xmax = 1.;
59  npar = 3 ;//Number of parameters
60 
61  // create TF1 using a free C function
62  f1 = new TF1("f1",PetrukhinFunc,xmin,xmax,npar);
63  //Setting parameters
64  f1->SetParameters(EMuon,A,Z);
66  //d = distance for several materials
67  //X0 = radLen
68  //d = radLengths * X0(for tracker,yoke,ECAL and HCAL)
69  d = radLengths * radLen ;
70  //Integration
71  bremProba = density * d *(NA/A)* (f1->Integral(0.,1.));
72 
73 
74  // Number of photons to be radiated.
75  unsigned int nPhotons = poisson(bremProba);
76  _theUpdatedState.reserve(nPhotons);
77 
78 
79  if ( !nPhotons ) return;
80 
81  //Rotate to the lab frame
82  double chi = Particle.theta();
83  double psi = Particle.phi();
84  RawParticle::RotationZ rotZ(psi);
85  RawParticle::RotationY rotY(chi);
86 
87 
88 
89  // Energy of these photons
90  for ( unsigned int i=0; i<nPhotons; ++i ) {
91 
92  // Check that there is enough energy left.
93  if ( Particle.e() < photonEnergy ) break;
94  LogDebug("MuonBremsstrahlungSimulator")<< "MuonBremsstrahlungSimulator parameters:"<< std::endl;
95  LogDebug("MuonBremsstrahlungSimulator")<< "xmin-> " << xmin << std::endl;
96  LogDebug("MuonBremsstrahlungSimulator")<< "Atomic Weight-> " << A << std::endl;
97  LogDebug("MuonBremsstrahlungSimulator")<< "Density-> " << density << std::endl;
98  LogDebug("MuonBremsstrahlungSimulator")<< "Distance-> " << d << std::endl;
99  LogDebug("MuonBremsstrahlungSimulator")<< "bremProba->" << bremProba << std::endl;
100  LogDebug("MuonBremsstrahlungSimulator")<< "nPhotons->" << nPhotons << std::endl;
101  LogDebug("MuonBremsstrahlungSimulator")<< " Muon_Energy-> " << EMuon << std::endl;
102  LogDebug("MuonBremsstrahlungSimulator")<< "X0-> "<< radLen << std::endl;
103  LogDebug("MuonBremsstrahlungSimulator")<< " radLengths-> " << radLengths << std::endl;
104 
105  // Add a photon
106  RawParticle thePhoton(22,brem(Particle));
107  if (thePhoton.E()>0.){
108 
109  thePhoton.rotate(rotY);
110  thePhoton.rotate(rotZ);
111 
112  _theUpdatedState.push_back(thePhoton);
113 
114  // Update the original mu +/-
115  deltaPMuon = Particle -= thePhoton.momentum();
116  // Information of brem photon
117  brem_photon.SetXYZT(thePhoton.Px(),thePhoton.Py(),thePhoton.Pz(),thePhoton.E());
118 
119  LogDebug("MuonBremsstrahlungSimulator")<< " Muon Bremsstrahlung: photon_energy-> " << thePhoton.E() << std::endl;
120  LogDebug("MuonBremsstrahlungSimulator")<< "photon_px->" << thePhoton.Px() << std::endl;
121  LogDebug("MuonBremsstrahlungSimulator")<< "photon_py->" << thePhoton.Py() << std::endl;
122  LogDebug("MuonBremsstrahlungSimulator")<< "photon_pz->" << thePhoton.Pz() << std::endl;
123 
124  }
125 
126  }
127 
128 }
132 
133  // This is a simple version of a Muon Brem using Petrukhin model .
134  //Ref: http://pdg.lbl.gov/2008/AtomicNuclearProperties/adndt.pdf
135  double mumass = 0.105658367;//mu mass (GeV/c^2)
136  double xp =0.;
137  //+++++++++++++++++++++++++++++++++++++++++++++++++++
138  // Ref: Framework
139  //https://hypernews.cern.ch/HyperNews/CMS/get/swDevelopment/1969/1.html
140  // http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/DQM/DTMonitorClient/src/DTResolutionTest.cc?view=markup
142  try {
143 
144  xp = f1->GetRandom();
145  } catch (...) {
146 
147  edm::LogInfo ("MuonBremsstrahlungSimulator") << "Exception when fitting...";
148  }
149  LogDebug("MuonBremsstrahlungSimulator")<< "MuonBremsstrahlungSimulator: xp->" << xp << std::endl;
150 
151 
152  // Have photon energy. Now generate angles with respect to the z axis
153  // defined by the incoming particle's momentum.
154 
155  // Isotropic in phi
156  const double phi = random->flatShoot()*2*M_PI;
157  // theta from universal distribution
158  const double theta = gbteth(pp.e(),mumass,xp)*mumass/pp.e();
159 
160  // Make momentum components
161  double stheta = std::sin(theta);
162  double ctheta = std::cos(theta);
163  double sphi = std::sin(phi);
164  double cphi = std::cos(phi);
165 
166  return xp * pp.e() * XYZTLorentzVector(stheta*cphi,stheta*sphi,ctheta,1.);
167 
168 }
170 double
172  const double partm,
173  const double efrac) const {
174  const double alfa = 0.625;
175 
176  const double d = 0.13*(0.8+1.3/theZ())*(100.0+(1.0/ener))*(1.0+efrac);
177  const double w1 = 9.0/(9.0+d);
178  const double umax = ener*M_PI/partm;
179  double u;
180 
181  do {
182  double beta = (random->flatShoot()<=w1) ? alfa : 3.0*alfa;
183  u = -std::log(random->flatShoot()*random->flatShoot())/beta;
184  } while (u>=umax);
185 
186  return u;
187 }
188 
190 unsigned int
192 
193  unsigned int n = 0;
194  double prob = std::exp(-ymu);
195  double proba = prob;
196  double x = random->flatShoot();
197 
198  while ( proba <= x ) {
199  prob *= ymu / double(++n);
200  proba += prob;
201  }
202 
203  return n;
204 
205 }
206 
#define LogDebug(id)
const double Z[kNumberCalorimeter]
const double beta
int i
Definition: DBlmapReader.cc:9
tuple pp
Definition: createTree.py:15
double PetrukhinFunc(double *x, double *p)
double xmin
The fractional photon energy cut (determined from the above two)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
MuonBremsstrahlungSimulator(const RandomEngine *engine, double A, double Z, double density, double radLen, double photonEnergyCut, double photonFractECut)
Constructor.
XYZTLorentzVector brem(ParticlePropagator &p) const
Compute Brem photon energy and angles, if any.
ROOT::Math::RotationZ RotationZ
Definition: RawParticle.h:38
std::map< std::string, int, std::less< std::string > > psi
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:285
const T & max(const T &a, const T &b)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
ROOT::Math::RotationY RotationY
Definition: RawParticle.h:37
void rotate(double rphi, const XYZVector &raxis)
Definition: RawParticle.cc:154
unsigned int poisson(double ymu)
Generate numbers according to a Poisson distribution of mean ymu.
double gbteth(const double ener, const double partm, const double efrac) const
A universal angular distribution - still from GEANT.
#define M_PI
Definition: BFit3D.cc:3
double flatShoot(double xmin=0.0, double xmax=1.0) const
Definition: RandomEngine.h:30
double photonFractE
The minimum photon fractional energy (wrt that of the electron)
std::vector< RawParticle > _theUpdatedState
void compute(ParticlePropagator &Particle)
Generate Bremsstrahlung photons.
double photonEnergy
The minimum photon energy to be radiated, in GeV.
Definition: DDAxes.h:10
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
Definition: DDAxes.h:10