CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
AnalyticalPropagator Class Reference

#include <AnalyticalPropagator.h>

Inheritance diagram for AnalyticalPropagator:
Propagator

Public Member Functions

 AnalyticalPropagator (const MagneticField *field, PropagationDirection dir=alongMomentum, float maxDPhi=1.6, bool isOld=true)
 
virtual AnalyticalPropagatorclone () const
 
virtual bool setMaxDirectionChange (float phiMax)
 
void setMaxRelativeChangeInBz (const float maxDBz)
 
 ~AnalyticalPropagator ()
 
- Public Member Functions inherited from Propagator
virtual FreeTrajectoryState propagate (const FreeTrajectoryState &ftsStart, const GlobalPoint &pDest) const final
 
virtual FreeTrajectoryState propagate (const FreeTrajectoryState &ftsStart, const GlobalPoint &pDest1, const GlobalPoint &pDest2) const final
 
virtual FreeTrajectoryState propagate (const FreeTrajectoryState &ftsStart, const reco::BeamSpot &beamSpot) const final
 
virtual std::pair
< TrajectoryStateOnSurface,
double > 
propagateWithPath (const FreeTrajectoryState &, const Surface &) const final
 
virtual std::pair
< TrajectoryStateOnSurface,
double > 
propagateWithPath (const TrajectoryStateOnSurface &tsos, const Surface &sur) const final
 
virtual std::pair
< TrajectoryStateOnSurface,
double > 
propagateWithPath (const TrajectoryStateOnSurface &tsos, const Plane &sur) const
 
virtual std::pair
< TrajectoryStateOnSurface,
double > 
propagateWithPath (const TrajectoryStateOnSurface &tsos, const Cylinder &sur) const
 
virtual std::pair
< FreeTrajectoryState, double > 
propagateWithPath (const FreeTrajectoryState &ftsStart, const GlobalPoint &pDest) const
 
virtual std::pair
< FreeTrajectoryState, double > 
propagateWithPath (const FreeTrajectoryState &ftsStart, const GlobalPoint &pDest1, const GlobalPoint &pDest2) const
 Propagate to PCA to a line (given by 2 points) given a starting point. More...
 
virtual std::pair
< FreeTrajectoryState, double > 
propagateWithPath (const FreeTrajectoryState &ftsStart, const reco::BeamSpot &beamSpot) const
 Propagate to PCA to a line (given by beamSpot position and slope) given a starting point. More...
 
virtual PropagationDirection propagationDirection () const
 
 Propagator (PropagationDirection dir=alongMomentum)
 
virtual void setPropagationDirection (PropagationDirection dir)
 
virtual ~Propagator ()
 

Private Types

typedef std::pair
< TrajectoryStateOnSurface,
double > 
TsosWP
 

Private Member Functions

virtual const MagneticFieldmagneticField () const
 
std::pair
< TrajectoryStateOnSurface,
double > 
propagatedStateWithPath (const FreeTrajectoryState &fts, const Surface &surface, const GlobalTrajectoryParameters &gtp, const double &s) const
 propagation of errors (if needed) and generation of a new TSOS More...
 
bool propagateParametersOnCylinder (const FreeTrajectoryState &fts, const Cylinder &cylinder, GlobalPoint &x, GlobalVector &p, double &s) const
 parameter propagation to cylinder (returns position, momentum and path length) More...
 
bool propagateParametersOnPlane (const FreeTrajectoryState &fts, const Plane &plane, GlobalPoint &x, GlobalVector &p, double &s) const
 parameter propagation to plane (returns position, momentum and path length) More...
 
bool propagateWithHelixCrossing (HelixPlaneCrossing &, const Plane &, const float, GlobalPoint &, GlobalVector &, double &s) const
 helix parameter propagation to a plane using HelixPlaneCrossing More...
 
bool propagateWithLineCrossing (const GlobalPoint &, const GlobalVector &, const Plane &, GlobalPoint &, double &) const
 straight line parameter propagation to a plane More...
 
bool propagateWithLineCrossing (const GlobalPoint &, const GlobalVector &, const Cylinder &, GlobalPoint &, double &) const
 straight line parameter propagation to a cylinder More...
 
