CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes
BaseParticlePropagator Class Reference

#include <BaseParticlePropagator.h>

Inheritance diagram for BaseParticlePropagator:
ParticlePropagator

Public Member Functions

bool backPropagate ()
 
 BaseParticlePropagator ()
 Default c'tor. More...
 
 BaseParticlePropagator (const RawParticle &myPart, double r, double z, double B)
 
 BaseParticlePropagator (const RawParticle &myPart, double r, double z, double B, double t)
 
double getMagneticField () const
 Get the magnetic field. More...
 
int getSuccess () const
 Has propagation been performed and was barrel or endcap reached ? More...
 
bool hasDecayed () const
 Has the particle decayed while propagated ? More...
 
double helixCentreDistToAxis () const
 The distance between the cylinder and the helix axes. More...
 
double helixCentreDistToAxis (double xC, double yC) const
 
double helixCentrePhi () const
 The azimuth if the vector joining the cylinder and the helix axes. More...
 
double helixCentrePhi (double xC, double yC) const
 
double helixCentreX () const
 The x coordinate of the helix axis. More...
 
double helixCentreX (double radius, double phi) const
 
double helixCentreY () const
 The y coordinate of the helix axis. More...
 
double helixCentreY (double radius, double phi) const
 
double helixRadius () const
 The helix Radius. More...
 
double helixRadius (double pT) const
 
double helixStartPhi () const
 The azimuth of the momentum at the vertex. More...
 
void increaseRCyl (double delta)
 Just an internal trick. More...
 
void init ()
 Initialize internal switches and quantities. More...
 
bool inside () const
 Is the vertex inside the cylinder ? (stricly inside : true) More...
 
bool inside (double rPos2) const
 
bool onBarrel () const
 Is the vertex already on the cylinder barrel ? More...
 
bool onBarrel (double rPos2) const
 
bool onEndcap () const
 Is the vertex already on the cylinder endcap ? More...
 
bool onEndcap (double rPos2) const
 
bool onFiducial () const
 Is the vertex on some material ? More...
 
bool onSurface () const
 Is the vertex already on the cylinder surface ? More...
 
bool onSurface (double rPos2) const
 
RawParticle const & particle () const
 The particle being propagated. More...
 
RawParticleparticle ()
 
bool propagate ()
 
BaseParticlePropagator propagated () const
 
bool propagateToBeamCylinder (const XYZTLorentzVector &v, double radius=0.)
 
bool propagateToClosestApproach (double x0=0., double y0=0, bool first=true)
 
bool propagateToEcal (bool first=true)
 
bool propagateToEcalEntrance (bool first=true)
 
bool propagateToHcalEntrance (bool first=true)
 
bool propagateToHcalExit (bool first=true)
 
bool propagateToHOLayer (bool first=true)
 
bool propagateToNominalVertex (const XYZTLorentzVector &hit2=XYZTLorentzVector(0., 0., 0., 0.))
 
bool propagateToPreshowerLayer1 (bool first=true)
 
bool propagateToPreshowerLayer2 (bool first=true)
 
bool propagateToVFcalEntrance (bool first=true)
 
void resetDebug ()
 
void setDebug ()
 Set the debug leve;. More...
 
void setMagneticField (double b)
 Set the magnetic field. More...
 
void setParticle (RawParticle const &iParticle)
 
void setPropagationConditions (double r, double z, bool firstLoop=true)
 Set the propagation characteristics (rCyl, zCyl and first loop only) More...
 
void setProperDecayTime (double t)
 Set the proper decay time. More...
 
double xyImpactParameter (double x0=0., double y0=0.) const
 Transverse impact parameter. More...
 
double zImpactParameter (double x0=0, double y0=0.) const
 Longitudinal impact parameter. More...
 

Protected Member Functions

double c_light () const
 The speed of light in mm/ns (!) without clhep (yeaaahhh!) More...
 

Protected Attributes

bool fiducial
 The particle traverses some real material. More...
 
int success
 0:propagation still be done, 1:reached 'barrel', 2:reached 'endcaps' More...
 

Private Attributes

double bField
 Magnetic field in the cylinder, oriented along the Z axis. More...
 
bool debug
 The debug level. More...
 
bool decayed
 The particle decayed while propagated ! More...
 
bool firstLoop
 Do only the first half-loop. More...
 
RawParticle particle_
 
int propDir
 The propagation direction. More...
 
double properDecayTime
 The proper decay time of the particle. More...
 
double properTime
 The proper time of the particle. More...
 
double rCyl
 Simulated particle that is to be resp has been propagated. More...
 
double rCyl2
 
double zCyl
 Half-height of the cylinder (centred at 0,0,0) to which propagation is done. More...
 

Detailed Description

This class is aimed at propagating charged and neutral particles (yet under the form of a RawParticle) from a given vertex to a cylinder, defined by a radius rCyl and a length 2*zCyl, centered in (0,0,0) and whose axis is parallel to the B field (B is oriented along z, by definition of the z axis).

The main method

bool propagate()

returns true if an intersection between the propagated RawParticle and the cylinder is found. The location of the intersection is given by the value of success:

The method

propagated()

returns a new RawParticle with the propagated coordinates, if overwriting is not considered an advantage by the user.

Particles can be propagated backward with the method

backPropagate()

Member functions

o propagateToPreshowerLayer1(),
o propagateToPreshowerLayer2(),
o propagateToEcalEntrance(),
o propagateToHcalEntrance(),
o propagateToHcalExit(),
o propagateToClosestApproach(),
only need a momentum, a vertex and an electric charge to operate. Radii, half-lengths and default B field (4T) are set therein by default.

As of today, no average loss of energy (dE/dx, Brems) is considered in the propagation. No uncertainty (multiple scattering, dE/dx, Brems) is yet implemented.

Author
Patrick Janot $Date : 19-Aug-2002, with subsequent modification for FAMOS
Version
15-Dec-2003

Definition at line 82 of file BaseParticlePropagator.h.

Constructor & Destructor Documentation

◆ BaseParticlePropagator() [1/3]

BaseParticlePropagator::BaseParticlePropagator ( )

Default c'tor.

Definition at line 5 of file BaseParticlePropagator.cc.

5 : rCyl(0.), rCyl2(0.), zCyl(0.), bField(0), properDecayTime(1E99) { ; }
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
double rCyl
Simulated particle that is to be resp has been propagated.
double properDecayTime
The proper decay time of the particle.
double bField
Magnetic field in the cylinder, oriented along the Z axis.

◆ BaseParticlePropagator() [2/3]

BaseParticlePropagator::BaseParticlePropagator ( const RawParticle myPart,
double  r,
double  z,
double  B 
)

Constructors taking as arguments a RawParticle, as well as the radius, half-height and magnetic field defining the cylinder for which propagation is to be performed, and optionally, the proper decay time

Definition at line 12 of file BaseParticlePropagator.cc.

References init().

13  : particle_(myPart), rCyl(R), rCyl2(R * R), zCyl(Z), bField(B), properDecayTime(1E99) {
14  init();
15 }
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
Definition: APVGainStruct.h:7
double rCyl
Simulated particle that is to be resp has been propagated.
void init()
Initialize internal switches and quantities.
double properDecayTime
The proper decay time of the particle.
double bField
Magnetic field in the cylinder, oriented along the Z axis.

◆ BaseParticlePropagator() [3/3]

BaseParticlePropagator::BaseParticlePropagator ( const RawParticle myPart,
double  r,
double  z,
double  B,
double  t 
)

Definition at line 7 of file BaseParticlePropagator.cc.

References init().

8  : particle_(myPart), rCyl(R), rCyl2(R * R), zCyl(Z), bField(B), properDecayTime(t) {
9  init();
10 }
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
Definition: APVGainStruct.h:7
double rCyl
Simulated particle that is to be resp has been propagated.
void init()
Initialize internal switches and quantities.
double properDecayTime
The proper decay time of the particle.
double bField
Magnetic field in the cylinder, oriented along the Z axis.

Member Function Documentation

◆ backPropagate()

bool BaseParticlePropagator::backPropagate ( )

Update the current instance, after the back-propagation of the particle to the surface of the cylinder

Definition at line 291 of file BaseParticlePropagator.cc.

References RawParticle::charge(), fileCollector::done, RawParticle::E(), particle_, propagate(), propDir, RawParticle::Px(), RawParticle::Py(), RawParticle::Pz(), RawParticle::setCharge(), and RawParticle::setMomentum().

Referenced by propagated(), and propagateToClosestApproach().

291  {
292  // Backpropagate
295  propDir = -1;
296  bool done = propagate();
299  propDir = +1;
300 
301  return done;
302 }
void setMomentum(const XYZTLorentzVector &vtx)
set the momentum
Definition: RawParticle.h:328
void setCharge(float q)
set the MEASURED charge
Definition: RawParticle.cc:32
double Pz() const
z of the momentum
Definition: RawParticle.h:303
double E() const
energy of the momentum
Definition: RawParticle.h:306
double charge() const
get the MEASURED charge
Definition: RawParticle.h:294
double Py() const
y of the momentum
Definition: RawParticle.h:300
int propDir
The propagation direction.
double Px() const
x of the momentum
Definition: RawParticle.h:297

