CMS 3D CMS Logo

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

Implementation of Bremsstrahlung from mu+/mu- in the tracker layers based on a Petrukhin Model (nuclear screening correction). More...

Inheritance diagram for fastsim::MuonBremsstrahlung:
fastsim::InteractionModel

Public Member Functions

void interact (Particle &particle, const SimplifiedGeometry &layer, std::vector< std::unique_ptr< Particle > > &secondaries, const RandomEngineAndDistribution &random) override
 Perform the interaction. More...
 
 MuonBremsstrahlung (const std::string &name, const edm::ParameterSet &cfg)
 Constructor. More...
 
 ~MuonBremsstrahlung () 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...
 

Static Private Member Functions

static double PetrukhinFunc (double *x, double *p)
 Petrukhin Function: Returns cross section using nuclear-electron screening correction from G4 style. More...
 

Private Attributes

double A_
 Atomic weight of material (usually silicon A=28.0855) More...
 
double density_
 Density of material (usually silicon rho=2.329) More...
 
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...
 
TF1 * Petrfunc
 The Petrukhin Function. More...
 
double radLenInCm_
 Radiation length of material (usually silicon X0=9.360) More...
 
double Z_
 Atomic number of material (usually silicon Z=14) More...
 

Detailed Description

Implementation of Bremsstrahlung from mu+/mu- in the tracker layers based on a Petrukhin Model (nuclear screening correction).

Computes the number, energy and angles of Bremsstrahlung photons emitted by muons and modifies mu+/mu- particle accordingly.

Definition at line 39 of file MuonBremsstrahlung.cc.

Constructor & Destructor Documentation

◆ MuonBremsstrahlung()

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

Constructor.

Definition at line 95 of file MuonBremsstrahlung.cc.

References A_, looper::cfg, density_, minPhotonEnergy_, minPhotonEnergyFraction_, radLenInCm_, and Z_.

97  // Set the minimal photon energy for a Brem from mu+/-
98  minPhotonEnergy_ = cfg.getParameter<double>("minPhotonEnergy");
99  minPhotonEnergyFraction_ = cfg.getParameter<double>("minPhotonEnergyFraction");
100  // Material properties
101  A_ = cfg.getParameter<double>("A");
102  Z_ = cfg.getParameter<double>("Z");
103  density_ = cfg.getParameter<double>("density");
104  radLenInCm_ = cfg.getParameter<double>("radLen");
105 }
double minPhotonEnergyFraction_
Cut on minimum fraction of particle&#39;s energy which has to be carried by photon.
double radLenInCm_
Radiation length of material (usually silicon X0=9.360)
Base class for any interaction model between a particle and a tracker layer.
double density_
Density of material (usually silicon rho=2.329)
double minPhotonEnergy_
Cut on minimum energy of bremsstrahlung photons.
double A_
Atomic weight of material (usually silicon A=28.0855)
double Z_
Atomic number of material (usually silicon Z=14)

◆ ~MuonBremsstrahlung()

fastsim::MuonBremsstrahlung::~MuonBremsstrahlung ( )
inlineoverride

Default destructor.

Definition at line 45 of file MuonBremsstrahlung.cc.

45 { ; };

Member Function Documentation

◆ brem()

math::XYZTLorentzVector fastsim::MuonBremsstrahlung::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 193 of file MuonBremsstrahlung.cc.

References funct::cos(), RandomEngineAndDistribution::flatShoot(), M_PI, fastsim::Particle::momentum(), fastsim::Constants::muMass, funct::sin(), and theta().