std::pair
< TrajectoryStateOnSurface,
double > 
propagateWithPath (const FreeTrajectoryState &fts, const Plane &plane) const override
 propagation to plane with path length More...
 
std::pair
< TrajectoryStateOnSurface,
double > 
propagateWithPath (const FreeTrajectoryState &fts, const Cylinder &cylinder) const override
 propagation to cylinder with path length More...
 

Private Attributes

bool isOldPropagationType
 
const MagneticFieldtheField
 
float theMaxDBzRatio
 
float theMaxDPhi2
 

Detailed Description

(Mostly) analytical helix propagation to cylindrical or planar surfaces. Based on GtfGeometricalPropagator with successive replacement of components (currently: propagation to arbitrary plane).

Definition at line 23 of file AnalyticalPropagator.h.

Member Typedef Documentation

typedef std::pair<TrajectoryStateOnSurface,double> AnalyticalPropagator::TsosWP
private

Definition at line 119 of file AnalyticalPropagator.h.

Constructor & Destructor Documentation

AnalyticalPropagator::AnalyticalPropagator ( const MagneticField field,
PropagationDirection  dir = alongMomentum,
float  maxDPhi = 1.6,
bool  isOld = true 
)
inline

Definition at line 27 of file AnalyticalPropagator.h.

Referenced by clone().

29  :
30  Propagator(dir),
31  theMaxDPhi2(maxDPhi*maxDPhi),
32  theMaxDBzRatio(0.5),
33  theField(field),
34  isOldPropagationType(isOld) {}
Propagator(PropagationDirection dir=alongMomentum)
Definition: Propagator.h:48
const MagneticField * theField
dbl *** dir
Definition: mlp_gen.cc:35
AnalyticalPropagator::~AnalyticalPropagator ( )
inline

Definition at line 36 of file AnalyticalPropagator.h.

36 {}

Member Function Documentation

virtual AnalyticalPropagator* AnalyticalPropagator::clone ( void  ) const
inlinevirtual

Implements Propagator.

Definition at line 70 of file AnalyticalPropagator.h.

References AnalyticalPropagator().

71  {
72  return new AnalyticalPropagator(*this);
73  }
AnalyticalPropagator(const MagneticField *field, PropagationDirection dir=alongMomentum, float maxDPhi=1.6, bool isOld=true)
virtual const MagneticField* AnalyticalPropagator::magneticField ( ) const
inlineprivatevirtual

Implements Propagator.

Definition at line 116 of file AnalyticalPropagator.h.

References theField.

116 {return theField;}
const MagneticField * theField
std::pair< TrajectoryStateOnSurface, double > AnalyticalPropagator::propagatedStateWithPath ( const FreeTrajectoryState fts,
const Surface surface,
const GlobalTrajectoryParameters gtp,
const double &  s 
) const
private

propagation of errors (if needed) and generation of a new TSOS

Definition at line 111 of file AnalyticalPropagator.cc.

References SurfaceSideDefinition::afterSurface, alongMomentum, SurfaceSideDefinition::beforeSurface, FreeTrajectoryState::curvilinearError(), FreeTrajectoryState::hasError(), AnalyticalCurvilinearJacobian::jacobian(), CurvilinearTrajectoryError::matrix(), GlobalTrajectoryParameters::momentum(), FreeTrajectoryState::parameters(), GlobalTrajectoryParameters::position(), and alignCSCRings::s.

115 {
116  //
117  // for forward propagation: state is before surface,
118  // for backward propagation: state is after surface
119  //
122  //
123  //
124  // error propagation (if needed) and conversion to a TrajectoryStateOnSurface
125  //
126  if (fts.hasError()) {
127  //
128  // compute jacobian
129  //
130  AnalyticalCurvilinearJacobian analyticalJacobian(fts.parameters(), gtp.position(), gtp.momentum(), s);
131  const AlgebraicMatrix55 &jacobian = analyticalJacobian.jacobian();
132  // CurvilinearTrajectoryError cte(ROOT::Math::Similarity(jacobian, fts.curvilinearError().matrix()));
133  return TsosWP(TrajectoryStateOnSurface(gtp,
134  ROOT::Math::Similarity(jacobian, fts.curvilinearError().matrix()),
135  surface,side),s);
136  }
137  else {
138  //
139  // return state without errors
140  //
141  return TsosWP(TrajectoryStateOnSurface(gtp,surface,side),s);
142  }
143 }
virtual PropagationDirection propagationDirection() const
Definition: Propagator.h:155
const GlobalTrajectoryParameters & parameters() const
const AlgebraicMatrix55 & jacobian() const
const CurvilinearTrajectoryError & curvilinearError() const
const AlgebraicSymMatrix55 & matrix() const
std::pair< TrajectoryStateOnSurface, double > TsosWP
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
bool AnalyticalPropagator::propagateParametersOnCylinder ( const FreeTrajectoryState fts,
const Cylinder cylinder,
GlobalPoint x,
GlobalVector p,
double &  s 
) const
private