◆ c_light()

double BaseParticlePropagator::c_light ( ) const
inlineprotected

The speed of light in mm/ns (!) without clhep (yeaaahhh!)

Definition at line 150 of file BaseParticlePropagator.h.

Referenced by helixRadius(), and propagateToBeamCylinder().

150 { return 299.792458; }

◆ getMagneticField()

double BaseParticlePropagator::getMagneticField ( ) const
inline

Get the magnetic field.

Definition at line 302 of file BaseParticlePropagator.h.

References bField.

Referenced by TrajectoryManager::createPSimHits(), and ConvBremSeedProducer::produce().

302 { return bField; }
double bField
Magnetic field in the cylinder, oriented along the Z axis.

◆ getSuccess()

int BaseParticlePropagator::getSuccess ( ) const
inline

Has propagation been performed and was barrel or endcap reached ?

Definition at line 296 of file BaseParticlePropagator.h.

References success.

Referenced by PFTrackTransformer::addPoints(), PFTrackTransformer::addPointsAndBrems(), reco::tau::atECALEntrance(), lowptgsfeleid::features_V0(), FBaseSimEvent::fill(), ConvBremSeedProducer::GoodCluster(), ConvBremSeedProducer::produce(), TrajectoryManager::propagateToCalorimeters(), and TrajectoryManager::propagateToLayer().

296 { return success; }
int success
0:propagation still be done, 1:reached &#39;barrel&#39;, 2:reached &#39;endcaps&#39;

◆ hasDecayed()

bool BaseParticlePropagator::hasDecayed ( ) const
inline

Has the particle decayed while propagated ?

Definition at line 293 of file BaseParticlePropagator.h.

References decayed.

Referenced by ParticlePropagator::propagateToBoundSurface(), and TrajectoryManager::propagateToCalorimeters().

293 { return decayed; }
bool decayed
The particle decayed while propagated !

◆ helixCentreDistToAxis() [1/2]

double BaseParticlePropagator::helixCentreDistToAxis ( ) const
inline

The distance between the cylinder and the helix axes.

Definition at line 235 of file BaseParticlePropagator.h.

References helixCentreX(), helixCentreY(), and mathSSE::sqrt().

Referenced by propagate(), propagateToBeamCylinder(), propagateToClosestApproach(), and xyImpactParameter().

235  {
236  // The distance between the cylinder and the helix axes
237  double xC = helixCentreX();
238  double yC = helixCentreY();
239  return std::sqrt(xC * xC + yC * yC);
240  }
T sqrt(T t)
Definition: SSEVec.h:19
double helixCentreX() const
The x coordinate of the helix axis.
double helixCentreY() const
The y coordinate of the helix axis.

◆ helixCentreDistToAxis() [2/2]

double BaseParticlePropagator::helixCentreDistToAxis ( double  xC,
double  yC 
) const
inline

Definition at line 242 of file BaseParticlePropagator.h.

References mathSSE::sqrt().

242  {
243  // Faster version of helixCentreDistToAxis
244  return std::sqrt(xC * xC + yC * yC);
245  }
T sqrt(T t)
Definition: SSEVec.h:19

◆ helixCentrePhi() [1/2]

double BaseParticlePropagator::helixCentrePhi ( ) const
inline

The azimuth if the vector joining the cylinder and the helix axes.

Definition at line 248 of file BaseParticlePropagator.h.

References helixCentreX(), and helixCentreY().

Referenced by propagate().

248  {
249  // The azimuth if the vector joining the cylinder and the helix axes
250  double xC = helixCentreX();
251  double yC = helixCentreY();
252  return xC == 0.0 && yC == 0.0 ? 0.0 : std::atan2(yC, xC);
253  }
double helixCentreX() const
The x coordinate of the helix axis.
double helixCentreY() const
The y coordinate of the helix axis.

◆ helixCentrePhi() [2/2]

double BaseParticlePropagator::helixCentrePhi ( double  xC,
double  yC 
) const
inline

Definition at line 255 of file BaseParticlePropagator.h.

255  {
256  // Faster version of helixCentrePhi()
257  return xC == 0.0 && yC == 0.0 ? 0.0 : std::atan2(yC, xC);
258  }

◆ helixCentreX() [1/2]

double BaseParticlePropagator::helixCentreX ( ) const
inline

The x coordinate of the helix axis.

Definition at line 213 of file BaseParticlePropagator.h.

References helixRadius(), helixStartPhi(), particle_, funct::sin(), and RawParticle::X().

Referenced by helixCentreDistToAxis(), helixCentrePhi(), propagate(), propagateToBeamCylinder(), propagateToClosestApproach(), and xyImpactParameter().

213  {
214  // The x coordinate of the helix axis
215  return particle_.X() - helixRadius() * std::sin(helixStartPhi());
216  }
double helixStartPhi() const
The azimuth of the momentum at the vertex.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double helixRadius() const
The helix Radius.
double X() const
x of vertex
Definition: RawParticle.h:286

◆ helixCentreX() [2/2]

double BaseParticlePropagator::helixCentreX ( double  radius,
double  phi 
) const
inline

Definition at line 218 of file BaseParticlePropagator.h.

References particle_, phi, CosmicsPD_Skims::radius, funct::sin(), and RawParticle::X().

218  {
219  // Fast version of helixCentreX()
220  return particle_.X() - radius * std::sin(phi);
221  }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double X() const
x of vertex
Definition: RawParticle.h:286

◆ helixCentreY() [1/2]

double BaseParticlePropagator::helixCentreY ( ) const
inline

The y coordinate of the helix axis.

Definition at line 224 of file BaseParticlePropagator.h.

References funct::cos(), helixRadius(), helixStartPhi(), particle_, and RawParticle::Y().

Referenced by helixCentreDistToAxis(), helixCentrePhi(), propagate(), propagateToBeamCylinder(), propagateToClosestApproach(), and xyImpactParameter().

224  {
225  // The y coordinate of the helix axis
226  return particle_.Y() + helixRadius() * std::cos(helixStartPhi());
227  }
double helixStartPhi() const
The azimuth of the momentum at the vertex.
double helixRadius() const
The helix Radius.
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double Y() const
y of vertex
Definition: RawParticle.h:287

◆ helixCentreY() [2/2]

double BaseParticlePropagator::helixCentreY ( double  radius,
double  phi 
) const
inline

Definition at line 229 of file BaseParticlePropagator.h.

References funct::cos(), particle_, phi, CosmicsPD_Skims::radius, and RawParticle::Y().

229  {
230  // Fast version of helixCentreX()
231  return particle_.Y() + radius * std::cos(phi);
232  }
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double Y() const
y of vertex
Definition: RawParticle.h:287

◆ helixRadius() [1/2]

double BaseParticlePropagator::helixRadius ( ) const
inline

The helix Radius.

Definition at line 189 of file BaseParticlePropagator.h.

References bField, c_light(), RawParticle::charge(), MillePedeFileConverter_cfg::e, particle_, and RawParticle::Pt().

Referenced by helixCentreX(), helixCentreY(), propagate(), propagateToClosestApproach(), propagateToNominalVertex(), and xyImpactParameter().

189  {
190  // The helix Radius
191  //
192  // The helix' Radius sign accounts for the orientation of the magnetic field
193  // (+ = along z axis) and the sign of the electric charge of the particle.
194  // It signs the rotation of the (charged) particle around the z axis:
195  // Positive means anti-clockwise, negative means clockwise rotation.
196  //
197  // The radius is returned in cm to match the units in RawParticle.
198  return particle_.charge() == 0 ? 0.0 : -particle_.Pt() / (c_light() * 1e-5 * bField * particle_.charge());
199  }
double c_light() const
The speed of light in mm/ns (!) without clhep (yeaaahhh!)
double charge() const
get the MEASURED charge
Definition: RawParticle.h:294
double Pt() const
transverse momentum
Definition: RawParticle.h:308
double bField
Magnetic field in the cylinder, oriented along the Z axis.

◆ helixRadius() [2/2]

double BaseParticlePropagator::helixRadius ( double  pT) const
inline

Definition at line 201 of file BaseParticlePropagator.h.

References bField, c_light(), RawParticle::charge(), MillePedeFileConverter_cfg::e, particle_, and pv::pT.

201  {
202  // a faster version of helixRadius, once Perp() has been computed
203  return particle_.charge() == 0 ? 0.0 : -pT / (c_light() * 1e-5 * bField * particle_.charge());
204  }
double c_light() const
The speed of light in mm/ns (!) without clhep (yeaaahhh!)
double charge() const
get the MEASURED charge
Definition: RawParticle.h:294
double bField
Magnetic field in the cylinder, oriented along the Z axis.

◆ helixStartPhi()

double BaseParticlePropagator::helixStartPhi ( ) const
inline

The azimuth of the momentum at the vertex.

Definition at line 207 of file BaseParticlePropagator.h.

References particle_, RawParticle::Px(), and RawParticle::Py().

Referenced by helixCentreX(), helixCentreY(), propagate(), propagateToClosestApproach(), and xyImpactParameter().

