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::ProducerBase &producer) 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 36 of file Bremsstrahlung.cc.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 83 of file Bremsstrahlung.cc.

References edm::ParameterSet::getParameter(), minPhotonEnergy_, minPhotonEnergyFraction_, and Z_.

85 {
86  // Set the minimal photon energy for a Brem from e+/-
87  minPhotonEnergy_ = cfg.getParameter<double>("minPhotonEnergy");
88  minPhotonEnergyFraction_ = cfg.getParameter<double>("minPhotonEnergyFraction");
89  // Material properties
90  Z_ = cfg.getParameter<double>("Z");
91 }
T getParameter(std::string const &) const
double minPhotonEnergy_
Cut on minimum energy of bremsstrahlung photons.
double Z_
Atomic number of material (usually silicon Z=14)
Base class for any interaction model between a particle and a tracker layer.
double minPhotonEnergyFraction_
Cut on minimum fraction of particle&#39;s energy which has to be carried by photon.
fastsim::Bremsstrahlung::~Bremsstrahlung ( )
inlineoverride

Default destructor.

Definition at line 43 of file Bremsstrahlung.cc.

References brem(), gbteth(), interact(), random, and TrackerOfflineValidation_Dqm_cff::xmin.

43 {;};

Member Function Documentation

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 172 of file Bremsstrahlung.cc.

References funct::cos(), fastsim::Constants::eMass, JetChargeProducer_cfi::exp, RandomEngineAndDistribution::flatShoot(), gbteth(), cmsBatch::log, M_PI, fastsim::Particle::momentum(), random, funct::sin(), and theta().

Referenced by interact(), and ~Bremsstrahlung().

173 {
174  double xp=0;
175  double weight = 0.;
176 
177  do{
178  xp = xmin * std::exp ( -std::log(xmin) * random.flatShoot() );
179  weight = 1. - xp + 3./4.*xp*xp;
180  }while(weight < random.flatShoot());
181 
182 
183  // Have photon energy. Now generate angles with respect to the z axis
184  // defined by the incoming particle's momentum.
185 
186  // Isotropic in phi
187  const double phi = random.flatShoot() * 2. * M_PI;
188  // theta from universal distribution
189  const double theta = gbteth(particle.momentum().E(), fastsim::Constants::eMass, xp, random)
190  * fastsim::Constants::eMass / particle.momentum().E();
191 
192  // Make momentum components
193  double stheta = std::sin(theta);
194  double ctheta = std::cos(theta);
195  double sphi = std::sin(phi);
196  double cphi = std::cos(phi);
197 
198  return xp * particle.momentum().E() * math::XYZTLorentzVector(stheta*cphi, stheta*sphi, ctheta, 1.);
199 }
double flatShoot(double xmin=0.0, double xmax=1.0) const
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
static double constexpr eMass
Electron mass[GeV].
Definition: Constants.h:17
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
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, const RandomEngineAndDistribution &random) const
A universal angular distribution.
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:145
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 202 of file Bremsstrahlung.cc.

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

Referenced by brem(), and ~Bremsstrahlung().

206 {
207  // Details on implementation here
208  // http://www.dnp.fmph.uniba.sk/cernlib/asdoc/geant_html3/node299.html#GBTETH
209  // http://svn.cern.ch/guest/AliRoot/tags/v3-07-03/GEANT321/gphys/gbteth.F
210 
211  const double alfa = 0.625;
212 
213  const double d = 0.13*(0.8+1.3/Z_)*(100.0+(1.0/ener))*(1.0+efrac);
214  const double w1 = 9.0/(9.0+d);
215  const double umax = ener*M_PI/partm;
216  double u;
217 
218  do
219  {
220  double beta = (random.flatShoot()<=w1) ? alfa : 3.0*alfa;
221  u = -std::log(random.flatShoot()*random.flatShoot())/beta;
222  }while (u >= umax);
223 
224  return u;
225 }
double flatShoot(double xmin=0.0, double xmax=1.0) const
double Z_
Atomic number of material (usually silicon Z=14)
#define M_PI
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 94 of file Bremsstrahlung.cc.

