CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
fastsim::Bremsstrahlung Class Reference

Implementation of Bremsstrahlung from e+/e- in the tracker layers. More...

Inheritance diagram for fastsim::Bremsstrahlung:
fastsim::InteractionModel

Public Member Functions

 Bremsstrahlung (const std::string &name, const edm::ParameterSet &cfg)
 Constructor. More...
 
void interact (Particle &particle, const SimplifiedGeometry &layer, std::vector< std::unique_ptr< Particle > > &secondaries, const RandomEngineAndDistribution &random) override
 Perform the interaction. More...
 
 ~Bremsstrahlung () override
 Default destructor. More...
 
- Public Member Functions inherited from fastsim::InteractionModel
const std::string getName ()
 Return (unique) name of this interaction. More...
 
 InteractionModel (std::string name)
 Constructor. More...
 
virtual void registerProducts (edm::ProducesCollector) const
 In case interaction produces and stores content in the event (e.g. TrackerSimHits). More...
 
virtual void storeProducts (edm::Event &iEvent)
 In case interaction produces and stores content in the event (e.g. TrackerSimHits). More...
 
virtual ~InteractionModel ()
 Default destructor. More...
 

Private Member Functions

math::XYZTLorentzVector brem (Particle &particle, double xmin, const RandomEngineAndDistribution &random) const
 Compute Brem photon energy and angles, if any. More...
 
double gbteth (const double ener, const double partm, const double efrac, const RandomEngineAndDistribution &random) const
 A universal angular distribution. More...
 

Private Attributes

double minPhotonEnergy_
 Cut on minimum energy of bremsstrahlung photons. More...
 
double minPhotonEnergyFraction_
 Cut on minimum fraction of particle's energy which has to be carried by photon. More...
 
double Z_
 Atomic number of material (usually silicon Z=14) More...
 

Detailed Description

Implementation of Bremsstrahlung from e+/e- in the tracker layers.

Computes the number, energy and angles of Bremsstrahlung photons emitted by electrons and positrons and modifies e+/e- particle accordingly.

Definition at line 33 of file Bremsstrahlung.cc.

Constructor & Destructor Documentation

◆ Bremsstrahlung()

fastsim::Bremsstrahlung::Bremsstrahlung ( const std::string &  name,
const edm::ParameterSet cfg 
)

Constructor.

Definition at line 82 of file Bremsstrahlung.cc.

84  // Set the minimal photon energy for a Brem from e+/-
85  minPhotonEnergy_ = cfg.getParameter<double>("minPhotonEnergy");
86  minPhotonEnergyFraction_ = cfg.getParameter<double>("minPhotonEnergyFraction");
87  // Material properties
88  Z_ = cfg.getParameter<double>("Z");
89 }

References looper::cfg, minPhotonEnergy_, minPhotonEnergyFraction_, and Z_.

◆ ~Bremsstrahlung()

fastsim::Bremsstrahlung::~Bremsstrahlung ( )
inlineoverride

Default destructor.

Definition at line 39 of file Bremsstrahlung.cc.

39 { ; };

Member Function Documentation

◆ brem()

math::XYZTLorentzVector fastsim::Bremsstrahlung::brem ( fastsim::Particle particle,
double  xmin,
const RandomEngineAndDistribution random 
) const
private

Compute Brem photon energy and angles, if any.

Parameters
particleThe particle that interacts with the matter.
xminMinimum fraction of the particle's energy that has to be converted to a photon.
randomThe Random Engine.
Returns
Momentum 4-vector of a bremsstrahlung photon.

Definition at line 161 of file Bremsstrahlung.cc.

163  {
164  double xp = 0;
165  double weight = 0.;
166 
167  do {
168  xp = xmin * std::exp(-std::log(xmin) * random.flatShoot());
169  weight = 1. - xp + 3. / 4. * xp * xp;
170  } while (weight < random.flatShoot());
171 
172  // Have photon energy. Now generate angles with respect to the z axis
173  // defined by the incoming particle's momentum.
174 
175  // Isotropic in phi
176  const double phi = random.flatShoot() * 2. * M_PI;
177  // theta from universal distribution
178  const double theta = gbteth(particle.momentum().E(), fastsim::Constants::eMass, xp, random) *
179  fastsim::Constants::eMass / particle.momentum().E();
180 
181  // Make momentum components
182  double stheta = std::sin(theta);
183  double ctheta = std::cos(theta);
184  double sphi = std::sin(phi);
185  double cphi = std::cos(phi);
186 
187  return xp * particle.momentum().E() * math::XYZTLorentzVector(stheta * cphi, stheta * sphi, ctheta, 1.);
188 }

