CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
MuonBremsstrahlungSimulator Class Reference

#include <MuonBremsstrahlungSimulator.h>

Inheritance diagram for MuonBremsstrahlungSimulator:
MaterialEffectsSimulator

Public Member Functions

const XYZTLorentzVectordeltaP_BremMuon () const
 
const XYZTLorentzVectordeltaP_BremPhoton () const
 
 MuonBremsstrahlungSimulator (double A, double Z, double density, double radLen, double photonEnergyCut, double photonFractECut)
 Constructor. More...
 
 ~MuonBremsstrahlungSimulator () override
 Default destructor. More...
 
- Public Member Functions inherited from MaterialEffectsSimulator
RHEP_const_iter beginDaughters () const
 Returns const iterator to the beginning of the daughters list. More...
 
int closestDaughterId ()
 The id of the closest charged daughter (filled for nuclear interactions only) More...
 
double eMass () const
 Electron mass in GeV/c2. More...
 
RHEP_const_iter endDaughters () const
 Returns const iterator to the end of the daughters list. More...
 
double excitE () const
 Mean excitation energy (in GeV) More...
 
 MaterialEffectsSimulator (double A=28.0855, double Z=14.0000, double density=2.329, double radLen=9.360)
 
unsigned nDaughters () const
 Returns the number of daughters. More...
 
XYZVector orthogonal (const XYZVector &) const
 A vector orthogonal to another one (because it's not in XYZTLorentzVector) More...
 
double radLenIncm () const
 One radiation length in cm. More...
 
double rho () const
 Density in g/cm3. More...
 
virtual void save ()
 Used by NuclearInteractionSimulator to save last sampled event. More...
 
void setNormalVector (const GlobalVector &normal)
 Sets the vector normal to the surface traversed. More...
 
double theA () const
 A. More...
 
double theZ () const
 Z. More...
 
void updateState (ParticlePropagator &myTrack, double radlen, RandomEngineAndDistribution const *)
 Compute the material effect (calls the sub class) More...
 
virtual ~MaterialEffectsSimulator ()
 

Private Member Functions

XYZTLorentzVector brem (ParticlePropagator &p, RandomEngineAndDistribution const *) const
 Compute Brem photon energy and angles, if any. More...
 
void compute (ParticlePropagator &Particle, RandomEngineAndDistribution const *) override
 Generate Bremsstrahlung photons. More...
 
double gbteth (const double ener, const double partm, const double efrac, RandomEngineAndDistribution const *) const
 A universal angular distribution - still from GEANT. More...
 
unsigned int poisson (double ymu)
 Generate numbers according to a Poisson distribution of mean ymu. More...
 

Private Attributes

XYZTLorentzVector brem_photon
 
double bremProba
 
double d
 
XYZTLorentzVector deltaPMuon
 
TF1 * f1
 
int npar
 
double photonEnergy
 The minimum photon energy to be radiated, in GeV. More...
 
double photonFractE
 The minimum photon fractional energy (wrt that of the electron) More...
 
double rand
 
double xmax
 
double xmin
 The fractional photon energy cut (determined from the above two) More...
 

Additional Inherited Members

- Public Types inherited from MaterialEffectsSimulator
typedef std::vector< RawParticle >::const_iterator RHEP_const_iter
 
- Protected Attributes inherited from MaterialEffectsSimulator
std::vector< RawParticle_theUpdatedState
 
double A
 
double density
 
double radLen
 
double radLengths
 
int theClosestChargedDaughterId
 
GlobalVector theNormalVector
 
double Z
 

Detailed Description

Definition at line 31 of file MuonBremsstrahlungSimulator.h.

Constructor & Destructor Documentation

MuonBremsstrahlungSimulator::MuonBremsstrahlungSimulator ( double  A,
double  Z,
double  density,
double  radLen,
double  photonEnergyCut,
double  photonFractECut 
)

Constructor.

Definition at line 23 of file MuonBremsstrahlungSimulator.cc.

References d, LogDebug, photonEnergy, and photonFractE.

24  :
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 }
#define LogDebug(id)
MaterialEffectsSimulator(double A=28.0855, double Z=14.0000, double density=2.329, double radLen=9.360)
double photonFractE
The minimum photon fractional energy (wrt that of the electron)
double photonEnergy
The minimum photon energy to be radiated, in GeV.
MuonBremsstrahlungSimulator::~MuonBremsstrahlungSimulator ( )
inlineoverride

Default destructor.

Definition at line 40 of file MuonBremsstrahlungSimulator.h.

40 {}

Member Function Documentation

XYZTLorentzVector MuonBremsstrahlungSimulator::brem ( ParticlePropagator p,
RandomEngineAndDistribution const *  random 
) const
private