207  {
208  // The azimuth of the momentum at the vertex
209  return particle_.Px() == 0.0 && particle_.Py() == 0.0 ? 0.0 : std::atan2(particle_.Py(), particle_.Px());
210  }
double Py() const
y of the momentum
Definition: RawParticle.h:300
double Px() const
x of the momentum
Definition: RawParticle.h:297

◆ increaseRCyl()

void BaseParticlePropagator::increaseRCyl ( double  delta)
inline

Just an internal trick.

Definition at line 172 of file BaseParticlePropagator.h.

References dumpMFGeometry_cfg::delta, rCyl, and rCyl2.

172  {
173  rCyl = rCyl + delta;
174  rCyl2 = rCyl * rCyl;
175  }
double rCyl
Simulated particle that is to be resp has been propagated.

◆ init()

void BaseParticlePropagator::init ( void  )

Initialize internal switches and quantities.

Definition at line 17 of file BaseParticlePropagator.cc.

References debug, decayed, fiducial, firstLoop, propDir, properTime, and success.

Referenced by BaseParticlePropagator().

17  {
18  // The status of the propagation
19  success = 0;
20  // Propagate only for the first half-loop
21  firstLoop = true;
22  //
23  fiducial = true;
24  // The particle has not yet decayed
25  decayed = false;
26  // The proper time is set to zero
27  properTime = 0.;
28  // The propagation direction is along the momentum
29  propDir = 1;
30  //
31  debug = false;
32 }
bool firstLoop
Do only the first half-loop.
int success
0:propagation still be done, 1:reached &#39;barrel&#39;, 2:reached &#39;endcaps&#39;
double properTime
The proper time of the particle.
int propDir
The propagation direction.
bool decayed
The particle decayed while propagated !
bool debug
The debug level.
bool fiducial
The particle traverses some real material.

◆ inside() [1/2]

bool BaseParticlePropagator::inside ( ) const
inline

Is the vertex inside the cylinder ? (stricly inside : true)

Definition at line 261 of file BaseParticlePropagator.h.

References particle_, RawParticle::R2(), rCyl, rCyl2, RawParticle::Z(), and zCyl.

Referenced by propagate().

261  {
262  return (particle_.R2() < rCyl2 - 0.00001 * rCyl && fabs(particle_.Z()) < zCyl - 0.00001);
263  }
double R2() const
vertex radius**2
Definition: RawParticle.h:291
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
double rCyl
Simulated particle that is to be resp has been propagated.
double Z() const
z of vertex
Definition: RawParticle.h:288

◆ inside() [2/2]

bool BaseParticlePropagator::inside ( double  rPos2) const
inline

Definition at line 265 of file BaseParticlePropagator.h.

References particle_, rCyl, rCyl2, RawParticle::Z(), and zCyl.

265  {
266  return (rPos2 < rCyl2 - 0.00001 * rCyl && fabs(particle_.Z()) < zCyl - 0.00001);
267  }
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
double rCyl
Simulated particle that is to be resp has been propagated.
double Z() const
z of vertex
Definition: RawParticle.h:288

◆ onBarrel() [1/2]

bool BaseParticlePropagator::onBarrel ( ) const
inline

Is the vertex already on the cylinder barrel ?

Definition at line 275 of file BaseParticlePropagator.h.

References particle_, RawParticle::R2(), rCyl, rCyl2, RawParticle::Z(), and zCyl.

Referenced by onSurface(), and propagate().

275  {
276  double rPos2 = particle_.R2();
277  return (fabs(rPos2 - rCyl2) < 0.00001 * rCyl && fabs(particle_.Z()) <= zCyl);
278  }
double R2() const
vertex radius**2
Definition: RawParticle.h:291
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
double rCyl
Simulated particle that is to be resp has been propagated.
double Z() const
z of vertex
Definition: RawParticle.h:288

◆ onBarrel() [2/2]

bool BaseParticlePropagator::onBarrel ( double  rPos2) const
inline

Definition at line 280 of file BaseParticlePropagator.h.

References particle_, rCyl, rCyl2, RawParticle::Z(), and zCyl.

280  {
281  return (fabs(rPos2 - rCyl2) < 0.00001 * rCyl && fabs(particle_.Z()) <= zCyl);
282  }
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
double rCyl
Simulated particle that is to be resp has been propagated.
double Z() const
z of vertex
Definition: RawParticle.h:288

◆ onEndcap() [1/2]

bool BaseParticlePropagator::onEndcap ( ) const
inline

Is the vertex already on the cylinder endcap ?

Definition at line 285 of file BaseParticlePropagator.h.

References particle_, RawParticle::R2(), rCyl2, RawParticle::Z(), and zCyl.

Referenced by onSurface(), and propagate().

285 { return (fabs(fabs(particle_.Z()) - zCyl) < 0.00001 && particle_.R2() <= rCyl2); }
double R2() const
vertex radius**2
Definition: RawParticle.h:291
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
double Z() const
z of vertex
Definition: RawParticle.h:288

◆ onEndcap() [2/2]

bool BaseParticlePropagator::onEndcap ( double  rPos2) const
inline

Definition at line 287 of file BaseParticlePropagator.h.

References particle_, rCyl2, RawParticle::Z(), and zCyl.

287 { return (fabs(fabs(particle_.Z()) - zCyl) < 0.00001 && rPos2 <= rCyl2); }
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
double Z() const
z of vertex
Definition: RawParticle.h:288

◆ onFiducial()

bool BaseParticlePropagator::onFiducial ( ) const
inline

Is the vertex on some material ?

Definition at line 290 of file BaseParticlePropagator.h.

References fiducial.

Referenced by TrajectoryManager::propagateToLayer().

290 { return fiducial; }
bool fiducial
The particle traverses some real material.

◆ onSurface() [1/2]

bool BaseParticlePropagator::onSurface ( ) const
inline

Is the vertex already on the cylinder surface ?

Definition at line 270 of file BaseParticlePropagator.h.

References onBarrel(), and onEndcap().

Referenced by propagate().

270 { return (onBarrel() || onEndcap()); }
bool onEndcap() const
Is the vertex already on the cylinder endcap ?
bool onBarrel() const
Is the vertex already on the cylinder barrel ?

◆ onSurface() [2/2]

bool BaseParticlePropagator::onSurface ( double  rPos2) const
inline

Definition at line 272 of file BaseParticlePropagator.h.

References onBarrel(), and onEndcap().

272 { return (onBarrel(rPos2) || onEndcap(rPos2)); }
bool onEndcap() const
Is the vertex already on the cylinder endcap ?
bool onBarrel() const
Is the vertex already on the cylinder barrel ?

◆ particle() [1/2]

RawParticle const& BaseParticlePropagator::particle ( ) const
inline

◆ particle() [2/2]

RawParticle& BaseParticlePropagator::particle ( )
inline

Definition at line 165 of file BaseParticlePropagator.h.

References particle_.

165 { return particle_; }

◆ propagate()

bool BaseParticlePropagator::propagate ( )

Update the current instance, after the propagation of the particle to the surface of the cylinder

Definition at line 34 of file BaseParticlePropagator.cc.

References b2, bField, RawParticle::charge(), funct::cos(), decayed, dumpMFGeometry_cfg::delta, SiPixelRawToDigiRegional_cfi::deltaPhi, RawParticle::E(), firstLoop, helixCentreDistToAxis(), helixCentrePhi(), helixCentreX(), helixCentreY(), helixRadius(), helixStartPhi(), inside(), M_PI, RawParticle::mass(), SiStripPI::min, onBarrel(), onEndcap(), onSurface(), particle_, RawParticle::Perp2(), propDir, properDecayTime, properTime, pv::pT, RawParticle::Pt(), HLT_2022v15_cff::pT2, RawParticle::Px(), RawParticle::Py(), RawParticle::Pz(), SiStripMonitorCluster_cfi::q0, RawParticle::R2(), CosmicsPD_Skims::radius, rCyl, rCyl2, RawParticle::setMomentum(), RawParticle::setVertex(), funct::sin(), mathSSE::sqrt(), success, RawParticle::T(), RawParticle::X(), RawParticle::Y(), RawParticle::Z(), and zCyl.

Referenced by PFTrackTransformer::addPoints(), PFTrackTransformer::addPointsAndBrems(), backPropagate(), ConvBremSeedProducer::produce(), propagated(), ParticlePropagator::propagateToBoundSurface(), propagateToEcal(), propagateToEcalEntrance(), propagateToHcalEntrance(), propagateToHcalExit(), propagateToHOLayer(), propagateToPreshowerLayer1(), propagateToPreshowerLayer2(), and propagateToVFcalEntrance().

