CMS 3D CMS Logo

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

#include <BremsstrahlungSimulator.h>

Inheritance diagram for BremsstrahlungSimulator:
MaterialEffectsSimulator

Public Member Functions

 BremsstrahlungSimulator (double photonEnergyCut, double photonFractECut)
 Constructor. More...
 
 ~BremsstrahlungSimulator () 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, RandomEngineAndDistribution const *)
 Generate numbers according to a Poisson distribution of mean ymu. More...
 

Private Attributes

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 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 26 of file BremsstrahlungSimulator.h.

Constructor & Destructor Documentation

BremsstrahlungSimulator::BremsstrahlungSimulator ( double  photonEnergyCut,
double  photonFractECut 
)

Constructor.

Definition at line 8 of file BremsstrahlungSimulator.cc.

References photonEnergy, and photonFractE.

10 {
11  // Set the minimal photon energy for a Brem from e+/-
12  photonEnergy = photonEnergyCut;
13  photonFractE = photonFractECut;
14 }
double photonEnergy
The minimum photon energy to be radiated, in GeV.
double photonFractE
The minimum photon fractional energy (wrt that of the electron)
BremsstrahlungSimulator::~BremsstrahlungSimulator ( )
inlineoverride

Default destructor.

Definition at line 35 of file BremsstrahlungSimulator.h.

35 {}

Member Function Documentation

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

Compute Brem photon energy and angles, if any.

Definition at line 69 of file BremsstrahlungSimulator.cc.

References funct::cos(), RawParticle::e(), JetChargeProducer_cfi::exp, RandomEngineAndDistribution::flatShoot(), gbteth(), cmsBatch::log, M_PI, BaseParticlePropagator::particle(), phi, random, funct::sin(), theta(), and xmin.

Referenced by compute().