Compute Brem photon energy and angles, if any.

Definition at line 131 of file MuonBremsstrahlungSimulator.cc.

References funct::cos(), gather_cfg::cout, RawParticle::e(), f1, RandomEngineAndDistribution::flatShoot(), gbteth(), LogDebug, M_PI, BaseParticlePropagator::particle(), phi, random, funct::sin(), and theta().

Referenced by compute().

131  {
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 
137  // RANDOM_NUMBER_ERROR
138  // Random number should be generated by the engines from the
139  // RandomNumberGeneratorService. This appears to use the global
140  // engine in ROOT. This is not thread safe unless the module using
141  // it is a one module and declares a shared resource and all
142  // other modules using it also declare the same shared resource.
143  // This also breaks replay.
144  // It might be that this is never used because of the std::cout
145  // statement 2 lines below. I've never seen or heard complaints
146  // about that vebosity ....
147  double xp = f1->GetRandom();
148  LogDebug("MuonBremsstrahlungSimulator")<< "MuonBremsstrahlungSimulator: xp->" << xp << std::endl;
149  std::cout << "MuonBremsstrahlungSimulator: xp->" << xp << std::endl;
150 
151  // Have photon energy. Now generate angles with respect to the z axis
152  // defined by the incoming particle's momentum.
153 
154  // Isotropic in phi
155  const double phi = random->flatShoot()*2*M_PI;
156  // theta from universal distribution
157  const double theta = gbteth(pp.particle().e(),mumass,xp,random)*mumass/pp.particle().e();
158 
159  // Make momentum components
160  double stheta = std::sin(theta);
161  double ctheta = std::cos(theta);
162  double sphi = std::sin(phi);
163  double cphi = std::cos(phi);
164 
165  return xp * pp.particle().e() * XYZTLorentzVector(stheta*cphi,stheta*sphi,ctheta,1.);
166 
167 }
#define LogDebug(id)
double gbteth(const double ener, const double partm, const double efrac, RandomEngineAndDistribution const *) const
A universal angular distribution - still from GEANT.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
TRandom random
Definition: MVATrainer.cc:138
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
#define M_PI
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:27
void MuonBremsstrahlungSimulator::compute ( ParticlePropagator Particle,
RandomEngineAndDistribution const *  random 
)
overrideprivatevirtual

Generate Bremsstrahlung photons.

Implements MaterialEffectsSimulator.

Definition at line 38 of file MuonBremsstrahlungSimulator.cc.

References MaterialEffectsSimulator::_theUpdatedState, MaterialEffectsSimulator::A, brem(), brem_photon, bremProba, d, deltaPMuon, MaterialEffectsSimulator::density, RawParticle::e(), RawParticle::E(), f1, mps_fire::i, LogDebug, makeParticle(), SiStripPI::max, RawParticle::momentum(), fastsim::Constants::NA, npar, BaseParticlePropagator::particle(), ParticlePropagator::particleDataTable(), PetrukhinFunc(), RawParticle::phi(), photonEnergy, photonFractE, RandomEngineAndDistribution::poissonShoot(), RawParticle::Px(), RawParticle::Py(), RawParticle::Pz(), MaterialEffectsSimulator::radLen, MaterialEffectsSimulator::radLengths, RawParticle::rotate(), RawParticle::setMomentum(), RawParticle::theta(), xmax, xmin, and MaterialEffectsSimulator::Z.