parameter propagation to cylinder (returns position, momentum and path length)

Definition at line 145 of file AnalyticalPropagator.cc.

References funct::abs(), constexpr, alignCSCRings::e, f, PV3DBase< T, PVType, FrameType >::mag(), FreeTrajectoryState::momentum(), AlCaHLTBitMon_ParallelJobs::p, GloballyPositioned< T >::position(), FreeTrajectoryState::position(), rho, alignCSCRings::s, FreeTrajectoryState::transverseCurvature(), Vector3DBase< T, FrameTag >::unit(), unlikely, x, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

148 {
149 
150  GlobalPoint const & sp = cylinder.position();
151  if unlikely(sp.x()!=0. || sp.y()!=0.) {
152  throw PropagationException("Cannot propagate to an arbitrary cylinder");
153  }
154  // preset output
155  x = fts.position();
156  p = fts.momentum();
157  s = 0;
158  // (transverse) curvature
159  auto rho = fts.transverseCurvature();
160  //
161  // Straight line approximation? |rho|<1.e-10 equivalent to ~ 1um
162  // difference in transversal position at 10m.
163  //
164  if unlikely( std::abs(rho)<1.e-10f )
165  return propagateWithLineCrossing(fts.position(),p,cylinder,x,s);
166  //
167  // Helix case
168  //
169  // check for possible intersection
170  constexpr float tolerance = 1.e-4; // 1 micron distance
171  auto rdiff = x.perp() - cylinder.radius();
172  if ( std::abs(rdiff) < tolerance ) return true;
173  //
174  // Instantiate HelixBarrelCylinderCrossing and get solutions
175  //
176  HelixBarrelCylinderCrossing cylinderCrossing(fts.position(),fts.momentum(),rho,
177  propagationDirection(),cylinder);
178  if unlikely( !cylinderCrossing.hasSolution() ) return false;
179  // path length
180  s = cylinderCrossing.pathLength();
181  // point
182  x = cylinderCrossing.position();
183  // direction (renormalised)
184  p = cylinderCrossing.direction().unit()*fts.momentum().mag();
185  return true;
186 }
virtual PropagationDirection propagationDirection() const
Definition: Propagator.h:155
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Definition: DDAxes.h:10
T y() const
Definition: PV3DBase.h:63
bool propagateWithLineCrossing(const GlobalPoint &, const GlobalVector &, const Plane &, GlobalPoint &, double &) const
straight line parameter propagation to a plane
#define constexpr
return((rh^lh)&mask)
#define unlikely(x)
string unit
Definition: csvLumiCalc.py:46
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
Common base class.
GlobalVector momentum() const
GlobalPoint position() const
double transverseCurvature() const
T perp() const
Magnitude of transverse component.
static int position[264][3]
Definition: ReadPGInfo.cc:509
volatile std::atomic< bool > shutdown_flag false
if(conf.exists("allCellsPositionCalc"))
T x() const
Definition: PV3DBase.h:62
const PositionType & position() const
bool AnalyticalPropagator::propagateParametersOnPlane ( const FreeTrajectoryState fts,
const Plane plane,
GlobalPoint x,
GlobalVector p,
double &  s 
) const
private

parameter propagation to plane (returns position, momentum and path length)

Definition at line 189 of file AnalyticalPropagator.cc.