34  {
35  //
36  // Check that the particle is not already on the cylinder surface
37  //
38  double rPos2 = particle_.R2();
39 
40  if (onBarrel(rPos2)) {
41  success = 1;
42  return true;
43  }
44  //
45  if (onEndcap(rPos2)) {
46  success = 2;
47  return true;
48  }
49  //
50  // Treat neutral particles (or no magnetic field) first
51  //
52  if (fabs(particle_.charge()) < 1E-12 || bField < 1E-12) {
53  //
54  // First check that the particle crosses the cylinder
55  //
56  double pT2 = particle_.Perp2();
57  // double pT2 = pT*pT;
58  // double b2 = rCyl * pT;
59  double b2 = rCyl2 * pT2;
60  double ac = fabs(particle_.X() * particle_.Py() - particle_.Y() * particle_.Px());
61  double ac2 = ac * ac;
62  //
63  // The particle never crosses (or never crossed) the cylinder
64  //
65  // if ( ac > b2 ) return false;
66  if (ac2 > b2)
67  return false;
68 
69  // double delta = std::sqrt(b2*b2 - ac*ac);
70  double delta = std::sqrt(b2 - ac2);
71  double tplus = -(particle_.X() * particle_.Px() + particle_.Y() * particle_.Py()) + delta;
72  double tminus = -(particle_.X() * particle_.Px() + particle_.Y() * particle_.Py()) - delta;
73  //
74  // Find the first (time-wise) intersection point with the cylinder
75  //
76 
77  double solution = tminus < 0 ? tplus / pT2 : tminus / pT2;
78  if (solution < 0.)
79  return false;
80  //
81  // Check that the particle does (did) not exit from the endcaps first.
82  //
83  double zProp = particle_.Z() + particle_.Pz() * solution;
84  if (fabs(zProp) > zCyl) {
85  tplus = (zCyl - particle_.Z()) / particle_.Pz();
86  tminus = (-zCyl - particle_.Z()) / particle_.Pz();
87  solution = tminus < 0 ? tplus : tminus;
88  if (solution < 0.)
89  return false;
90  success = 2;
91  } else {
92  success = 1;
93  }
94 
95  //
96  // Check the decay time
97  //
98  double delTime = propDir * particle_.mass() * solution;
99  double factor = 1.;
100  properTime += delTime;
101  if (properTime > properDecayTime) {
102  factor = 1. - (properTime - properDecayTime) / delTime;
104  decayed = true;
105  }
106 
107  //
108  // Compute the coordinates of the RawParticle after propagation
109  //
110  double xProp = particle_.X() + particle_.Px() * solution * factor;
111  double yProp = particle_.Y() + particle_.Py() * solution * factor;
112  zProp = particle_.Z() + particle_.Pz() * solution * factor;
113  double tProp = particle_.T() + particle_.E() * solution * factor;
114 
115  //
116  // Last check : Although propagated to the endcaps, R could still be
117  // larger than rCyl because the cylinder was too short...
118  //
119  if (success == 2 && xProp * xProp + yProp * yProp > rCyl2) {
120  success = -1;
121  return true;
122  }
123 
124  //
125  // ... and update the particle with its propagated coordinates
126  //
127  particle_.setVertex(XYZTLorentzVector(xProp, yProp, zProp, tProp));
128 
129  return true;
130  }
131  //
132  // Treat charged particles with magnetic field next.
133  //
134  else {
135  //
136  // First check that the particle can cross the cylinder
137  //
138  double pT = particle_.Pt();
139  double radius = helixRadius(pT);
140  double phi0 = helixStartPhi();
141  double xC = helixCentreX(radius, phi0);
142  double yC = helixCentreY(radius, phi0);
143  double dist = helixCentreDistToAxis(xC, yC);
144  //
145  // The helix either is outside or includes the cylinder -> no intersection
146  //
147  if (fabs(fabs(radius) - dist) > rCyl)
148  return false;
149  //
150  // The particle is already away from the endcaps, and will never come back
151  // Could be back-propagated, so warn the user that it is so.
152  //
153  if (particle_.Z() * particle_.Pz() > zCyl * fabs(particle_.Pz())) {
154  success = -2;
155  return true;
156  }
157 
158  double pZ = particle_.Pz();
159  double phiProp, zProp;
160  //
161  // Does the helix cross the cylinder barrel ? If not, do the first half-loop
162  //
163  double rProp = std::min(fabs(radius) + dist - 0.000001, rCyl);
164 
165  // Check for rounding errors in the ArcSin.
166  double sinPhiProp = (rProp * rProp - radius * radius - dist * dist) / (2. * dist * radius);
167 
168  //
169  double deltaPhi = 1E99;
170 
171  // Taylor development up to third order for large momenta
172  if (1. - fabs(sinPhiProp) < 1E-9) {
173  double cphi0 = std::cos(phi0);
174  double sphi0 = std::sin(phi0);
175  double r0 = (particle_.X() * cphi0 + particle_.Y() * sphi0) / radius;
176  double q0 = (particle_.X() * sphi0 - particle_.Y() * cphi0) / radius;
177  double rcyl2 = (rCyl2 - particle_.X() * particle_.X() - particle_.Y() * particle_.Y()) / (radius * radius);
178  double delta = r0 * r0 + rcyl2 * (1. - q0);
179 
180  // This is a solution of a second order equation, assuming phi = phi0 + epsilon
181  // and epsilon is small.
182  deltaPhi = radius > 0 ? (-r0 + std::sqrt(delta)) / (1. - q0) : (-r0 - std::sqrt(delta)) / (1. - q0);
183  }
184 
185  // Use complete calculation otherwise, or if the delta phi is large anyway.
186  if (fabs(deltaPhi) > 1E-3) {
187  // phiProp = fabs(sinPhiProp) < 0.99999999 ?
188  // std::asin(sinPhiProp) : M_PI/2. * sinPhiProp;
189  phiProp = std::asin(sinPhiProp);
190 
191  //
192  // Compute phi after propagation (two solutions, but the asin definition
193  // (between -pi/2 and +pi/2) takes care of the wrong one.
194  //
195  phiProp = helixCentrePhi(xC, yC) + (inside(rPos2) || onSurface(rPos2) ? phiProp : M_PI - phiProp);
196 
197  //
198  // Solve the obvious two-pi ambiguities: more than one turn!
199  //
200  if (fabs(phiProp - phi0) > 2. * M_PI)
201  phiProp += (phiProp > phi0 ? -2. * M_PI : +2. * M_PI);
202 
203  //
204  // Check that the phi difference has the right sign, according to Q*B
205  // (Another two pi ambuiguity, slightly more subtle)
206  //
207  if ((phiProp - phi0) * radius < 0.0)
208  radius > 0.0 ? phiProp += 2 * M_PI : phiProp -= 2 * M_PI;
209 
210  } else {
211  // Use Taylor
212  phiProp = phi0 + deltaPhi;
213  }
214 
215  //
216  // Check that the particle does not exit from the endcaps first.
217  //
218  zProp = particle_.Z() + (phiProp - phi0) * pZ * radius / pT;
219  if (fabs(zProp) > zCyl) {
220  zProp = zCyl * fabs(pZ) / pZ;
221  phiProp = phi0 + (zProp - particle_.Z()) / radius * pT / pZ;
222  success = 2;
223 
224  } else {
225  //
226  // If the particle does not cross the barrel, either process the
227  // first-half loop only or propagate all the way down to the endcaps
228  //
229  if (rProp < rCyl) {
230  if (firstLoop || fabs(pZ) / pT < 1E-10) {
231  success = 0;
232 
233  } else {
234  zProp = zCyl * fabs(pZ) / pZ;
235  phiProp = phi0 + (zProp - particle_.Z()) / radius * pT / pZ;
236  success = 2;
237  }
238  //
239  // The particle crossed the barrel
240  //
241  } else {
242  success = 1;
243  }
244  }
245  //
246  // Compute the coordinates of the RawParticle after propagation
247  //
248  //
249  // Check the decay time
250  //
251  double delTime = propDir * (phiProp - phi0) * radius * particle_.mass() / pT;
252  double factor = 1.;
253  properTime += delTime;
254  if (properTime > properDecayTime) {
255  factor = 1. - (properTime - properDecayTime) / delTime;
257  decayed = true;
258  }
259 
260  zProp = particle_.Z() + (zProp - particle_.Z()) * factor;
261 
262  phiProp = phi0 + (phiProp - phi0) * factor;
263 
264  double sProp = std::sin(phiProp);
265  double cProp = std::cos(phiProp);
266  double xProp = xC + radius * sProp;
267  double yProp = yC - radius * cProp;
268  double tProp = particle_.T() + (phiProp - phi0) * radius * particle_.E() / pT;
269  double pxProp = pT * cProp;
270  double pyProp = pT * sProp;
271 
272  //
273  // Last check : Although propagated to the endcaps, R could still be
274  // larger than rCyl because the cylinder was too short...
275  //
276  if (success == 2 && xProp * xProp + yProp * yProp > rCyl2) {
277  success = -1;
278  return true;
279  }
280  //
281  // ... and update the particle vertex and momentum
282  //
283  particle_.setVertex(XYZTLorentzVector(xProp, yProp, zProp, tProp));
284  particle_.setMomentum(pxProp, pyProp, pZ, particle_.E());
285  // particle_.SetPx(pxProp);
286  // particle_.SetPy(pyProp);
287  return true;
288  }
289 }
double R2() const
vertex radius**2
Definition: RawParticle.h:291
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
void setMomentum(const XYZTLorentzVector &vtx)
set the momentum
Definition: RawParticle.h:328
double helixCentrePhi() const
The azimuth if the vector joining the cylinder and the helix axes.
double Pz() const
z of the momentum
Definition: RawParticle.h:303
double rCyl
Simulated particle that is to be resp has been propagated.
double helixStartPhi() const
The azimuth of the momentum at the vertex.
bool firstLoop
Do only the first half-loop.
double Z() const
z of vertex
Definition: RawParticle.h:288
int success
0:propagation still be done, 1:reached &#39;barrel&#39;, 2:reached &#39;endcaps&#39;
double properTime
The proper time of the particle.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double helixRadius() const
The helix Radius.
double E() const
energy of the momentum
Definition: RawParticle.h:306
double charge() const
get the MEASURED charge
Definition: RawParticle.h:294
double helixCentreDistToAxis() const
The distance between the cylinder and the helix axes.
double Pt() const
transverse momentum
Definition: RawParticle.h:308
double properDecayTime
The proper decay time of the particle.
T sqrt(T t)
Definition: SSEVec.h:19
double Py() const
y of the momentum
Definition: RawParticle.h:300
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
weight_default_t b2[10]
Definition: b2.h:9
bool onEndcap() const
Is the vertex already on the cylinder endcap ?
double helixCentreX() const
The x coordinate of the helix axis.
int propDir
The propagation direction.
double Perp2() const
perpendicular momentum squared
Definition: RawParticle.h:311
#define M_PI
bool decayed
The particle decayed while propagated !
double T() const
vertex time
Definition: RawParticle.h:289
bool onSurface() const
Is the vertex already on the cylinder surface ?
double mass() const
get the MEASURED mass
Definition: RawParticle.h:295
bool onBarrel() const
Is the vertex already on the cylinder barrel ?
double Y() const
y of vertex
Definition: RawParticle.h:287
double X() const
x of vertex
Definition: RawParticle.h:286
double bField
Magnetic field in the cylinder, oriented along the Z axis.
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
Definition: RawParticle.h:325
bool inside() const
Is the vertex inside the cylinder ? (stricly inside : true)
double helixCentreY() const
The y coordinate of the helix axis.
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25
double Px() const
x of the momentum
Definition: RawParticle.h:297

