CMS 3D CMS Logo

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

Implementation of multiple scattering in the tracker layers. More...

Inheritance diagram for fastsim::MultipleScattering:
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...
 
 MultipleScattering (const std::string &name, const edm::ParameterSet &cfg)
 Constructor. More...
 
 ~MultipleScattering () 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

XYZVector orthogonal (const XYZVector &aVector) const
 Return an orthogonal vector. More...
 

Private Attributes

double minPt_
 Cut on minimum pT of particle. More...
 
double radLenInCm_
 Radiation length of material (usually silicon X0=9.360) More...
 

Detailed Description

Implementation of multiple scattering in the tracker layers.

Computes the direction change by multiple scattering of a charged particle (assumes constant properties of material).

Definition at line 34 of file MultipleScattering.cc.

Constructor & Destructor Documentation

◆ MultipleScattering()

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

Constructor.

Definition at line 63 of file MultipleScattering.cc.

65  // Set the minimal pT for interaction
66  minPt_ = cfg.getParameter<double>("minPt");
67  // Material properties
68  radLenInCm_ = cfg.getParameter<double>("radLen");
69 }

References looper::cfg, minPt_, and radLenInCm_.

◆ ~MultipleScattering()

fastsim::MultipleScattering::~MultipleScattering ( )
inlineoverride

Default destructor.

Definition at line 40 of file MultipleScattering.cc.

40 { ; };

Member Function Documentation

◆ interact()

void fastsim::MultipleScattering::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 71 of file MultipleScattering.cc.

74  {
75  //
76  // only charged particles
77  //
78  if (particle.charge() == 0) {
79  return;
80  }
81 
82  double radLengths = layer.getThickness(particle.position(), particle.momentum());
83  //
84  // no material
85  //
86  if (radLengths < 1E-10) {
87  return;
88  }
89 
90  // particle must have minimum pT
91  if (particle.momentum().Pt() < minPt_) {
92  return;
93  }
94 
95  double p2 = particle.momentum().Vect().Mag2();
96  double m2 = particle.momentum().mass2();
97  double e = std::sqrt(p2 + m2);
98 
99  // This is p*beta
100  double pbeta = p2 / e;
101 
102  // Average multiple scattering angle from Moliere radius
103  // The sqrt(2) factor is because of the *space* angle
104  double theta0 = 0.0136 / pbeta * particle.charge() * std::sqrt(2. * radLengths) * (1. + 0.038 * std::log(radLengths));
105 
106  // Generate multiple scattering space angle perpendicular to the particle motion
107  double theta = random.gaussShoot(0., theta0);
108  // Plus a random rotation angle around the particle motion
109  double phi = 2. * M_PI * random.flatShoot();
110 
111  // The two rotations
112  ROOT::Math::AxisAngle rotation1(orthogonal(particle.momentum().Vect()), theta);
113  ROOT::Math::AxisAngle rotation2(particle.momentum().Vect(), phi);
114  // Rotate!
115  XYZVector rotated = rotation2((rotation1(particle.momentum().Vect())));
116  particle.momentum().SetXYZT(rotated.X(), rotated.Y(), rotated.Z(), particle.momentum().E());
117 
118  // Generate mutiple scattering displacements in cm (assuming the detectors
119  // are silicon only to determine the thickness) in the directions orthogonal
120  // to the vector normal to the surface
121  double xp = (cos(phi) * theta / 2. + random.gaussShoot(0., theta0) / sqrt(12.)) * radLengths * radLenInCm_;
122  double yp = (sin(phi) * theta / 2. + random.gaussShoot(0., theta0) / sqrt(12.)) * radLengths * radLenInCm_;
123 
124  // Determine a unitary vector tangent to the surface
125  XYZVector normal;
126  if (layer.isForward()) {
127  normal = XYZVector(0., 0., particle.momentum().Pz() > 0. ? 1. : -1.);
128  } else {
129  double norm = particle.position().Rho();
130  normal = XYZVector(particle.position().X() / norm, particle.position().Y() / norm, 0.);
131  }
132  XYZVector tangent = orthogonal(normal);
133  // The total displacement
134  XYZVector delta = xp * tangent + yp * normal.Cross(tangent);
135  // Translate!
136  particle.position().SetXYZT(particle.position().X() + delta.X(),
137  particle.position().Y() + delta.Y(),
138  particle.position().Z() + delta.Z(),
139  particle.position().T());
140 
141  // Make sure particle is still physically on the layer (particle moved on a tangent but a barrel layer is bent)
142  // Happens only in very very few cases
143  if (!layer.isForward() && std::abs(layer.getGeomProperty() - particle.position().Rho()) > 1e-2) {
144  // Scale positions so the radii agree
145  double scalePos = layer.getGeomProperty() / particle.position().Rho();
146  particle.position().SetXYZT(particle.position().X() * scalePos,
147  particle.position().Y() * scalePos,
148  particle.position().Z(),
149  particle.position().T());
150 
151  // Add a protection in case something goes wrong
152  if (std::abs(layer.getGeomProperty() - particle.position().Rho()) > 1e-2) {
153  throw cms::Exception("fastsim::MultipleScattering") << "particle no longer on layer's surface";
154  }
155  }
156 }