References funct::abs(), constexpr, HelixArbitraryPlaneCrossing::direction(), alignCSCRings::e, f, likely, LogDebug, PV3DBase< T, PVType, FrameType >::mag(), Basic3DVector< T >::mag(), FreeTrajectoryState::momentum(), oppositeToMomentum, AlCaHLTBitMon_ParallelJobs::p, perp(), Geom::Phi< T >::phi(), PV3DBase< T, PVType, FrameType >::phi(), HelixArbitraryPlaneCrossing::position(), FreeTrajectoryState::position(), alignCSCRings::s, indexGen::s2, FreeTrajectoryState::transverseCurvature(), unlikely, and x.

194 {
195  // initialisation of position, momentum and path length
196  x = fts.position();
197  p = fts.momentum();
198  s = 0;
199  // (transverse) curvature
200  auto rho = fts.transverseCurvature();
201  //
202  // Straight line approximation? |rho|<1.e-10 equivalent to ~ 1um
203  // difference in transversal position at 10m.
204  //
205  if unlikely( std::abs(rho)<1.e-10f )
206  return propagateWithLineCrossing(fts.position(),p,plane,x,s);
207  //
208  // Helix case
209  //
210  GlobalVector u = plane.normalVector();
211  constexpr float small = 1.e-6; // for orientation of planes
212  //
213  // Frame-independant point and vector are created explicitely to
214  // avoid confusing gcc (refuses to compile with temporary objects
215  // in the constructor).
216  //
217  HelixPlaneCrossing::PositionType helixPos(x);
218  HelixPlaneCrossing::DirectionType helixDir(p);
220  if (std::abs(u.z()) < small) {
221  // barrel plane:
222  // instantiate HelixBarrelPlaneCrossing, get vector of solutions and check for existance
223  HelixBarrelPlaneCrossingByCircle planeCrossing(helixPos,helixDir,rho,propagationDirection());
224  return propagateWithHelixCrossing(planeCrossing,plane,fts.momentum().mag(),x,p,s);
225  }
226  if ( (std::abs(u.x()) < small) & (std::abs(u.y()) < small) ) {
227  // forward plane:
228  // instantiate HelixForwardPlaneCrossing, get vector of solutions and check for existance
229  HelixForwardPlaneCrossing planeCrossing(helixPos,helixDir,rho,propagationDirection());
230  return propagateWithHelixCrossing(planeCrossing,plane,fts.momentum().mag(),x,p,s);
231  }
232  // arbitrary plane:
233  // instantiate HelixArbitraryPlaneCrossing, get vector of solutions and check for existance
234  HelixArbitraryPlaneCrossing planeCrossing(helixPos,helixDir,rho,propagationDirection());
235  return propagateWithHelixCrossing(planeCrossing,plane,fts.momentum().mag(),x,p,s);
236  }
237 
238 
239 
240 
241  //--- Alternative implementation to be used for the propagation of the parameters of looping
242  // particles that cross twice the (infinite) surface of the plane. It is not trivial to determine
243  // which of the two intersections has to be returned.
244 
245  //---- FIXME: WHAT FOLLOWS HAS TO BE REWRITTEN IN A CLEANER (AND CPU-OPTIMIZED) WAY ---------
246  LogDebug("AnalyticalPropagator") << "In AnaliticalProp, calling HAPC " << "\n"
247  << "plane is centered in xyz: "
248  << plane.position().x() << " , "
249  << plane.position().y() << " , "
250  << plane.position().z() << "\n";
251 
252 
253  GlobalPoint gp1 = fts.position();
254  GlobalVector gm1 = fts.momentum();
255  double s1 = 0;
256  double rho1 = fts.transverseCurvature();
257  HelixPlaneCrossing::PositionType helixPos1(gp1);
258  HelixPlaneCrossing::DirectionType helixDir1(gm1);
259  LogDebug("AnalyticalPropagator") << "gp1 before calling planeCrossing1: " << gp1 << "\n";
260  HelixArbitraryPlaneCrossing planeCrossing1(helixPos1,helixDir1,rho1,propagationDirection());
261 
264 
265  double tolerance(0.0050);
267  tolerance *=-1;
268 
269  bool check1 = propagateWithHelixCrossing(planeCrossing1,plane,fts.momentum().mag(),gp1,gm1,s1);
270  double dphi1 = fabs(fts.momentum().phi()-gm1.phi());
271  LogDebug("AnalyticalPropagator") << "check1, s1, dphi, gp1: "
272  << check1 << " , "
273  << s1 << " , "
274  << dphi1 << " , "
275  << gp1 << "\n";
276 
277  //move forward a bit to avoid that the propagator doesn't propagate because the state is already on surface.
278  //we want to go to the other point of intersection between the helix and the plane
279  xGen = planeCrossing1.position(s1+tolerance);
280  pGen = planeCrossing1.direction(s1+tolerance);
281 
282  /*
283  if(!check1 || s1>170 ){
284  //PropagationDirection newDir = (propagationDirection() == alongMomentum) ? oppositeToMomentum : alongMomentum;
285  PropagationDirection newDir = anyDirection;
286  HelixArbitraryPlaneCrossing planeCrossing1B(helixPos1,helixDir1,rho1,newDir);
287  check1 = propagateWithHelixCrossing(planeCrossing1B,plane,fts.momentum().mag(),gp1,gm1,s1);
288  LogDebug("AnalyticalPropagator") << "after second attempt, check1, s1,gp1: "
289  << check1 << " , "
290  << s1 << " , " << gp1 << "\n";
291 
292  xGen = planeCrossing1B.position(s1+tolerance);
293  pGen = planeCrossing1B.direction(s1+tolerance);
294  }
295  */
296 
297  if(!check1){
298  LogDebug("AnalyticalPropagator") << "failed also second attempt. No idea what to do, then bailout" << "\n";
299  }
300 
301 
302  pGen *= gm1.mag()/pGen.mag();
303  GlobalPoint gp2(xGen);
304  GlobalVector gm2(pGen);
305  double s2 = 0;
306  double rho2 = rho1;
307  HelixPlaneCrossing::PositionType helixPos2(gp2);
308  HelixPlaneCrossing::DirectionType helixDir2(gm2);
309  HelixArbitraryPlaneCrossing planeCrossing2(helixPos2,helixDir2,rho2,propagationDirection());
310 
311  bool check2 = propagateWithHelixCrossing(planeCrossing2,plane,gm2.mag(),gp2,gm2,s2);
312 
313  if(!check2){
314  x = gp1;
315  p = gm1;
316  s = s1;
317  return check1;
318  }
319 
320  if(!check1){
321  edm::LogError("AnalyticalPropagator") << "LOGIC ERROR: I should not have entered here!" << "\n";
322  return false;
323  }
324 
325 
326  LogDebug("AnalyticalPropagator") << "check2, s2, gp2: "
327  << check2 << " , "
328  << s2 << " , " << gp2 << "\n";
329 
330 
331  double dist1 = (plane.position()-gp1).perp();
332  double dist2 = (plane.position()-gp2).perp();
333 
334 
335  LogDebug("AnalyticalPropagator") << "propDir, dist1, dist2: "
336  << propagationDirection() << " , "
337  << dist1 << " , "
338  << dist2 << "\n";
339 
340  //If there are two solutions, the one which is the closest to the module's center is chosen
341  if(dist1<2*dist2){
342  x = gp1;
343  p = gm1;
344  s = s1;
345  return check1;
346  }else if(dist2<2*dist1){
347  x = gp2;
348  p = gm2;
349  s = s1+s2+tolerance;
350  return check2;
351  }else{
352  if(fabs(s1)<fabs(s2)){
353  x = gp1;
354  p = gm1;
355  s = s1;
356  return check1;
357  }else{
358  x = gp2;
359  p = gm2;
360  s = s1+s2+tolerance;
361  return check2;
362  }
363  }
364 
365  //-------- END of ugly piece of code ---------------
366 
367 }
#define LogDebug(id)
virtual PropagationDirection propagationDirection() const
Definition: Propagator.h:155
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Definition: DDAxes.h:10
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
bool propagateWithLineCrossing(const GlobalPoint &, const GlobalVector &, const Plane &, GlobalPoint &, double &) const
straight line parameter propagation to a plane
#define constexpr
tuple s2
Definition: indexGen.py:106
return((rh^lh)&mask)
Point3DBase< Scalar, GlobalTag > PositionType
Definition: Definitions.h:30
#define unlikely(x)
#define likely(x)
T mag() const
Definition: PV3DBase.h:67
T phi() const
Definition: Phi.h:41
bool propagateWithHelixCrossing(HelixPlaneCrossing &, const Plane &, const float, GlobalPoint &, GlobalVector &, double &s) const
helix parameter propagation to a plane using HelixPlaneCrossing
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
GlobalVector momentum() const
GlobalPoint position() const
double transverseCurvature() const
T perp() const
Magnitude of transverse component.
static int position[264][3]
Definition: ReadPGInfo.cc:509
Definition: DDAxes.h:10
if(conf.exists("allCellsPositionCalc"))
bool AnalyticalPropagator::propagateWithHelixCrossing ( HelixPlaneCrossing planeCrossing,
const Plane plane,
const float  pmag,
GlobalPoint x,
GlobalVector p,
double &  s 
) const
private