69  {
70 
71  // This is a simple version (a la PDG) of a Brem generator.
72  // It replaces the buggy GEANT3 -> C++ former version.
73  // Author : Patrick Janot - 25-Dec-2003
74  double emass = 0.0005109990615;
75  double xp=0;
76  double weight = 0.;
77 
78  do {
79  xp = xmin * std::exp ( -std::log(xmin) * random->flatShoot() );
80  weight = 1. - xp + 3./4.*xp*xp;
81  } while ( weight < random->flatShoot() );
82 
83 
84  // Have photon energy. Now generate angles with respect to the z axis
85  // defined by the incoming particle's momentum.
86 
87  // Isotropic in phi
88  const double phi = random->flatShoot()*2*M_PI;
89  // theta from universal distribution
90  const double theta = gbteth(pp.particle().e(),emass,xp,random)*emass/pp.particle().e();
91 
92  // Make momentum components
93  double stheta = std::sin(theta);
94  double ctheta = std::cos(theta);
95  double sphi = std::sin(phi);
96  double cphi = std::cos(phi);
97 
98  return xp * pp.particle().e() * XYZTLorentzVector(stheta*cphi,stheta*sphi,ctheta,1.);
99 
100 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
Definition: weight.py:1
TRandom random
Definition: MVATrainer.cc:138
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
#define M_PI
double gbteth(const double ener, const double partm, const double efrac, RandomEngineAndDistribution const *) const
A universal angular distribution - still from GEANT.
double xmin
The fractional photon energy cut (determined from the above two)
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:27
void BremsstrahlungSimulator::compute ( ParticlePropagator Particle,
RandomEngineAndDistribution const *  random 
)
overrideprivatevirtual

Generate Bremsstrahlung photons.

Implements MaterialEffectsSimulator.

Definition at line 18 of file BremsstrahlungSimulator.cc.

References MaterialEffectsSimulator::_theUpdatedState, brem(), RawParticle::e(), mps_fire::i, cmsBatch::log, makeParticle(), SiStripPI::max, RawParticle::momentum(), BaseParticlePropagator::particle(), ParticlePropagator::particleDataTable(), RawParticle::phi(), photonEnergy, photonFractE, poisson(), MaterialEffectsSimulator::radLengths, RawParticle::rotate(), RawParticle::setMomentum(), RawParticle::theta(), and xmin.

19 {
20 
21  // Protection : Just stop the electron if more than 1 radiation lengths.
22  // This case corresponds to an electron entering the layer parallel to
23  // the layer axis - no reliable simulation can be done in that case...
24  // 08/02/06 - pv: increase protection from 1 to 4 X0 for eta>4.8 region
25  // if ( radLengths > 1. ) Particle.particle().setMomentum(0.,0.,0.,0.);
26  if ( radLengths > 4. ) Particle.particle().setMomentum(0.,0.,0.,0.);
27 
28  // Hard brem probability with a photon Energy above photonEnergy.
29  if (Particle.particle().e()<photonEnergy) return;
31  if ( xmin >=1. || xmin <=0. ) return;
32 
33  double bremProba = radLengths * ( 4./3. * std::log(1./xmin)
34  - 4./3. * (1.-xmin)
35  + 1./2. * (1.-xmin*xmin) );
36 
37 
38  // Number of photons to be radiated.
39  unsigned int nPhotons = poisson(bremProba, random);
40  _theUpdatedState.reserve(nPhotons);
41 
42  if ( !nPhotons ) return;
43 
44  //Rotate to the lab frame
45  double chi = Particle.particle().theta();
46  double psi = Particle.particle().phi();
47  RawParticle::RotationZ rotZ(psi);
48  RawParticle::RotationY rotY(chi);
49 
50  // Energy of these photons
51  for ( unsigned int i=0; i<nPhotons; ++i ) {
52 
53  // Check that there is enough energy left.
54  if ( Particle.particle().e() < photonEnergy ) break;
55 
56  // Add a photon
57  RawParticle thePhoton = makeParticle(Particle.particleDataTable(),22,brem(Particle, random));
58  thePhoton.rotate(rotY);
59  thePhoton.rotate(rotZ);
60  _theUpdatedState.push_back(thePhoton);
61 
62  // Update the original e+/-
63  Particle.particle().momentum() -= thePhoton.momentum();
64 
65  }
66 }
void setMomentum(const XYZTLorentzVector &vtx)
set the momentum
Definition: RawParticle.h:347
const HepPDT::ParticleDataTable * particleDataTable() const
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
unsigned int poisson(double ymu, RandomEngineAndDistribution const *)
Generate numbers according to a Poisson distribution of mean ymu.
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:340
double e() const
energy of the momentum
Definition: RawParticle.h:324
ROOT::Math::RotationY RotationY
Definition: RawParticle.h:50
void rotate(double rphi, const XYZVector &raxis)
Definition: RawParticle.cc:83
double theta() const
theta of momentum vector
Definition: RawParticle.h:334
double xmin
The fractional photon energy cut (determined from the above two)
std::vector< RawParticle > _theUpdatedState
double photonEnergy
The minimum photon energy to be radiated, in GeV.
XYZTLorentzVector brem(ParticlePropagator &p, RandomEngineAndDistribution const *) const
Compute Brem photon energy and angles, if any.
double photonFractE
The minimum photon fractional energy (wrt that of the electron)
double BremsstrahlungSimulator::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 103 of file BremsstrahlungSimulator.cc.

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

Referenced by brem().

106  {
107  const double alfa = 0.625;
108 
109  const double d = 0.13*(0.8+1.3/theZ())*(100.0+(1.0/ener))*(1.0+efrac);
110  const double w1 = 9.0/(9.0+d);
111  const double umax = ener*M_PI/partm;
112  double u;
113 
114  do {
115  double beta = (random->flatShoot()<=w1) ? alfa : 3.0*alfa;
116  u = -std::log(random->flatShoot()*random->flatShoot())/beta;
117  } while (u>=umax);
118 
119  return u;
120 }
TRandom random
Definition: MVATrainer.cc:138
#define M_PI
unsigned int BremsstrahlungSimulator::poisson ( double  ymu,
RandomEngineAndDistribution const *  random 
)
private

Generate numbers according to a Poisson distribution of mean ymu.

Definition at line 124 of file BremsstrahlungSimulator.cc.

References JetChargeProducer_cfi::exp, RandomEngineAndDistribution::flatShoot(), gen::n, TtFullHadEvtBuilder_cfi::prob, and x.

Referenced by compute().

124  {
125 
126  unsigned int n = 0;
127  double prob = std::exp(-ymu);
128  double proba = prob;
129  double x = random->flatShoot();
130 
131  while ( proba <= x ) {
132  prob *= ymu / double(++n);
133  proba += prob;
134  }
135 
136  return n;
137 
138 }
TRandom random
Definition: MVATrainer.cc:138

Member Data Documentation

double BremsstrahlungSimulator::photonEnergy
private

The minimum photon energy to be radiated, in GeV.

Definition at line 40 of file BremsstrahlungSimulator.h.

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

double BremsstrahlungSimulator::photonFractE
private

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

Definition at line 43 of file BremsstrahlungSimulator.h.

Referenced by BremsstrahlungSimulator(), and compute().

double BremsstrahlungSimulator::xmin
private

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

Definition at line 46 of file BremsstrahlungSimulator.h.

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