◆ propagated()

BaseParticlePropagator BaseParticlePropagator::propagated ( ) const

Return a new instance, corresponding to the particle propagated to the surface of the cylinder

Definition at line 304 of file BaseParticlePropagator.cc.

References backPropagate(), firstLoop, propagate(), and success.

Referenced by ParticlePropagator::propagated().

304  {
305  // Copy the input particle in a new instance
306  BaseParticlePropagator myPropPart(*this);
307  // Allow for many loops
308  myPropPart.firstLoop = false;
309  // Propagate the particle !
310  myPropPart.propagate();
311  // Back propagate if the forward propagation was not a success
312  if (myPropPart.success < 0)
313  myPropPart.backPropagate();
314  // Return the new instance
315  return myPropPart;
316 }

◆ propagateToBeamCylinder()

bool BaseParticlePropagator::propagateToBeamCylinder ( const XYZTLorentzVector v,
double  radius = 0. 
)

Definition at line 586 of file BaseParticlePropagator.cc.

References a, b, bField, HltBtagPostValidation_cff::c, c_light(), RawParticle::charge(), funct::cos(), gather_cfg::cout, PVValHelper::dx, PVValHelper::dy, PVValHelper::dz, MillePedeFileConverter_cfg::e, helixCentreDistToAxis(), helixCentreX(), helixCentreY(), mps_fire::i, RawParticle::momentum(), particle_, propagateToClosestApproach(), pv::pT, RawParticle::Pt(), RawParticle::Px(), RawParticle::Py(), alignCSCRings::r, diffTwoXMLs::r2, CosmicsPD_Skims::radius, RawParticle::setCharge(), RawParticle::SetE(), RawParticle::setMomentum(), funct::sin(), mathSSE::sqrt(), findQualityFiles::v, RawParticle::X(), RawParticle::Y(), and RawParticle::Z().

586  {
587  // For neutral (or BField = 0, simply check that the track passes through the cylinder
588  if (particle_.charge() == 0. || bField == 0.)
589  return fabs(particle_.Px() * particle_.Y() - particle_.Py() * particle_.X()) / particle_.Pt() < radius;
590 
591  // Now go to the charged particles
592 
593  // A few needed intermediate variables (to make the code understandable)
594  // The inner cylinder
595  double r = radius;
596  double r2 = r * r;
597  double r4 = r2 * r2;
598  // The two hits
599  double dx = particle_.X() - v.X();
600  double dy = particle_.Y() - v.Y();
601  double dz = particle_.Z() - v.Z();
602  double Sxy = particle_.X() * v.X() + particle_.Y() * v.Y();
603  double Dxy = particle_.Y() * v.X() - particle_.X() * v.Y();
604  double Dxy2 = Dxy * Dxy;
605  double Sxy2 = Sxy * Sxy;
606  double SDxy = dx * dx + dy * dy;
607  double SSDxy = std::sqrt(SDxy);
608  double ATxy = std::atan2(dy, dx);
609  // The coefficients of the second order equation for the trajectory radius
610  double a = r2 - Dxy2 / SDxy;
611  double b = r * (r2 - Sxy);
612  double c = r4 - 2. * Sxy * r2 + Sxy2 + Dxy2;
613 
614  // Here are the four possible solutions for
615  // 1) The trajectory radius
616  std::array<double, 4> helixRs = {{0.}};
617  helixRs[0] = (b - std::sqrt(b * b - a * c)) / (2. * a);
618  helixRs[1] = (b + std::sqrt(b * b - a * c)) / (2. * a);
619  helixRs[2] = -helixRs[0];
620  helixRs[3] = -helixRs[1];
621  // 2) The azimuthal direction at the second point
622  std::array<double, 4> helixPhis = {{0.}};
623  helixPhis[0] = std::asin(SSDxy / (2. * helixRs[0]));
624  helixPhis[1] = std::asin(SSDxy / (2. * helixRs[1]));
625  helixPhis[2] = -helixPhis[0];
626  helixPhis[3] = -helixPhis[1];
627  // Only two solutions are valid though
628  double solution1 = 0.;
629  double solution2 = 0.;
630  double phi1 = 0.;
631  double phi2 = 0.;
632  // Loop over the four possibilities
633  for (unsigned int i = 0; i < 4; ++i) {
634  helixPhis[i] = ATxy + helixPhis[i];
635  // The centre of the helix
636  double xC = helixCentreX(helixRs[i], helixPhis[i]);
637  double yC = helixCentreY(helixRs[i], helixPhis[i]);
638  // The distance of that centre to the origin
639  double dist = helixCentreDistToAxis(xC, yC);
640  /*
641  std::cout << "Solution "<< i << " = " << helixRs[i] << " accepted ? "
642  << fabs(fabs(helixRs[i]) - dist ) << " - " << radius << " = "
643  << fabs(fabs(helixRs[i]) - dist ) - fabs(radius) << " "
644  << ( fabs( fabs(fabs(helixRs[i]) - dist) -fabs(radius) ) < 1e-6) << std::endl;
645  */
646  // Check that the propagation will lead to a trajectroy tangent to the inner cylinder
647  if (fabs(fabs(fabs(helixRs[i]) - dist) - fabs(radius)) < 1e-6) {
648  // Fill first solution
649  if (solution1 == 0.) {
650  solution1 = helixRs[i];
651  phi1 = helixPhis[i];
652  // Fill second solution (ordered)
653  } else if (solution2 == 0.) {
654  if (helixRs[i] < solution1) {
655  solution2 = solution1;
656  solution1 = helixRs[i];
657  phi2 = phi1;
658  phi1 = helixPhis[i];
659  } else {
660  solution2 = helixRs[i];
661  phi2 = helixPhis[i];
662  }
663  // Must not happen!
664  } else {
665  std::cout << "warning !!! More than two solutions for this track !!! " << std::endl;
666  }
667  }
668  }
669 
670  // Find the largest possible pT compatible with coming from within the inne cylinder
671  double pT = 0.;
672  double helixPhi = 0.;
673  double helixR = 0.;
674  // This should not happen
675  if (solution1 * solution2 == 0.) {
676  std::cout << "warning !!! Less than two solution for this track! " << solution1 << " " << solution2 << std::endl;
677  return false;
678  // If the two solutions have opposite sign, it means that pT = infinity is a possibility.
679  } else if (solution1 * solution2 < 0.) {
680  pT = 1000.;
681  double norm = pT / SSDxy;
682  particle_.setCharge(+1.);
683  particle_.setMomentum(dx * norm, dy * norm, dz * norm, 0.);
684  particle_.SetE(std::sqrt(particle_.momentum().Vect().Mag2()));
685  // Otherwise take the solution that gives the largest transverse momentum
686  } else {
687  if (solution1 < 0.) {
688  helixR = solution1;
689  helixPhi = phi1;
690  particle_.setCharge(+1.);
691  } else {
692  helixR = solution2;
693  helixPhi = phi2;
694  particle_.setCharge(-1.);
695  }
696  pT = fabs(helixR) * 1e-5 * c_light() * bField;
697  double norm = pT / SSDxy;
698  particle_.setMomentum(pT * std::cos(helixPhi), pT * std::sin(helixPhi), dz * norm, 0.);
699  particle_.SetE(std::sqrt(particle_.momentum().Vect().Mag2()));
700  }
701 
702  // Propagate to closest approach to get the Z value (a bit of an overkill)
704 }
void setMomentum(const XYZTLorentzVector &vtx)
set the momentum
Definition: RawParticle.h:328
void setCharge(float q)
set the MEASURED charge
Definition: RawParticle.cc:32
double c_light() const
The speed of light in mm/ns (!) without clhep (yeaaahhh!)
double Z() const
z of vertex
Definition: RawParticle.h:288
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void SetE(double)
Definition: RawParticle.h:334
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:321
double charge() const
get the MEASURED charge
Definition: RawParticle.h:294
double helixCentreDistToAxis() const
The distance between the cylinder and the helix axes.
double Pt() const
transverse momentum
Definition: RawParticle.h:308
bool propagateToClosestApproach(double x0=0., double y0=0, bool first=true)
T sqrt(T t)
Definition: SSEVec.h:19
double Py() const
y of the momentum
Definition: RawParticle.h:300
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double helixCentreX() const
The x coordinate of the helix axis.
double b
Definition: hdecay.h:118
double Y() const
y of vertex
Definition: RawParticle.h:287
double a
Definition: hdecay.h:119
double X() const
x of vertex
Definition: RawParticle.h:286
double bField
Magnetic field in the cylinder, oriented along the Z axis.
double helixCentreY() const
The y coordinate of the helix axis.
double Px() const
x of the momentum
Definition: RawParticle.h:297