References funct::cos(), fastsim::Constants::eMass, JetChargeProducer_cfi::exp, RandomEngineAndDistribution::flatShoot(), dqm-mbProfile::log, M_PI, fastsim::Particle::momentum(), funct::sin(), theta(), and TrackerOfflineValidation_Dqm_cff::xmin.

◆ gbteth()

double fastsim::Bremsstrahlung::gbteth ( const double  ener,
const double  partm,
const double  efrac,
const RandomEngineAndDistribution random 
) const
private

A universal angular distribution.

Parameters
ener
partm
efrac
randomThe Random Engine.
Returns
Theta from universal distribution

Definition at line 190 of file Bremsstrahlung.cc.

193  {
194  // Details on implementation here
195  // http://www.dnp.fmph.uniba.sk/cernlib/asdoc/geant_html3/node299.html#GBTETH
196  // http://svn.cern.ch/guest/AliRoot/tags/v3-07-03/GEANT321/gphys/gbteth.F
197 
198  const double alfa = 0.625;
199 
200  const double d = 0.13 * (0.8 + 1.3 / Z_) * (100.0 + (1.0 / ener)) * (1.0 + efrac);
201  const double w1 = 9.0 / (9.0 + d);
202  const double umax = ener * M_PI / partm;
203  double u;
204 
205  do {
206  double beta = (random.flatShoot() <= w1) ? alfa : 3.0 * alfa;
207  u = -std::log(random.flatShoot() * random.flatShoot()) / beta;
208  } while (u >= umax);
209 
210  return u;
211 }

References HLT_FULL_cff::beta, ztail::d, RandomEngineAndDistribution::flatShoot(), dqm-mbProfile::log, and M_PI.

◆ interact()

void fastsim::Bremsstrahlung::interact ( fastsim::Particle particle,
const SimplifiedGeometry layer,
std::vector< std::unique_ptr< Particle > > &  secondaries,
const RandomEngineAndDistribution random 
)
overridevirtual

Perform the interaction.

Parameters
particleThe particle that interacts with the matter.
layerThe detector layer that interacts with the particle.
secondariesParticles that are produced in the interaction (if any).
randomThe Random Engine.

Implements fastsim::InteractionModel.

Definition at line 91 of file Bremsstrahlung.cc.

94  {
95  // only consider electrons and positrons
96  if (std::abs(particle.pdgId()) != 11) {
97  return;
98  }
99 
100  double radLengths = layer.getThickness(particle.position(), particle.momentum());
101  //
102  // no material
103  //
104  if (radLengths < 1E-10) {
105  return;
106  }
107 
108  // Protection : Just stop the electron if more than 1 radiation lengths.
109  // This case corresponds to an electron entering the layer parallel to
110  // the layer axis - no reliable simulation can be done in that case...
111  if (radLengths > 4.) {
112  particle.momentum().SetXYZT(0., 0., 0., 0.);
113  return;
114  }
115 
116  // electron must have more energy than minimum photon energy
117  if (particle.momentum().E() - particle.momentum().mass() < minPhotonEnergy_) {
118  return;
119  }
120 
121  // Hard brem probability with a photon Energy above threshold.
122  double xmin = std::max(minPhotonEnergy_ / particle.momentum().E(), minPhotonEnergyFraction_);
123  if (xmin >= 1. || xmin <= 0.) {
124  return;
125  }
126 
127  // probability to radiate a photon
128  double bremProba =
129  radLengths * (4. / 3. * std::log(1. / xmin) - 4. / 3. * (1. - xmin) + 1. / 2. * (1. - xmin * xmin));
130 
131  // Number of photons to be radiated.
132  unsigned int nPhotons = random.poissonShoot(bremProba);
133  if (nPhotons == 0) {
134  return;
135  }
136 
137  // Needed to rotate photons to the lab frame
138  double theta = particle.momentum().Theta();
139  double phi = particle.momentum().Phi();
140 
141  // Calculate energy of these photons and add them to the event
142  for (unsigned int i = 0; i < nPhotons; ++i) {
143  // Throw momentum of the photon
144  math::XYZTLorentzVector photonMom = brem(particle, xmin, random);
145 
146  // Check that there is enough energy left.
147  if (particle.momentum().E() - particle.momentum().mass() < photonMom.E())
148  break;
149 
150  // Rotate to the lab frame
151  photonMom = ROOT::Math::RotationZ(phi) * (ROOT::Math::RotationY(theta) * photonMom);
152 
153  // Add a photon
154  secondaries.emplace_back(new fastsim::Particle(22, particle.position(), photonMom));
155 
156  // Update the original e+/-
157  particle.momentum() -= photonMom;
158  }
159 }