References funct::abs(), brem(), fastsim::SimplifiedGeometry::getThickness(), mps_fire::i, cmsBatch::log, SiStripPI::max, minPhotonEnergy_, minPhotonEnergyFraction_, fastsim::Particle::momentum(), fastsim::Particle::pdgId(), RandomEngineAndDistribution::poissonShoot(), fastsim::Particle::position(), theta(), and TrackerOfflineValidation_Dqm_cff::xmin.

Referenced by ~Bremsstrahlung().

95 {
96  // only consider electrons and positrons
97  if(std::abs(particle.pdgId())!=11)
98  {
99  return;
100  }
101 
102  double radLengths = layer.getThickness(particle.position(),particle.momentum());
103  //
104  // no material
105  //
106  if(radLengths < 1E-10)
107  {
108  return;
109  }
110 
111  // Protection : Just stop the electron if more than 1 radiation lengths.
112  // This case corresponds to an electron entering the layer parallel to
113  // the layer axis - no reliable simulation can be done in that case...
114  if(radLengths > 4.)
115  {
116  particle.momentum().SetXYZT(0.,0.,0.,0.);
117  return;
118  }
119 
120  // electron must have more energy than minimum photon energy
121  if(particle.momentum().E() - particle.momentum().mass() < minPhotonEnergy_)
122  {
123  return;
124  }
125 
126  // Hard brem probability with a photon Energy above threshold.
128  if(xmin >=1. || xmin <=0.)
129  {
130  return;
131  }
132 
133  // probability to radiate a photon
134  double bremProba = radLengths * (4./3. * std::log(1./xmin)
135  - 4./3. * (1.-xmin)
136  + 1./2. * (1.-xmin*xmin));
137 
138 
139  // Number of photons to be radiated.
140  unsigned int nPhotons = random.poissonShoot(bremProba);
141  if(nPhotons == 0)
142  {
143  return;
144  }
145 
146  // Needed to rotate photons to the lab frame
147  double theta = particle.momentum().Theta();
148  double phi = particle.momentum().Phi();
149 
150  // Calculate energy of these photons and add them to the event
151  for(unsigned int i=0; i<nPhotons; ++i)
152  {
153  // Throw momentum of the photon
154  math::XYZTLorentzVector photonMom = brem(particle, xmin, random);
155 
156  // Check that there is enough energy left.
157  if(particle.momentum().E() - particle.momentum().mass() < photonMom.E()) break;
158 
159  // Rotate to the lab frame
160  photonMom = ROOT::Math::RotationZ(phi) * (ROOT::Math::RotationY(theta) * photonMom);
161 
162  // Add a photon
163  secondaries.emplace_back(new fastsim::Particle(22, particle.position(), photonMom));
164 
165  // Update the original e+/-
166  particle.momentum() -= photonMom;
167  }
168 }
const math::XYZTLorentzVector & position() const
Return position of the particle.
Definition: Particle.h:142
double minPhotonEnergy_
Cut on minimum energy of bremsstrahlung photons.
Geom::Theta< T > theta() const
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:136
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int poissonShoot(double mean) const
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:145
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:19
double minPhotonEnergyFraction_
Cut on minimum fraction of particle&#39;s energy which has to be carried by photon.
math::XYZTLorentzVector brem(Particle &particle, double xmin, const RandomEngineAndDistribution &random) const
Compute Brem photon energy and angles, if any.

Member Data Documentation

double fastsim::Bremsstrahlung::minPhotonEnergy_
private

Cut on minimum energy of bremsstrahlung photons.

Definition at line 77 of file Bremsstrahlung.cc.

Referenced by Bremsstrahlung(), and interact().

double fastsim::Bremsstrahlung::minPhotonEnergyFraction_
private

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

Definition at line 78 of file Bremsstrahlung.cc.

Referenced by Bremsstrahlung(), and interact().

double fastsim::Bremsstrahlung::Z_
private

Atomic number of material (usually silicon Z=14)

Definition at line 79 of file Bremsstrahlung.cc.

Referenced by Bremsstrahlung(), and gbteth().