◆ propagateToClosestApproach()

bool BaseParticlePropagator::propagateToClosestApproach ( double  x0 = 0.,
double  y0 = 0,
bool  first = true 
)

Update the particle after propagation to the closest approach from Z axis, to the preshower layer 1 & 2, to the ECAL entrance, to the HCAL entrance, the HCAL 2nd and 3rd layer (not coded yet), the VFCAL entrance, or any BoundSurface(disk or cylinder)

Definition at line 325 of file BaseParticlePropagator.cc.

References backPropagate(), bField, RawParticle::charge(), fileCollector::done, dqmdumpme::first, helixCentreDistToAxis(), helixCentreX(), helixCentreY(), helixRadius(), helixStartPhi(), RawParticle::momentum(), particle_, pv::pT, RawParticle::Pt(), RawParticle::Px(), RawParticle::Py(), CosmicsPD_Skims::radius, RawParticle::setMomentum(), setPropagationConditions(), RawParticle::setVertex(), mathSSE::sqrt(), RawParticle::vertex(), RawParticle::X(), RawParticle::Y(), and z.

Referenced by propagateToBeamCylinder(), ParticlePropagator::propagateToClosestApproach(), and propagateToNominalVertex().

325  {
326  // Pre-computed quantities
327  double pT = particle_.Pt();
328  double radius = helixRadius(pT);
329  double phi0 = helixStartPhi();
330 
331  // The Helix centre coordinates are computed wrt to the z axis
332  // Actually the z axis might be centred on 0,0: it's the beam spot position (x0,y0)!
333  double xC = helixCentreX(radius, phi0);
334  double yC = helixCentreY(radius, phi0);
335  double distz = helixCentreDistToAxis(xC - x0, yC - y0);
336  double dist0 = helixCentreDistToAxis(xC, yC);
337 
338  //
339  // Propagation to point of clostest approach to z axis
340  //
341  double rz, r0, z;
342  if (particle_.charge() != 0.0 && bField != 0.0) {
343  rz = fabs(fabs(radius) - distz) + std::sqrt(x0 * x0 + y0 * y0) + 0.0000001;
344  r0 = fabs(fabs(radius) - dist0) + 0.0000001;
345  } else {
346  rz = fabs(particle_.Px() * (particle_.Y() - y0) - particle_.Py() * (particle_.X() - x0)) / particle_.Pt();
347  r0 = fabs(particle_.Px() * particle_.Y() - particle_.Py() * particle_.X()) / particle_.Pt();
348  }
349 
350  z = 999999.;
351 
352  // Propagate to the first interesection point
353  // with cylinder of radius sqrt(x0**2+y0**2)
355  bool done = backPropagate();
356 
357  // Unsuccessful propagation - should not happen!
358  if (!done)
359  return done;
360 
361  // The z axis is (0,0) - no need to go further
362  if (fabs(rz - r0) < 1E-10)
363  return done;
364  double dist1 = (particle_.X() - x0) * (particle_.X() - x0) + (particle_.Y() - y0) * (particle_.Y() - y0);
365 
366  // We are already at closest approach - no need to go further
367  if (dist1 < 1E-10)
368  return done;
369 
370  // Keep for later if it happens to be the right solution
371  XYZTLorentzVector vertex1 = particle_.vertex();
372  XYZTLorentzVector momentum1 = particle_.momentum();
373 
374  // Propagate to the distance of closest approach to (0,0)
376  done = backPropagate();
377  if (!done)
378  return done;
379 
380  // Propagate to the first interesection point
381  // with cylinder of radius sqrt(x0**2+y0**2)
383  done = backPropagate();
384  if (!done)
385  return done;
386  double dist2 = (particle_.X() - x0) * (particle_.X() - x0) + (particle_.Y() - y0) * (particle_.Y() - y0);
387 
388  // Keep the good solution.
389  if (dist2 > dist1) {
390  particle_.setVertex(vertex1);
391  particle_.setMomentum(momentum1.X(), momentum1.Y(), momentum1.Z(), momentum1.E());
392  }
393 
394  // Done
395  return done;
396 }
void setMomentum(const XYZTLorentzVector &vtx)
set the momentum
Definition: RawParticle.h:328
double helixStartPhi() const
The azimuth of the momentum at the vertex.
double helixRadius() const
The helix Radius.
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:321
double charge() const
get the MEASURED charge
Definition: RawParticle.h:294
double helixCentreDistToAxis() const
The distance between the cylinder and the helix axes.
double Pt() const
transverse momentum
Definition: RawParticle.h:308
T sqrt(T t)
Definition: SSEVec.h:19
double Py() const
y of the momentum
Definition: RawParticle.h:300
double helixCentreX() const
The x coordinate of the helix axis.
double Y() const
y of vertex
Definition: RawParticle.h:287
double X() const
x of vertex
Definition: RawParticle.h:286
double bField
Magnetic field in the cylinder, oriented along the Z axis.
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
Definition: RawParticle.h:325
double helixCentreY() const
The y coordinate of the helix axis.
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:320
double Px() const
x of the momentum
Definition: RawParticle.h:297

◆ propagateToEcal()

bool BaseParticlePropagator::propagateToEcal ( bool  first = true)

Definition at line 398 of file BaseParticlePropagator.cc.

References dqmdumpme::first, propagate(), and setPropagationConditions().

398  {
399  //
400  // Propagation to Ecal (including preshower) after the
401  // last Tracker layer
402  // TODO: include proper geometry
403  // Geometry taken from CMS ECAL TDR
404  //
405  // First propagate to global barrel / endcap cylinder
406  // setPropagationConditions(1290. , 3045 , first);
407  // New position of the preshower as of CMSSW_1_7_X
408  setPropagationConditions(129.0, 303.353, first);
409  return propagate();
410 }
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)

◆ propagateToEcalEntrance()

bool BaseParticlePropagator::propagateToEcalEntrance ( bool  first = true)

Definition at line 450 of file BaseParticlePropagator.cc.

References RawParticle::cos2ThetaV(), fileCollector::done, dqmdumpme::first, particle_, propagate(), setPropagationConditions(), and success.

Referenced by PFTrackTransformer::addPointsAndBrems(), reco::tau::atECALEntrance(), lowptgsfeleid::features_V0(), FBaseSimEvent::fill(), ConvBremSeedProducer::GoodCluster(), reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::operator()(), GoodSeedProducer::produce(), l1tpf::propagateToCalo(), TrajectoryManager::propagateToCalorimeters(), and LowPtGsfElectronSeedProducer::propagateTrackToCalo().

450  {
451  //
452  // Propagation to ECAL entrance
453  // TODO: include proper geometry
454  // Geometry taken from CMS ECAL TDR
455  //
456  // First propagate to global barrel / endcap cylinder
457  setPropagationConditions(129.0, 320.9, first);
458  bool done = propagate();
459 
460  // Go to endcap cylinder in the "barrel cut corner"
461  // eta = 1.479 -> cos^2(theta) = 0.81230
462  // if ( done && eta > 1.479 && success == 1 ) {
463  if (done && particle_.cos2ThetaV() > 0.81230 && success == 1) {
464  setPropagationConditions(152.6, 320.9, first);
465  done = propagate();
466  }
467 
468  // We are not in the ECAL acceptance
469  // eta = 3.0 -> cos^2(theta) = 0.99013
470  if (particle_.cos2ThetaV() > 0.99014)
471  success = 0;
472 
473  return done;
474 }
int success
0:propagation still be done, 1:reached &#39;barrel&#39;, 2:reached &#39;endcaps&#39;
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)
double cos2ThetaV() const
Definition: RawParticle.h:281