39 {
40 
41  double NA = 6.022e+23; //Avogadro's number
42 
43  if ( radLengths > 4. ) {
44  Particle.particle().setMomentum(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.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.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 = random->poissonShoot(bremProba);
76  _theUpdatedState.reserve(nPhotons);
77 
78 
79  if ( !nPhotons ) return;
80 
81  //Rotate to the lab frame
82  double chi = Particle.particle().theta();
83  double psi = Particle.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.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=makeParticle(Particle.particleDataTable(),22,brem(Particle, random));
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.particle().momentum() -= 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 }
#define LogDebug(id)
static double constexpr NA
Avogadro&#39;s number.
Definition: Constants.h:20
void setMomentum(const XYZTLorentzVector &vtx)
set the momentum
Definition: RawParticle.h:347
const HepPDT::ParticleDataTable * particleDataTable() const
double PetrukhinFunc(double *x, double *p)
double xmin
The fractional photon energy cut (determined from the above two)
TRandom random
Definition: MVATrainer.cc:138
RawParticle const & particle() const
The particle being propagated.
std::map< std::string, int, std::less< std::string > > psi
ROOT::Math::RotationZ RotationZ
Definition: RawParticle.h:51
RawParticle makeParticle(HepPDT::ParticleDataTable const *, int id, const math::XYZTLorentzVector &p)
Definition: makeParticle.cc:29
double phi() const
phi of momentum vector
Definition: RawParticle.h:335
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:340
double e() const
energy of the momentum
Definition: RawParticle.h:324
double Py() const
y of the momentum
Definition: RawParticle.h:319
ROOT::Math::RotationY RotationY
Definition: RawParticle.h:50
void rotate(double rphi, const XYZVector &raxis)
Definition: RawParticle.cc:83
double Pz() const
z of the momentum
Definition: RawParticle.h:322
double photonFractE
The minimum photon fractional energy (wrt that of the electron)
double theta() const
theta of momentum vector
Definition: RawParticle.h:334
std::vector< RawParticle > _theUpdatedState
double Px() const
x of the momentum
Definition: RawParticle.h:316
double E() const
energy of the momentum
Definition: RawParticle.h:325
XYZTLorentzVector brem(ParticlePropagator &p, RandomEngineAndDistribution const *) const
Compute Brem photon energy and angles, if any.
double photonEnergy
The minimum photon energy to be radiated, in GeV.
const XYZTLorentzVector& MuonBremsstrahlungSimulator::deltaP_BremMuon ( ) const
inline

Definition at line 44 of file MuonBremsstrahlungSimulator.h.

References deltaPMuon.

44 { return deltaPMuon ; }
const XYZTLorentzVector& MuonBremsstrahlungSimulator::deltaP_BremPhoton ( ) const
inline

Definition at line 47 of file MuonBremsstrahlungSimulator.h.

References brem_photon.

47 { return brem_photon ; }
double MuonBremsstrahlungSimulator::gbteth ( const double  ener,
const double  partm,
const double  efrac,
RandomEngineAndDistribution const *  random 
) const
private

A universal angular distribution - still from GEANT.

Definition at line 170 of file MuonBremsstrahlungSimulator.cc.

References pfBoostedDoubleSVAK8TagInfos_cfi::beta, d, RandomEngineAndDistribution::flatShoot(), cmsBatch::log, M_PI, and MaterialEffectsSimulator::theZ().

Referenced by brem().

173  {
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 }
TRandom random
Definition: MVATrainer.cc:138
#define M_PI
unsigned int MuonBremsstrahlungSimulator::poisson ( double  ymu)
private

Generate numbers according to a Poisson distribution of mean ymu.

Member Data Documentation

XYZTLorentzVector MuonBremsstrahlungSimulator::brem_photon
private

Definition at line 84 of file MuonBremsstrahlungSimulator.h.

Referenced by compute(), and deltaP_BremPhoton().

double MuonBremsstrahlungSimulator::bremProba
private

Definition at line 58 of file MuonBremsstrahlungSimulator.h.

Referenced by compute().

double MuonBremsstrahlungSimulator::d
private

Definition at line 64 of file MuonBremsstrahlungSimulator.h.

Referenced by compute(), gbteth(), and MuonBremsstrahlungSimulator().

XYZTLorentzVector MuonBremsstrahlungSimulator::deltaPMuon
private

Definition at line 82 of file MuonBremsstrahlungSimulator.h.

Referenced by compute(), and deltaP_BremMuon().

TF1* MuonBremsstrahlungSimulator::f1
private

Definition at line 52 of file MuonBremsstrahlungSimulator.h.

Referenced by brem(), and compute().

int MuonBremsstrahlungSimulator::npar
private

Definition at line 54 of file MuonBremsstrahlungSimulator.h.

Referenced by compute().

double MuonBremsstrahlungSimulator::photonEnergy
private

The minimum photon energy to be radiated, in GeV.

Definition at line 57 of file MuonBremsstrahlungSimulator.h.

Referenced by compute(), MuonBremsstrahlungSimulator(), and Jet.Jet::photonEnergyFraction().

double MuonBremsstrahlungSimulator::photonFractE
private

The minimum photon fractional energy (wrt that of the electron)

Definition at line 60 of file MuonBremsstrahlungSimulator.h.

Referenced by compute(), and MuonBremsstrahlungSimulator().

double MuonBremsstrahlungSimulator::rand
private

Definition at line 63 of file MuonBremsstrahlungSimulator.h.

double MuonBremsstrahlungSimulator::xmax
private

Definition at line 63 of file MuonBremsstrahlungSimulator.h.

Referenced by svgfig.XAxis::__repr__(), and compute().

double MuonBremsstrahlungSimulator::xmin
private

The fractional photon energy cut (determined from the above two)

Definition at line 63 of file MuonBremsstrahlungSimulator.h.

Referenced by svgfig.XAxis::__repr__(), svgfig.Axes::__repr__(), svgfig.HGrid::__repr__(), svgfig.Grid::__repr__(), compute(), and svgfig.Axes::SVG().