References funct::abs(), mps_fire::i, phase1PixelTopology::layer, dqm-mbProfile::log, SiStripPI::max, fastsim::Particle::momentum(), fastsim::Particle::pdgId(), RandomEngineAndDistribution::poissonShoot(), fastsim::Particle::position(), theta(), and TrackerOfflineValidation_Dqm_cff::xmin.

Member Data Documentation

◆ minPhotonEnergy_

double fastsim::Bremsstrahlung::minPhotonEnergy_
private

Cut on minimum energy of bremsstrahlung photons.

Definition at line 76 of file Bremsstrahlung.cc.

Referenced by Bremsstrahlung().

◆ minPhotonEnergyFraction_

double fastsim::Bremsstrahlung::minPhotonEnergyFraction_
private

Cut on minimum fraction of particle's energy which has to be carried by photon.

Definition at line 77 of file Bremsstrahlung.cc.

Referenced by Bremsstrahlung().

◆ Z_

double fastsim::Bremsstrahlung::Z_
private

Atomic number of material (usually silicon Z=14)

Definition at line 78 of file Bremsstrahlung.cc.

Referenced by Bremsstrahlung().

fastsim::Particle::momentum
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:143
mps_fire.i
i
Definition: mps_fire.py:428
fastsim::Bremsstrahlung::Z_
double Z_
Atomic number of material (usually silicon Z=14)
Definition: Bremsstrahlung.cc:78
HLT_FULL_cff.beta
beta
Definition: HLT_FULL_cff.py:8651
fastsim::Bremsstrahlung::minPhotonEnergy_
double minPhotonEnergy_
Cut on minimum energy of bremsstrahlung photons.
Definition: Bremsstrahlung.cc:76
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
fastsim::InteractionModel
Base class for any interaction model between a particle and a tracker layer.
Definition: InteractionModel.h:29
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
fastsim::Bremsstrahlung::brem
math::XYZTLorentzVector brem(Particle &particle, double xmin, const RandomEngineAndDistribution &random) const
Compute Brem photon energy and angles, if any.
Definition: Bremsstrahlung.cc:161
theta
Geom::Theta< T > theta() const
Definition: Basic3DVectorLD.h:150
fastsim::Particle::position
const math::XYZTLorentzVector & position() const
Return position of the particle.
Definition: Particle.h:140
fastsim::Bremsstrahlung::gbteth
double gbteth(const double ener, const double partm, const double efrac, const RandomEngineAndDistribution &random) const
A universal angular distribution.
Definition: Bremsstrahlung.cc:190
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
fastsim::Constants::eMass
static constexpr double eMass
Electron mass[GeV].
Definition: Constants.h:13
PVValHelper::phi
Definition: PVValidationHelpers.h:69
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
fastsim::Particle::pdgId
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:134
fastsim::Particle
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:16
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
RandomEngineAndDistribution::flatShoot
double flatShoot(double xmin=0.0, double xmax=1.0) const
Definition: RandomEngineAndDistribution.h:27
looper.cfg
cfg
Definition: looper.py:296
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
fastsim::Bremsstrahlung::minPhotonEnergyFraction_
double minPhotonEnergyFraction_
Cut on minimum fraction of particle's energy which has to be carried by photon.
Definition: Bremsstrahlung.cc:77
ztail.d
d
Definition: ztail.py:151
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TrackerOfflineValidation_Dqm_cff.xmin
xmin
Definition: TrackerOfflineValidation_Dqm_cff.py:10
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
RandomEngineAndDistribution::poissonShoot
unsigned int poissonShoot(double mean) const
Definition: RandomEngineAndDistribution.h:33
weight
Definition: weight.py:1