◆ propagateToHcalEntrance()

bool BaseParticlePropagator::propagateToHcalEntrance ( bool  first = true)

Definition at line 476 of file BaseParticlePropagator.cc.

References RawParticle::cos2ThetaV(), fileCollector::done, dqmdumpme::first, particle_, propagate(), propDir, setPropagationConditions(), and success.

Referenced by PFTrackTransformer::addPointsAndBrems(), FBaseSimEvent::fill(), and TrajectoryManager::propagateToCalorimeters().

476  {
477  //
478  // Propagation to HCAL entrance
479  // TODO: include proper geometry
480  // Geometry taken from CMSSW_3_1_X xml (Sunanda)
481  //
482 
483  // First propagate to global barrel / endcap cylinder
484  setPropagationConditions(177.5, 335.0, first);
485  propDir = 0;
486  bool done = propagate();
487  propDir = 1;
488 
489  // If went through the bottom of HB cylinder -> re-propagate to HE surface
490  if (done && success == 2) {
491  setPropagationConditions(300.0, 400.458, first);
492  propDir = 0;
493  done = propagate();
494  propDir = 1;
495  }
496 
497  // out of the HB/HE acceptance
498  // eta = 3.0 -> cos^2(theta) = 0.99014
499  if (done && particle_.cos2ThetaV() > 0.99014)
500  success = 0;
501 
502  return done;
503 }
int success
0:propagation still be done, 1:reached &#39;barrel&#39;, 2:reached &#39;endcaps&#39;
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)
double cos2ThetaV() const
Definition: RawParticle.h:281
int propDir
The propagation direction.

◆ propagateToHcalExit()

bool BaseParticlePropagator::propagateToHcalExit ( bool  first = true)

Definition at line 528 of file BaseParticlePropagator.cc.

References fileCollector::done, dqmdumpme::first, propagate(), propDir, and setPropagationConditions().

Referenced by PFTrackTransformer::addPointsAndBrems(), and FBaseSimEvent::fill().

528  {
529  //
530  // Propagation to HCAL exit
531  // TODO: include proper geometry
532  // Geometry taken from CMSSW_3_1_X xml (Sunanda)
533  //
534 
535  // Approximate it to a single cylinder as it is not that crucial.
536  setPropagationConditions(263.9, 554.1, first);
537  // this->rawPart().particle_.setCharge(0.0); ?? Shower Propagation ??
538  propDir = 0;
539  bool done = propagate();
540  propDir = 1;
541 
542  return done;
543 }
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)
int propDir
The propagation direction.

◆ propagateToHOLayer()

bool BaseParticlePropagator::propagateToHOLayer ( bool  first = true)

Definition at line 545 of file BaseParticlePropagator.cc.

References funct::abs(), fileCollector::done, dqmdumpme::first, particle_, propagate(), propDir, setPropagationConditions(), success, and RawParticle::Z().

Referenced by PFTrackTransformer::addPointsAndBrems(), and FBaseSimEvent::fill().

545  {
546  //
547  // Propagation to Layer0 of HO (Ring-0)
548  // TODO: include proper geometry
549  // Geometry taken from CMSSW_3_1_X xml (Sunanda)
550  //
551 
552  // Approximate it to a single cylinder as it is not that crucial.
553  setPropagationConditions(387.6, 1000.0, first); //Dia is the middle of HO \sqrt{384.8**2+46.2**2}
554  // setPropagationConditions(410.2, 1000.0, first); //sqrt{406.6**2+54.2**2} //for outer layer
555  // this->rawPart().setCharge(0.0); ?? Shower Propagation ??
556  propDir = 0;
557  bool done = propagate();
558  propDir = 1;
559 
560  if (done && std::abs(particle_.Z()) > 700.25)
561  success = 0;
562 
563  return done;
564 }
double Z() const
z of vertex
Definition: RawParticle.h:288
int success
0:propagation still be done, 1:reached &#39;barrel&#39;, 2:reached &#39;endcaps&#39;
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int propDir
The propagation direction.

◆ propagateToNominalVertex()

bool BaseParticlePropagator::propagateToNominalVertex ( const XYZTLorentzVector hit2 = XYZTLorentzVector(0., 0., 0., 0.))

Definition at line 566 of file BaseParticlePropagator.cc.

References bField, RawParticle::charge(), funct::cos(), PVValHelper::dx, PVValHelper::dy, helixRadius(), particle_, phi, propagateToClosestApproach(), pv::pT, RawParticle::pt(), RawParticle::SetPx(), RawParticle::SetPy(), funct::sin(), mathSSE::sqrt(), findQualityFiles::v, RawParticle::X(), and RawParticle::Y().

Referenced by ParticlePropagator::propagateToNominalVertex().

566  {
567  // Not implemented for neutrals (used for electrons only)
568  if (particle_.charge() == 0. || bField == 0.)
569  return false;
570 
571  // Define the proper pT direction to meet the point (vx,vy) in (x,y)
572  double dx = particle_.X() - v.X();
573  double dy = particle_.Y() - v.Y();
574  double phi = std::atan2(dy, dx) + std::asin(std::sqrt(dx * dx + dy * dy) / (2. * helixRadius()));
575 
576  // The absolute pT (kept to the original value)
577  double pT = particle_.pt();
578 
579  // Set the new pT
582 
583  return propagateToClosestApproach(v.X(), v.Y());
584 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double helixRadius() const
The helix Radius.
double charge() const
get the MEASURED charge
Definition: RawParticle.h:294
double pt() const
transverse momentum
Definition: RawParticle.h:309
bool propagateToClosestApproach(double x0=0., double y0=0, bool first=true)
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
void SetPx(double)
Definition: RawParticle.h:331
void SetPy(double)
Definition: RawParticle.h:332
double Y() const
y of vertex
Definition: RawParticle.h:287
double X() const
x of vertex
Definition: RawParticle.h:286
double bField
Magnetic field in the cylinder, oriented along the Z axis.

◆ propagateToPreshowerLayer1()

bool BaseParticlePropagator::propagateToPreshowerLayer1 ( bool  first = true)

Definition at line 412 of file BaseParticlePropagator.cc.

References fileCollector::done, dqmdumpme::first, particle_, propagate(), RawParticle::R2(), setPropagationConditions(), and success.

Referenced by PFTrackTransformer::addPointsAndBrems(), FBaseSimEvent::fill(), and TrajectoryManager::propagateToCalorimeters().

412  {
413  //
414  // Propagation to Preshower Layer 1
415  // TODO: include proper geometry
416  // Geometry taken from CMS ECAL TDR
417  //
418  // First propagate to global barrel / endcap cylinder
419  // setPropagationConditions(1290., 3045 , first);
420  // New position of the preshower as of CMSSW_1_7_X
421  setPropagationConditions(129.0, 303.353, first);
422  bool done = propagate();
423 
424  // Check that were are on the Layer 1
425  if (done && (particle_.R2() > 125.0 * 125.0 || particle_.R2() < 45.0 * 45.0))
426  success = 0;
427 
428  return done;
429 }
double R2() const
vertex radius**2
Definition: RawParticle.h:291
int success
0:propagation still be done, 1:reached &#39;barrel&#39;, 2:reached &#39;endcaps&#39;
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)

◆ propagateToPreshowerLayer2()

bool BaseParticlePropagator::propagateToPreshowerLayer2 ( bool  first = true)

Definition at line 431 of file BaseParticlePropagator.cc.

References fileCollector::done, dqmdumpme::first, particle_, propagate(), RawParticle::R2(), setPropagationConditions(), and success.

Referenced by PFTrackTransformer::addPointsAndBrems(), FBaseSimEvent::fill(), and TrajectoryManager::propagateToCalorimeters().

431  {
432  //
433  // Propagation to Preshower Layer 2
434  // TODO: include proper geometry
435  // Geometry taken from CMS ECAL TDR
436  //
437  // First propagate to global barrel / endcap cylinder
438  // setPropagationConditions(1290. , 3090 , first);
439  // New position of the preshower as of CMSSW_1_7_X
440  setPropagationConditions(129.0, 307.838, first);
441  bool done = propagate();
442 
443  // Check that we are on Layer 2
444  if (done && (particle_.R2() > 125.0 * 125.0 || particle_.R2() < 45.0 * 45.0))
445  success = 0;
446 
447  return done;
448 }
double R2() const
vertex radius**2
Definition: RawParticle.h:291
int success
0:propagation still be done, 1:reached &#39;barrel&#39;, 2:reached &#39;endcaps&#39;
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)

◆ propagateToVFcalEntrance()

bool BaseParticlePropagator::propagateToVFcalEntrance ( bool  first = true)

Definition at line 505 of file BaseParticlePropagator.cc.

References RawParticle::cos2ThetaV(), fileCollector::done, dqmdumpme::first, particle_, propagate(), propDir, setPropagationConditions(), and success.