helix parameter propagation to a plane using HelixPlaneCrossing

Definition at line 420 of file AnalyticalPropagator.cc.

References HelixPlaneCrossing::direction(), PV3DBase< T, PVType, FrameType >::mag(), HelixPlaneCrossing::pathLength(), HelixPlaneCrossing::position(), and unlikely.

425  {
426  // get solution
427  std::pair<bool,double> propResult = planeCrossing.pathLength(plane);
428  if unlikely( !propResult.first ) return false;
429 
430  s = propResult.second;
431  x = GlobalPoint(planeCrossing.position(s));
432  // direction (reconverted to GlobalVector, renormalised)
433  GlobalVector pGen = GlobalVector(planeCrossing.direction(s));
434  pGen *= pmag/pGen.mag();
435  p = pGen;
436  //
437  return true;
438 }
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
return((rh^lh)&mask)
#define unlikely(x)
U second(std::pair< T, U > const &p)
virtual std::pair< bool, double > pathLength(const Plane &)=0
static int position[264][3]
Definition: ReadPGInfo.cc:509
volatile std::atomic< bool > shutdown_flag false
bool AnalyticalPropagator::propagateWithLineCrossing ( const GlobalPoint ,
const GlobalVector ,
const Plane ,
GlobalPoint ,
double &   
) const
private