References funct::abs(), fastsim::Particle::charge(), funct::cos(), dumpMFGeometry_cfg::delta, MillePedeFileConverter_cfg::e, Exception, RandomEngineAndDistribution::flatShoot(), RandomEngineAndDistribution::gaussShoot(), phase1PixelTopology::layer, dqm-mbProfile::log, M_PI, fastsim::Particle::momentum(), p2, fastsim::Particle::position(), funct::sin(), mathSSE::sqrt(), and theta().

◆ orthogonal()

XYZVector fastsim::MultipleScattering::orthogonal ( const XYZVector aVector) const
private

Return an orthogonal vector.

Definition at line 158 of file MultipleScattering.cc.

158  {
159  double x = fabs(aVector.X());
160  double y = fabs(aVector.Y());
161  double z = fabs(aVector.Z());
162 
163  if (x < y)
164  return x < z ? XYZVector(0., aVector.Z(), -aVector.Y()) : XYZVector(aVector.Y(), -aVector.X(), 0.);
165  else
166  return y < z ? XYZVector(-aVector.Z(), 0., aVector.X()) : XYZVector(aVector.Y(), -aVector.X(), 0.);
167 }

Member Data Documentation

◆ minPt_

double fastsim::MultipleScattering::minPt_
private

Cut on minimum pT of particle.

Definition at line 58 of file MultipleScattering.cc.

Referenced by MultipleScattering().

◆ radLenInCm_

double fastsim::MultipleScattering::radLenInCm_
private

Radiation length of material (usually silicon X0=9.360)

Definition at line 59 of file MultipleScattering.cc.

Referenced by MultipleScattering().

fastsim::Particle::momentum
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:143
detailsBasic3DVector::z
float float float z
Definition: extBasic3DVector.h:14
XYZVector
math::XYZVector XYZVector
Definition: RawParticle.h:26
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
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
p2
double p2[4]
Definition: TauolaWrapper.h:90
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::MultipleScattering::orthogonal
XYZVector orthogonal(const XYZVector &aVector) const
Return an orthogonal vector.
Definition: MultipleScattering.cc:158
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
PVValHelper::phi
Definition: PVValidationHelpers.h:69
dumpMFGeometry_cfg.delta
delta
Definition: dumpMFGeometry_cfg.py:25
fastsim::Particle::charge
double charge() const
Return charge of the particle.
Definition: Particle.h:137
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:297
Exception
Definition: hltDiff.cc:245
genVertex_cff.x
x
Definition: genVertex_cff.py:12
detailsBasic3DVector::y
float float y
Definition: extBasic3DVector.h:14
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
XYZVector
math::XYZVector XYZVector
Definition: MultipleScattering.cc:27
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
RandomEngineAndDistribution::gaussShoot
double gaussShoot(double mean=0.0, double sigma=1.0) const
Definition: RandomEngineAndDistribution.h:29
fastsim::MultipleScattering::minPt_
double minPt_
Cut on minimum pT of particle.
Definition: MultipleScattering.cc:58
fastsim::MultipleScattering::radLenInCm_
double radLenInCm_
Radiation length of material (usually silicon X0=9.360)
Definition: MultipleScattering.cc:59
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37