Referenced by FBaseSimEvent::fill(), and TrajectoryManager::propagateToCalorimeters().

505  {
506  //
507  // Propagation to VFCAL (HF) entrance
508  // Geometry taken from xml: Geometry/CMSCommonData/data/cms.xml
509 
510  setPropagationConditions(400.0, 1114.95, first);
511  propDir = 0;
512  bool done = propagate();
513  propDir = 1;
514 
515  if (!done)
516  success = 0;
517 
518  // We are not in the VFCAL acceptance
519  // eta = 3.0 -> cos^2(theta) = 0.99014
520  // eta = 5.2 -> cos^2(theta) = 0.9998755
521  double c2teta = particle_.cos2ThetaV();
522  if (done && (c2teta < 0.99014 || c2teta > 0.9998755))
523  success = 0;
524 
525  return done;
526 }
int success
0:propagation still be done, 1:reached &#39;barrel&#39;, 2:reached &#39;endcaps&#39;
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)
double cos2ThetaV() const
Definition: RawParticle.h:281
int propDir
The propagation direction.

◆ resetDebug()

void BaseParticlePropagator::resetDebug ( )
inline

Definition at line 306 of file BaseParticlePropagator.h.

References debug.

306 { debug = false; }
bool debug
The debug level.

◆ setDebug()

void BaseParticlePropagator::setDebug ( )
inline

Set the debug leve;.

Definition at line 305 of file BaseParticlePropagator.h.

References debug.

305 { debug = true; }
bool debug
The debug level.

◆ setMagneticField()

void BaseParticlePropagator::setMagneticField ( double  b)
inline

◆ setParticle()

void BaseParticlePropagator::setParticle ( RawParticle const &  iParticle)
inline

Definition at line 166 of file BaseParticlePropagator.h.

References particle_.

166 { particle_ = iParticle; }

◆ setPropagationConditions()

void BaseParticlePropagator::setPropagationConditions ( double  r,
double  z,
bool  firstLoop = true 
)

Set the propagation characteristics (rCyl, zCyl and first loop only)

Definition at line 318 of file BaseParticlePropagator.cc.

References f, firstLoop, dttmaxenums::R, rCyl, rCyl2, BeamSpotPI::Z, and zCyl.

Referenced by PFTrackTransformer::addPoints(), PFTrackTransformer::addPointsAndBrems(), propagateToClosestApproach(), propagateToEcal(), propagateToEcalEntrance(), propagateToHcalEntrance(), propagateToHcalExit(), propagateToHOLayer(), propagateToPreshowerLayer1(), propagateToPreshowerLayer2(), propagateToVFcalEntrance(), and ParticlePropagator::setPropagationConditions().

318  {
319  rCyl = R;
320  rCyl2 = R * R;
321  zCyl = Z;
322  firstLoop = f;
323 }
double zCyl
Half-height of the cylinder (centred at 0,0,0) to which propagation is done.
double rCyl
Simulated particle that is to be resp has been propagated.
bool firstLoop
Do only the first half-loop.
double f[11][100]

◆ setProperDecayTime()

void BaseParticlePropagator::setProperDecayTime ( double  t)
inline

Set the proper decay time.

Definition at line 169 of file BaseParticlePropagator.h.

References properDecayTime, and submitPVValidationJobs::t.

Referenced by ParticlePropagator::initProperDecayTime(), and ParticlePropagator::ParticlePropagator().

169 { properDecayTime = t; }
double properDecayTime
The proper decay time of the particle.

◆ xyImpactParameter()

double BaseParticlePropagator::xyImpactParameter ( double  x0 = 0.,
double  y0 = 0. 
) const

Transverse impact parameter.

Definition at line 706 of file BaseParticlePropagator.cc.

References bField, RawParticle::charge(), helixCentreDistToAxis(), helixCentreX(), helixCentreY(), helixRadius(), helixStartPhi(), particle_, pv::pT, RawParticle::Pt(), RawParticle::Px(), RawParticle::Py(), CosmicsPD_Skims::radius, RawParticle::X(), and RawParticle::Y().

706  {
707  double ip = 0.;
708  double pT = particle_.Pt();
709 
710  if (particle_.charge() != 0.0 && bField != 0.0) {
711  double radius = helixRadius(pT);
712  double phi0 = helixStartPhi();
713 
714  // The Helix centre coordinates are computed wrt to the z axis
715  double xC = helixCentreX(radius, phi0);
716  double yC = helixCentreY(radius, phi0);
717  double distz = helixCentreDistToAxis(xC - x0, yC - y0);
718  ip = distz - fabs(radius);
719  } else {
720  ip = fabs(particle_.Px() * (particle_.Y() - y0) - particle_.Py() * (particle_.X() - x0)) / pT;
721  }
722 
723  return ip;
724 }
double helixStartPhi() const
The azimuth of the momentum at the vertex.
double helixRadius() const
The helix Radius.
double charge() const
get the MEASURED charge
Definition: RawParticle.h:294
double helixCentreDistToAxis() const
The distance between the cylinder and the helix axes.
double Pt() const
transverse momentum
Definition: RawParticle.h:308
double Py() const
y of the momentum
Definition: RawParticle.h:300
double helixCentreX() const
The x coordinate of the helix axis.
double Y() const
y of vertex
Definition: RawParticle.h:287
double X() const
x of vertex
Definition: RawParticle.h:286
double bField
Magnetic field in the cylinder, oriented along the Z axis.
double helixCentreY() const
The y coordinate of the helix axis.
double Px() const
x of the momentum
Definition: RawParticle.h:297

◆ zImpactParameter()

double BaseParticlePropagator::zImpactParameter ( double  x0 = 0,
double  y0 = 0. 
) const
inline

Longitudinal impact parameter.

Definition at line 181 of file BaseParticlePropagator.h.

References particle_, RawParticle::Perp2(), RawParticle::Pz(), mathSSE::sqrt(), RawParticle::X(), RawParticle::Y(), and RawParticle::Z().

181  {
182  // Longitudinal impact parameter
183  return particle_.Z() - particle_.Pz() * std::sqrt(((particle_.X() - x0) * (particle_.X() - x0) +
184  (particle_.Y() - y0) * (particle_.Y() - y0)) /
185  particle_.Perp2());
186  }
double Pz() const
z of the momentum
Definition: RawParticle.h:303
double Z() const
z of vertex
Definition: RawParticle.h:288
T sqrt(T t)
Definition: SSEVec.h:19
double Perp2() const
perpendicular momentum squared
Definition: RawParticle.h:311
double Y() const
y of vertex
Definition: RawParticle.h:287
double X() const
x of vertex
Definition: RawParticle.h:286

Member Data Documentation

◆ bField

double BaseParticlePropagator::bField
private

Magnetic field in the cylinder, oriented along the Z axis.

Definition at line 137 of file BaseParticlePropagator.h.

Referenced by getMagneticField(), helixRadius(), propagate(), propagateToBeamCylinder(), propagateToClosestApproach(), propagateToNominalVertex(), setMagneticField(), and xyImpactParameter().

◆ debug

bool BaseParticlePropagator::debug
private

◆ decayed

bool BaseParticlePropagator::decayed
private

The particle decayed while propagated !

Definition at line 156 of file BaseParticlePropagator.h.

Referenced by hasDecayed(), init(), and propagate().

◆ fiducial

bool BaseParticlePropagator::fiducial
protected

The particle traverses some real material.

Definition at line 147 of file BaseParticlePropagator.h.

Referenced by init(), onFiducial(), and ParticlePropagator::propagateToBoundSurface().

◆ firstLoop

bool BaseParticlePropagator::firstLoop
private

Do only the first half-loop.

Definition at line 154 of file BaseParticlePropagator.h.

Referenced by init(), propagate(), propagated(), ParticlePropagator::setPropagationConditions(), and setPropagationConditions().

◆ particle_

RawParticle BaseParticlePropagator::particle_
private

◆ propDir

int BaseParticlePropagator::propDir
private

◆ properDecayTime

double BaseParticlePropagator::properDecayTime
private

The proper decay time of the particle.

Definition at line 139 of file BaseParticlePropagator.h.

Referenced by ParticlePropagator::initProperDecayTime(), propagate(), and setProperDecayTime().

◆ properTime

double BaseParticlePropagator::properTime
private

The proper time of the particle.

Definition at line 158 of file BaseParticlePropagator.h.

Referenced by init(), and propagate().

◆ rCyl

double BaseParticlePropagator::rCyl
private

Simulated particle that is to be resp has been propagated.

Radius of the cylinder (centred at 0,0,0) to which propagation is done

Definition at line 132 of file BaseParticlePropagator.h.

Referenced by increaseRCyl(), inside(), onBarrel(), propagate(), and setPropagationConditions().

◆ rCyl2

double BaseParticlePropagator::rCyl2
private

◆ success

int BaseParticlePropagator::success
protected

◆ zCyl

double BaseParticlePropagator::zCyl
private

Half-height of the cylinder (centred at 0,0,0) to which propagation is done.

Definition at line 135 of file BaseParticlePropagator.h.

Referenced by inside(), onBarrel(), onEndcap(), propagate(), and setPropagationConditions().