195  {
196  // This is a simple version of a Muon Brem using Petrukhin model.
197  // Ref: http://pdg.lbl.gov/2008/AtomicNuclearProperties/adndt.pdf
198  double xp = Petrfunc->GetRandom();
199 
200  // Have photon energy. Now generate angles with respect to the z axis
201  // defined by the incoming particle's momentum.
202 
203  // Isotropic in phi
204  const double phi = random.flatShoot() * 2. * M_PI;
205  // theta from universal distribution
206  const double theta = gbteth(particle.momentum().E(), fastsim::Constants::muMass, xp, random) *
207  fastsim::Constants::muMass / particle.momentum().E();
208 
209  // Make momentum components
210  double stheta = std::sin(theta);
211  double ctheta = std::cos(theta);
212  double sphi = std::sin(phi);
213  double cphi = std::cos(phi);
214 
215  return xp * particle.momentum().E() * math::XYZTLorentzVector(stheta * cphi, stheta * sphi, ctheta, 1.);
216 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
static double constexpr muMass
Muon mass [GeV].
Definition: Constants.h:14
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double gbteth(const double ener, const double partm, const double efrac, const RandomEngineAndDistribution &random) const
A universal angular distribution.
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
TF1 * Petrfunc
The Petrukhin Function.
#define M_PI
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:143
double flatShoot(double xmin=0.0, double xmax=1.0) const
Geom::Theta< T > theta() const

◆ gbteth()

double fastsim::MuonBremsstrahlung::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 218 of file MuonBremsstrahlung.cc.

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

221  {
222  // Details on implementation here
223  // http://www.dnp.fmph.uniba.sk/cernlib/asdoc/geant_html3/node299.html#GBTETH
224  // http://svn.cern.ch/guest/AliRoot/tags/v3-07-03/GEANT321/gphys/gbteth.F
225 
226  const double alfa = 0.625;
227 
228  const double d = 0.13 * (0.8 + 1.3 / Z_) * (100.0 + (1.0 / ener)) * (1.0 + efrac);
229  const double w1 = 9.0 / (9.0 + d);
230  const double umax = ener * M_PI / partm;
231  double u;
232 
233  do {
234  double beta = (random.flatShoot() <= w1) ? alfa : 3.0 * alfa;
235  u = -std::log(random.flatShoot() * random.flatShoot()) / beta;
236  } while (u >= umax);
237 
238  return u;
239 }
d
Definition: ztail.py:151
#define M_PI
double Z_
Atomic number of material (usually silicon Z=14)
double flatShoot(double xmin=0.0, double xmax=1.0) const

◆ interact()

void fastsim::MuonBremsstrahlung::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 107 of file MuonBremsstrahlung.cc.

References funct::abs(), HLT_2022v12_cff::distance, mps_fire::i, phase1PixelTopology::layer, SiStripPI::max, fastsim::Particle::momentum(), fastsim::Constants::NA, fastsim::Particle::pdgId(), PetrukhinFunc(), RandomEngineAndDistribution::poissonShoot(), fastsim::Particle::position(), theta(), TrackerOfflineValidation_Dqm_cff::xmax, and TrackerOfflineValidation_Dqm_cff::xmin.

110  {
111  // only consider muons
112  if (std::abs(particle.pdgId()) != 13) {
113  return;
114  }
115 
116  double radLengths = layer.getThickness(particle.position(), particle.momentum());
117  //
118  // no material
119  //
120  if (radLengths < 1E-10) {
121  return;
122  }
123 
124  // Protection : Just stop the electron if more than 1 radiation lengths.
125  // This case corresponds to an electron entering the layer parallel to
126  // the layer axis - no reliable simulation can be done in that case...
127  if (radLengths > 4.) {
128  particle.momentum().SetXYZT(0., 0., 0., 0.);
129  return;
130  }
131 
132  // muon must have more energy than minimum photon energy
133  if (particle.momentum().E() - particle.momentum().mass() < minPhotonEnergy_) {
134  return;
135  }
136 
137  // Min fraction of muon's energy transferred to the photon
138  double xmin = std::max(minPhotonEnergy_ / particle.momentum().E(), minPhotonEnergyFraction_);
139  // Hard brem probability with a photon Energy above threshold.
140  if (xmin >= 1. || xmin <= 0.) {
141  return;
142  }
143 
144  // Max fraction of muon's energy transferred to the photon
145  double xmax = 1.;
146 
147  // create TF1 using a free C function
148  Petrfunc = new TF1("Petrfunc", PetrukhinFunc, xmin, xmax, 3);
149  // Setting parameters
150  Petrfunc->SetParameters(particle.momentum().E(), A_, Z_);
151  // d = distance for several materials
152  // X0 = radLen
153  // d = radLengths * X0 (for tracker, yoke, ECAL and HCAL)
154  double distance = radLengths * radLenInCm_;
155 
156  // Integration
157  // Fixed previous version which used Petrfunc->Integral(0.,1.) -> does not make sense
158  double bremProba = density_ * distance * (fastsim::Constants::NA / A_) * (Petrfunc->Integral(xmin, xmax));
159  if (bremProba < 1E-10) {
160  return;
161  }
162 
163  // Number of photons to be radiated.
164  unsigned int nPhotons = random.poissonShoot(bremProba);
165  if (nPhotons == 0) {
166  return;
167  }
168 
169  //Rotate to the lab frame
170  double theta = particle.momentum().Theta();
171  double phi = particle.momentum().Phi();
172 
173  // Energy of these photons
174  for (unsigned int i = 0; i < nPhotons; ++i) {
175  // Throw momentum of the photon
176  math::XYZTLorentzVector photonMom = brem(particle, xmin, random);
177 
178  // Check that there is enough energy left.
179  if (particle.momentum().E() - particle.momentum().mass() < photonMom.E())
180  break;
181 
182  // Rotate to the lab frame
183  photonMom = ROOT::Math::RotationZ(phi) * (ROOT::Math::RotationY(theta) * photonMom);
184 
185  // Add a photon
186  secondaries.emplace_back(new fastsim::Particle(22, particle.position(), photonMom));
187 
188  // Update the original e+/-
189  particle.momentum() -= photonMom;
190  }
191 }
static double constexpr NA
Avogadro&#39;s number.
Definition: Constants.h:16
const math::XYZTLorentzVector & position() const
Return position of the particle.
Definition: Particle.h:140
double minPhotonEnergyFraction_
Cut on minimum fraction of particle&#39;s energy which has to be carried by photon.
double radLenInCm_
Radiation length of material (usually silicon X0=9.360)
constexpr std::array< uint8_t, layerIndexSize > layer
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
math::XYZTLorentzVector brem(Particle &particle, double xmin, const RandomEngineAndDistribution &random) const
Compute Brem photon energy and angles, if any.
double density_
Density of material (usually silicon rho=2.329)
double minPhotonEnergy_
Cut on minimum energy of bremsstrahlung photons.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TF1 * Petrfunc
The Petrukhin Function.
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:134
static double PetrukhinFunc(double *x, double *p)
Petrukhin Function: Returns cross section using nuclear-electron screening correction from G4 style...
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:143
unsigned int poissonShoot(double mean) const
double A_
Atomic weight of material (usually silicon A=28.0855)
double Z_
Atomic number of material (usually silicon Z=14)
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:16
Geom::Theta< T > theta() const

◆ PetrukhinFunc()

double fastsim::MuonBremsstrahlung::PetrukhinFunc ( double *  x,
double *  p 
)
staticprivate

Petrukhin Function: Returns cross section using nuclear-electron screening correction from G4 style.

Definition at line 241 of file MuonBremsstrahlung.cc.

References A, alpha, B, dumpMFGeometry_cfg::delta, triggerObjects_cff::Dn, f, AlCaHLTBitMon_ParallelJobs::p, funct::pow(), x, and BeamSpotPI::Z.

241  {
242  // Function independent variable
243  double nu = x[0]; // fraction of muon's energy transferred to the photon
244 
245  // Parameters
246  double E = p[0]; // Muon Energy (in GeV)
247  double A = p[1]; // Atomic weight
248  double Z = p[2]; // Atomic number
249 
250  /*
252  // Function of Muom Brem using nuclear screening correction
253  // Ref: http://pdg.lbl.gov/2008/AtomicNuclearProperties/adndt.pdf
254  // http://geant4.cern.ch/G4UsersDocuments/UsersGuides/PhysicsReferenceManual/html/node48.html
255 
256  // Physical constants
257  double B = 182.7;
258  double ee = sqrt(2.7181) ; // sqrt(e)
259  double ZZ= pow( Z,-1./3.); // Z^-1/3
261  double emass = 0.0005109990615; // electron mass (GeV/c^2)
262  double mumass = 0.105658367;//mu mass (GeV/c^2)
263 
264  double re = 2.817940285e-13;// Classical electron radius (Units: cm)
265  double alpha = 1./137.03599976; // fine structure constant
266  double Dn = 1.54* (pow(A,0.27));
267  double constant = pow((2.0 * Z * emass/mumass * re ),2.0);
269 
270  double delta = (mumass * mumass * nu) /(2.* E * (1.- nu));
271 
272  double Delta_n = TMath::Log(Dn / (1.+ delta *( Dn * ee -2.)/ mumass)); //nuclear screening correction
273 
274  double Phi = TMath::Log((B * mumass * ZZ / emass)/ (1.+ delta * ee * B * ZZ / emass)) - Delta_n;//phi(delta)
275 
276  // Diff. Cross Section for Muon Brem from a screened nuclear (Equation 16: REF: LBNL-44742)
277  double f = alpha * constant *(4./3.-4./3.*nu + nu*nu)*Phi/nu;
278  */
279 
281  // Function for Muon Brem Xsec from G4
283 
284  // Physical constants
285  double B = 183.;
286  double Bl = 1429.;
287  double ee = 1.64872; // sqrt(e)
288  double Z13 = pow(Z, -1. / 3.); // Z^-1/3
289  double Z23 = pow(Z, -2. / 3.); // Z^-2/3
290 
291  // Original values of paper
292  double emass = 0.0005109990615; // electron mass (GeV/c^2)
293  double mumass = 0.105658367; // muon mass (GeV/c^2)
294  // double re = 2.817940285e-13; // Classical electron radius (Units: cm)
295  double alpha = 0.00729735; // 1./137.03599976; // fine structure constant
296  double constant = 1.85736e-30; // pow( ( emass / mumass * re ) , 2.0);
297 
298  // Use nomenclature from reference -> Follow those formula step by step
299  if (nu * E >= E - mumass)
300  return 0;
301 
302  double Dn = 1.54 * (pow(A, 0.27));
303  double Dnl = pow(Dn, (1. - 1. / Z));
304 
305  double delta = (mumass * mumass * nu) / (2. * E * (1. - nu));
306 
307  double Phi_n = TMath::Log(B * Z13 * (mumass + delta * (Dnl * ee - 2)) / (Dnl * (emass + delta * ee * B * Z13)));
308  if (Phi_n < 0)
309  Phi_n = 0;
310 
311  double Phi_e =
312  TMath::Log((Bl * Z23 * mumass) / (1. + delta * mumass / (emass * emass * ee)) / (emass + delta * ee * Bl * Z23));
313  if (Phi_e < 0 || nu * E >= E / (1. + (mumass * mumass / (2. * emass * E))))
314  Phi_e = 0;
315 
316  // Diff. Cross Section for Muon Brem from G4 (without NA/A factor)
317  double f = 16. / 3. * alpha * constant * Z * (Z * Phi_n + Phi_e) * (1. / nu) * (1. - nu + 0.75 * nu * nu);
318 
319  return f;
320 }
float alpha
Definition: AMPTWrapper.h:105
Definition: APVGainStruct.h:7
double f[11][100]
float x
Definition: APVGainStruct.h:7
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

Member Data Documentation

◆ A_

double fastsim::MuonBremsstrahlung::A_
private

Atomic weight of material (usually silicon A=28.0855)

Definition at line 90 of file MuonBremsstrahlung.cc.

Referenced by MuonBremsstrahlung().

◆ density_

double fastsim::MuonBremsstrahlung::density_
private

Density of material (usually silicon rho=2.329)

Definition at line 88 of file MuonBremsstrahlung.cc.

Referenced by MuonBremsstrahlung().

◆ minPhotonEnergy_

double fastsim::MuonBremsstrahlung::minPhotonEnergy_
private

Cut on minimum energy of bremsstrahlung photons.

Definition at line 86 of file MuonBremsstrahlung.cc.

Referenced by MuonBremsstrahlung().

◆ minPhotonEnergyFraction_

double fastsim::MuonBremsstrahlung::minPhotonEnergyFraction_
private

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

Definition at line 87 of file MuonBremsstrahlung.cc.

Referenced by MuonBremsstrahlung().

◆ Petrfunc

TF1* fastsim::MuonBremsstrahlung::Petrfunc
private

The Petrukhin Function.

Definition at line 85 of file MuonBremsstrahlung.cc.

◆ radLenInCm_

double fastsim::MuonBremsstrahlung::radLenInCm_
private

Radiation length of material (usually silicon X0=9.360)

Definition at line 89 of file MuonBremsstrahlung.cc.

Referenced by MuonBremsstrahlung().

◆ Z_

double fastsim::MuonBremsstrahlung::Z_
private

Atomic number of material (usually silicon Z=14)

Definition at line 91 of file MuonBremsstrahlung.cc.

Referenced by MuonBremsstrahlung().