CMS 3D CMS Logo

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

Computes the probability for photons to convert into an e+e- pair in the tracker layer. More...

Inheritance diagram for fastsim::PairProduction:
fastsim::InteractionModel

Public Member Functions

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

double gbteth (double ener, double partm, double efrac, const RandomEngineAndDistribution &random) const
 A universal angular distribution. More...
 

Private Attributes

double minPhotonEnergy_
 Cut on minimum energy of photons. More...
 
double Z_
 Atomic number of material (usually silicon Z=14) More...
 

Detailed Description

Computes the probability for photons to convert into an e+e- pair in the tracker layer.

In case, it returns a list of two Secondaries (e+ and e-).

Definition at line 34 of file PairProduction.cc.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 68 of file PairProduction.cc.

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

70 {
71  // Set the minimal photon energy for possible conversion
72  minPhotonEnergy_ = cfg.getParameter<double>("photonEnergyCut");
73  // Material properties
74  Z_ = cfg.getParameter<double>("Z");
75 }
T getParameter(std::string const &) const
double Z_
Atomic number of material (usually silicon Z=14)
Base class for any interaction model between a particle and a tracker layer.
double minPhotonEnergy_
Cut on minimum energy of photons.
fastsim::PairProduction::~PairProduction ( )
inlineoverride

Default destructor.

Definition at line 41 of file PairProduction.cc.

References gbteth(), interact(), and random.

41 {;};

Member Function Documentation

double fastsim::PairProduction::gbteth ( double  ener,
double  partm,
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 182 of file PairProduction.cc.

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

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

186 {
187  // Details on implementation here
188  // http://www.dnp.fmph.uniba.sk/cernlib/asdoc/geant_html3/node299.html#GBTETH
189  // http://svn.cern.ch/guest/AliRoot/tags/v3-07-03/GEANT321/gphys/gbteth.F
190 
191  const double alfa = 0.625;
192 
193  const double d = 0.13*(0.8+1.3/Z_)*(100.0+(1.0/ener))*(1.0+efrac);
194  const double w1 = 9.0/(9.0+d);
195  const double umax = ener*M_PI/partm;
196  double u;
197 
198  do
199  {
200  double beta = (random.flatShoot()<=w1) ? alfa : 3.0*alfa;
201  u = -std::log(random.flatShoot()*random.flatShoot())/beta;
202  }while (u>=umax);
203 
204  return u;
205 }
const double beta
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::PairProduction::interact ( fastsim::Particle particle,
const SimplifiedGeometry layer,
std::vector< std::unique_ptr< fastsim::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 78 of file PairProduction.cc.

References funct::cos(), fastsim::Constants::eMass, RandomEngineAndDistribution::flatShoot(), gbteth(), fastsim::SimplifiedGeometry::getThickness(), cmsBatch::log, M_PI, SiStripPI::max, minPhotonEnergy_, fastsim::Particle::momentum(), fastsim::Particle::pdgId(), fastsim::Particle::position(), funct::sin(), and mathSSE::sqrt().

Referenced by ~PairProduction().

79 {
80 
81  double eGamma = particle.momentum().e();
82  //
83  // only consider photons
84  //
85  if(particle.pdgId()!=22)
86  {
87  return;
88  }
89 
90  double radLengths = layer.getThickness(particle.position(),particle.momentum());
91  //
92  // no material
93  //
94  if(radLengths < 1E-10)
95  {
96  return;
97  }
98 
99  //
100  // The photon has enough energy to create a pair
101  //
102  if (eGamma < minPhotonEnergy_)
103  {
104  return;
105  }
106 
107  //
108  // Probability to convert is 7/9*(dx/X0)
109  //
110  if(-std::log(random.flatShoot()) > (7./9.)*radLengths)
111  {
112  return;
113  }
114 
115  double xe=0;
117  double xm=eMass/eGamma;
118  double weight = 0.;
119 
120  // Generate electron energy between emass and eGamma-emass
121  do{
122  xe = random.flatShoot()*(1.-2.*xm) + xm;
123  weight = 1. - 4./3.*xe*(1.-xe);
124  }while(weight < random.flatShoot());
125 
126  // the electron
127  double eElectron = xe * eGamma;
128  double tElectron = eElectron-eMass;
129  double pElectron = std::sqrt(std::max((eElectron+eMass)*tElectron,0.));
130 
131  // the positron
132  double ePositron = eGamma-eElectron;
133  double tPositron = ePositron-eMass;
134  double pPositron = std::sqrt((ePositron+eMass)*tPositron);
135 
136  // Generate angles
137  double phi = random.flatShoot()*2.*M_PI;
138  double sphi = std::sin(phi);
139  double cphi = std::cos(phi);
140 
141  double stheta1, stheta2, ctheta1, ctheta2;
142 
143  if(eElectron > ePositron){
144  double theta1 = gbteth(eElectron,eMass,xe,random)*eMass/eElectron;
145  stheta1 = std::sin(theta1);
146  ctheta1 = std::cos(theta1);
147  stheta2 = stheta1*pElectron/pPositron;
148  ctheta2 = std::sqrt(std::max(0.,1.0-(stheta2*stheta2)));
149  }else{
150  double theta2 = gbteth(ePositron,eMass,xe,random)*eMass/ePositron;
151  stheta2 = std::sin(theta2);
152  ctheta2 = std::cos(theta2);
153  stheta1 = stheta2*pPositron/pElectron;
154  ctheta1 = std::sqrt(std::max(0.,1.0-(stheta1*stheta1)));
155  }
156 
157  //Rotate to the lab frame
158  double thetaLab = particle.momentum().Theta();
159  double phiLab = particle.momentum().Phi();
160 
161  // Add a electron
162  secondaries.emplace_back(new fastsim::Particle(11,particle.position(),
163  math::XYZTLorentzVector(pElectron*stheta1*cphi,
164  pElectron*stheta1*sphi,
165  pElectron*ctheta1,
166  eElectron)));
167  secondaries.back()->momentum() = ROOT::Math::RotationZ(phiLab)*(ROOT::Math::RotationY(thetaLab)*secondaries.back()->momentum());
168 
169  // Add a positron
170  secondaries.emplace_back(new fastsim::Particle(-11,particle.position(),
171  math::XYZTLorentzVector(-pPositron*stheta2*cphi,
172  -pPositron*stheta2*sphi,
173  pPositron*ctheta2,
174  ePositron)));
175  secondaries.back()->momentum() = ROOT::Math::RotationZ(phiLab)*(ROOT::Math::RotationY(thetaLab)*secondaries.back()->momentum());
176 
177  // The photon converted
178  particle.momentum().SetXYZT(0.,0.,0.,0.);
179 }
const math::XYZTLorentzVector & position() const
Return position of the particle.
Definition: Particle.h:142
double flatShoot(double xmin=0.0, double xmax=1.0) const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Definition: weight.py:1
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:136
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
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
#define M_PI
double minPhotonEnergy_
Cut on minimum energy of photons.
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 gbteth(double ener, double partm, double efrac, const RandomEngineAndDistribution &random) const
A universal angular distribution.

Member Data Documentation

double fastsim::PairProduction::minPhotonEnergy_
private

Cut on minimum energy of photons.

Definition at line 63 of file PairProduction.cc.

Referenced by interact(), and PairProduction().

double fastsim::PairProduction::Z_
private

Atomic number of material (usually silicon Z=14)

Definition at line 64 of file PairProduction.cc.

Referenced by gbteth(), and PairProduction().