straight line parameter propagation to a plane

bool AnalyticalPropagator::propagateWithLineCrossing ( const GlobalPoint ,
const GlobalVector ,
const Cylinder ,
GlobalPoint ,
double &   
) const
private

straight line parameter propagation to a cylinder

std::pair<TrajectoryStateOnSurface,double> AnalyticalPropagator::propagateWithPath ( const FreeTrajectoryState fts,
const Plane plane 
) const
overrideprivatevirtual

propagation to plane with path length

Implements Propagator.

std::pair<TrajectoryStateOnSurface,double> AnalyticalPropagator::propagateWithPath ( const FreeTrajectoryState fts,
const Cylinder cylinder 
) const
overrideprivatevirtual

propagation to cylinder with path length

Implements Propagator.

virtual bool AnalyticalPropagator::setMaxDirectionChange ( float  phiMax)
inlinevirtual

limitation of change in transverse direction (to avoid loops).

Reimplemented from Propagator.

Definition at line 65 of file AnalyticalPropagator.h.

References theMaxDPhi2.

65  {
66  theMaxDPhi2 = phiMax*phiMax;
67  return true;
68  }
void AnalyticalPropagator::setMaxRelativeChangeInBz ( const float  maxDBz)
inline

Set the maximum relative change in Bz (Bz_at_end-Bz_at_start)/Bz_at_start for a single propagation. The default is no limit. NB: this propagator assumes constant, non-zero magnetic field parallel to the z-axis!

Definition at line 79 of file AnalyticalPropagator.h.

References theMaxDBzRatio.

79  {
80  theMaxDBzRatio = maxDBz;
81  }

Member Data Documentation

bool AnalyticalPropagator::isOldPropagationType
private

Definition at line 123 of file AnalyticalPropagator.h.

const MagneticField* AnalyticalPropagator::theField
private

Definition at line 122 of file AnalyticalPropagator.h.

Referenced by magneticField().

float AnalyticalPropagator::theMaxDBzRatio
private

Definition at line 121 of file AnalyticalPropagator.h.

Referenced by setMaxRelativeChangeInBz().

float AnalyticalPropagator::theMaxDPhi2
private

Definition at line 120 of file AnalyticalPropagator.h.

Referenced by setMaxDirectionChange().