Loading [MathJax]/extensions/tex2jax.js
ATLAS Offline Software
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
RungeKuttaIntersector.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // Runge-Kutta method for tracking a particle through a magnetic field
7 // (c) ATLAS Detector software
9 
10 #ifndef TRKEXRUNGEKUTTAINTERSECTOR_RUNGEKUTTAINTERSECTOR_H
11 #define TRKEXRUNGEKUTTAINTERSECTOR_RUNGEKUTTAINTERSECTOR_H
12 
13 #include <atomic>
14 
17 #include "GaudiKernel/PhysicalConstants.h"
18 #include "GaudiKernel/ToolHandle.h"
25 
26 namespace Trk {
27 
28 class RungeKuttaIntersector final : public extends<AthAlgTool, IIntersector> {
29 
30  public:
31  RungeKuttaIntersector(const std::string& type, const std::string& name,
32  const IInterface* parent);
33  virtual ~RungeKuttaIntersector() = default;
34 
35  virtual StatusCode initialize() override;
36  virtual StatusCode finalize() override;
37 
39  virtual std::optional<TrackSurfaceIntersection> intersectSurface(
40  const Surface& surface, const TrackSurfaceIntersection& trackIntersection,
41  const double qOverP) const override;
42 
44  virtual std::optional<TrackSurfaceIntersection> approachPerigeeSurface(
45  const PerigeeSurface& surface,
46  const TrackSurfaceIntersection& trackIntersection,
47  const double qOverP) const override;
48 
51  virtual std::optional<TrackSurfaceIntersection> approachStraightLineSurface(
52  const StraightLineSurface& surface,
53  const TrackSurfaceIntersection& trackIntersection,
54  const double qOverP) const override;
55 
58  virtual std::optional<TrackSurfaceIntersection> intersectCylinderSurface(
59  const CylinderSurface& surface,
60  const TrackSurfaceIntersection& trackIntersection,
61  const double qOverP) const override;
62 
64  virtual std::optional<TrackSurfaceIntersection> intersectDiscSurface(
65  const DiscSurface& surface,
66  const TrackSurfaceIntersection& trackIntersection,
67  const double qOverP) const override;
68 
70  virtual std::optional<TrackSurfaceIntersection> intersectPlaneSurface(
71  const PlaneSurface& surface,
72  const TrackSurfaceIntersection& trackIntersection,
73  const double qOverP) const override;
74 
77  virtual bool isValid(Amg::Vector3D /*startPosition*/,
78  Amg::Vector3D /*endPosition*/) const override {
79  return true;
80  }
81 
82  private:
83  // private methods
85  double& stepLength) const;
87  const Surface& surface, const double qOverP,
88  const double rStart, const double zStart,
89  const bool trapped) const;
90  double distanceToCylinder(const TrackSurfaceIntersection& isect,
91  const double cylinderRadius,
92  const double offsetRadius,
93  const Amg::Vector3D& offset,
94  double& stepLength) const;
95  double distanceToDisc(const TrackSurfaceIntersection& isect,
96  const double discZ, double& stepLength) const;
97 
98  double distanceToLine(const TrackSurfaceIntersection& isect,
99  const Amg::Vector3D& linePosition,
100  const Amg::Vector3D& lineDirection,
101  double& stepLength) const;
102 
103  double distanceToPlane(const TrackSurfaceIntersection& isect,
104  const Amg::Vector3D& planePosition,
105  const Amg::Vector3D& planeNormal,
106  double& stepLength) const;
107 
108  Amg::Vector3D field(const Amg::Vector3D& point,
109  MagField::AtlasFieldCache& fieldCache) const;
110 
111  void initializeFieldCache(MagField::AtlasFieldCache& fieldCache) const;
112 
113  std::optional<TrackSurfaceIntersection> newIntersection(
114  TrackSurfaceIntersection&& isect, const Surface& surface,
115  const double qOverP, const double rStart, const double zStart) const;
116 
118  const Amg::Vector3D& fieldValue, const double stepLength,
119  const double qOverP) const;
120 
122  double& stepLength, const double qOverP,
123  MagField::AtlasFieldCache& fieldCache) const;
124 
126  this, "AtlasFieldCacheCondObj", "fieldCondObj",
127  "Name of the Magnetic Field conditions object key"};
128 
129  // additional configuration
130  BooleanProperty m_productionMode{this, "ProductionMode", true};
131 
132  // some precalculated constants:
133  // r min for calo high field gradient region
134  const double m_caloR0 = 1900.*Gaudi::Units::mm;
135  // r max for calo high field gradient region
136  const double m_caloR1 = 2500.*Gaudi::Units::mm;
137  // r min for calo medium field gradient region
138  const double m_caloR2 = 3500.*Gaudi::Units::mm;
139  // r min for calo outer flux return region
140  const double m_caloR3 = 3700.*Gaudi::Units::mm;
141  // r max for calo medium field gradient region
142  const double m_caloR4 = 3800.*Gaudi::Units::mm;
143  // z min for calo medium field gradient region
144  const double m_caloZ0 = 2350.*Gaudi::Units::mm;
145  // z min for calo high field gradient region
146  const double m_caloZ1 = 2600.*Gaudi::Units::mm;
147  // z max for calo high field gradient region
148  const double m_caloZ2 = 3600.*Gaudi::Units::mm;
149  // z min for calo outer flux return region
150  const double m_caloZ3 = 6000.*Gaudi::Units::mm;
151  // end of central barrel near constant field region
152  const double m_inDetR0 = 400.*Gaudi::Units::mm;
153  // inner radius of middle/outer transition region
154  const double m_inDetR1 = 350.*Gaudi::Units::mm;
155  // outer radius of low field gradient field region
156  const double m_inDetR2 = 800.*Gaudi::Units::mm;
157  // end of central barrel near constant field region
158  const double m_inDetZ0 = 350.*Gaudi::Units::mm;
159  // start of well behaved transition region
160  const double m_inDetZ1 = 420.*Gaudi::Units::mm;
161  // start of endcap region
162  const double m_inDetZ2 = 700.*Gaudi::Units::mm;
163  // protection against loopers
164  const double m_momentumThreshold = 1./20.*Gaudi::Units::MeV;
165  // warning threshold for intersection failure
167  // inner radius of barrel toroid region
168  const double m_muonR0 = 4300.*Gaudi::Units::mm;
169  // start of endcap toroid region
170  const double m_muonZ0 = 6600.*Gaudi::Units::mm;
173  // r max after coil (will take small steps near coil)
174  const double m_solenoidR = 1300.*Gaudi::Units::mm;
175  // z end of InDet region
176  const double m_solenoidZ = 3500.*Gaudi::Units::mm;
180  double m_stepMax3 = 160.0*Gaudi::Units::mm;
181  double m_stepMax4 = 320.0*Gaudi::Units::mm;
183  const double m_third = 1./3.;
184  // endcap toroid central field - inner radius
185  const double m_toroidR0 = 1850.0*Gaudi::Units::mm;
186  // endcap toroid central field - outer radius
187  const double m_toroidR1 = 3500.0*Gaudi::Units::mm;
188  // after inner barrel or outer endcap coil
189  const double m_toroidR2 = 6000.0*Gaudi::Units::mm;
190  // toroid region - radius above which long steps OK
191  const double m_toroidR3 = 6500.0*Gaudi::Units::mm;
192  // endcap - near iron structure
193  const double m_toroidZ0 = 7000.0*Gaudi::Units::mm;
194  // endcap - high field gradient
195  const double m_toroidZ1 = 8000.0*Gaudi::Units::mm;
196  // barrel - before coil end loop
197  const double m_toroidZ2 = 8700.0*Gaudi::Units::mm;
198  // barrel - nearing coil end loop
199  const double m_toroidZ3 = 9100.0*Gaudi::Units::mm;
200  // endcap toroid full field and barrel coil regions
201  const double m_toroidZ4 = 9500.0*Gaudi::Units::mm;
202  // endcap toroid central field - inner z
203  const double m_toroidZ5 = 9800.0*Gaudi::Units::mm;
204  // endcap toroid central field - outer z
205  const double m_toroidZ6 = 11400.0*Gaudi::Units::mm;
206  // toroid exit fringe fields
207  const double m_toroidZ7 = 12900.0*Gaudi::Units::mm;
208  // essentially out of any toroid fields
209  const double m_toroidZ8 = 14000.0*Gaudi::Units::mm;
210 
211  // counters
212  mutable std::atomic<unsigned long long> m_countExtrapolations = 0;
213  mutable std::atomic<unsigned long long> m_countShortStep = 0;
214  mutable std::atomic<unsigned long long> m_countStep = 0;
215  mutable std::atomic<unsigned long long> m_countStepReduction = 0;
216 };
217 
218 //<<<<<< INLINE PRIVATE MEMBER FUNCTIONS >>>>>>
219 
221  const TrackSurfaceIntersection& isect, const double cylinderRadius,
222  const double offsetRadius, const Amg::Vector3D& offset,
223  double& stepLength) const {
224  const Amg::Vector3D& dir = isect.direction();
225 
226  double sinThsqinv = 1. / dir.perp2();
227  stepLength = (offset.x() * dir.x() + offset.y() * dir.y()) * sinThsqinv;
228  double deltaRSq = (cylinderRadius - offsetRadius) *
229  (cylinderRadius + offsetRadius) * sinThsqinv +
230  stepLength * stepLength;
231  if (deltaRSq > 0.)
232  stepLength += sqrt(deltaRSq);
233  return std::abs(stepLength);
234 }
235 
237  const TrackSurfaceIntersection& isect, const double discZ,
238  double& stepLength) const {
239  const Amg::Vector3D& pos = isect.position();
240  const Amg::Vector3D& dir = isect.direction();
241 
242  double distance = discZ - pos.z();
243  stepLength = distance / dir.z();
244  return std::abs(distance);
245 }
246 
248  const TrackSurfaceIntersection& isect, const Amg::Vector3D& linePosition,
249  const Amg::Vector3D& lineDirection, double& stepLength) const {
250  // offset joining track to line is given by
251  // offset = linePosition + a*lineDirection - trackPosition -
252  // b*trackDirection
253  //
254  // offset is perpendicular to both line and track at solution i.e.
255  // lineDirection.dot(offset) = 0
256  // trackDirection.dot(offset) = 0
257  const Amg::Vector3D& pos = isect.position();
258  const Amg::Vector3D& dir = isect.direction();
259 
260  double cosAngle = lineDirection.dot(dir);
261  stepLength = (linePosition - pos).dot(dir - lineDirection * cosAngle) /
262  (1. - cosAngle * cosAngle);
263  return std::abs(stepLength);
264 }
265 
267  const TrackSurfaceIntersection& isect, const Amg::Vector3D& planePosition,
268  const Amg::Vector3D& planeNormal, double& stepLength) const {
269  // take the normal component of the offset from track position to plane
270  // position this is equal to the normal component of the required distance
271  // along the track direction
272  const Amg::Vector3D& pos = isect.position();
273  const Amg::Vector3D& dir = isect.direction();
274 
275  double distance = planeNormal.dot(planePosition - pos);
276  stepLength = distance / planeNormal.dot(dir);
277  return std::abs(distance);
278 }
279 
281  MagField::AtlasFieldCache& fieldCache) const {
284  if (!fieldCondObj.isValid()) {
285  ATH_MSG_FATAL("Failed to get magnetic field conditions data "
287  }
288  fieldCondObj->getInitializedCache(fieldCache);
289 }
290 
292  const Amg::Vector3D& position,
293  MagField::AtlasFieldCache& fieldCache) const {
295  fieldCache.getField(position.data(), fieldValue.data());
296  return fieldValue;
297 }
298 
299 inline std::optional<TrackSurfaceIntersection> RungeKuttaIntersector::newIntersection(
300  TrackSurfaceIntersection&& isect, const Surface& surface,
301  const double qOverP, const double rStart, const double zStart) const {
302  // ensure intersection is valid (ie. on surface)
303  Intersection SLIntersect = surface.straightLineIntersection(
304  isect.position(), isect.direction(), false, false);
305  if (SLIntersect.valid) {
306  isect.position() = SLIntersect.position;
307  return std::move(isect);
308  }
309 
310  // invalid: take care to invalidate cache!
311  if (msgLvl(MSG::DEBUG)){
312  debugFailure(std::move(isect), surface, qOverP, rStart, zStart, false);
313  }
314 
315  return std::nullopt;
316 }
317 
318 } // namespace Trk
319 
320 #endif // TRKEXRUNGEKUTTAINTERSECTOR_RUNGEKUTTAINTERSECTOR_H
Trk::RungeKuttaIntersector::m_toroidZ0
const double m_toroidZ0
Definition: RungeKuttaIntersector.h:193
Trk::RungeKuttaIntersector::m_muonZ0
const double m_muonZ0
Definition: RungeKuttaIntersector.h:170
Trk::RungeKuttaIntersector::distanceToPlane
double distanceToPlane(const TrackSurfaceIntersection &isect, const Amg::Vector3D &planePosition, const Amg::Vector3D &planeNormal, double &stepLength) const
Definition: RungeKuttaIntersector.h:266
Trk::RungeKuttaIntersector::m_toroidR2
const double m_toroidR2
Definition: RungeKuttaIntersector.h:189
Trk::RungeKuttaIntersector
Definition: RungeKuttaIntersector.h:28
Trk::RungeKuttaIntersector::m_momentumWarnThreshold
const double m_momentumWarnThreshold
Definition: RungeKuttaIntersector.h:166
ATH_MSG_FATAL
#define ATH_MSG_FATAL(x)
Definition: AthMsgStreamMacros.h:34
Trk::RungeKuttaIntersector::m_countStepReduction
std::atomic< unsigned long long > m_countStepReduction
Definition: RungeKuttaIntersector.h:215
Trk::Intersection
Definition: Intersection.h:24
Trk::RungeKuttaIntersector::m_inDetR1
const double m_inDetR1
Definition: RungeKuttaIntersector.h:154
SG::ReadCondHandle
Definition: ReadCondHandle.h:44
Trk::RungeKuttaIntersector::m_third
const double m_third
Definition: RungeKuttaIntersector.h:183
Trk::Surface::straightLineIntersection
Intersection straightLineIntersection(const T &pars, bool forceDir=false, const Trk::BoundaryCheck &bchk=false) const
fst straight line intersection schema - templated for charged and neutral parameters
Definition: Tracking/TrkDetDescr/TrkSurfaces/TrkSurfaces/Surface.h:351
Trk::RungeKuttaIntersector::m_stepMax1
double m_stepMax1
Definition: RungeKuttaIntersector.h:178
Trk::PerigeeSurface
Definition: PerigeeSurface.h:43
Trk::RungeKuttaIntersector::m_caloR2
const double m_caloR2
Definition: RungeKuttaIntersector.h:138
AtlasFieldCacheCondObj.h
Trk::RungeKuttaIntersector::m_countExtrapolations
std::atomic< unsigned long long > m_countExtrapolations
Definition: RungeKuttaIntersector.h:212
Trk::RungeKuttaIntersector::m_productionMode
BooleanProperty m_productionMode
Definition: RungeKuttaIntersector.h:130
Trk::RungeKuttaIntersector::finalize
virtual StatusCode finalize() override
Definition: RungeKuttaIntersector.cxx:79
Trk::RungeKuttaIntersector::m_toroidZ6
const double m_toroidZ6
Definition: RungeKuttaIntersector.h:205
Trk::RungeKuttaIntersector::~RungeKuttaIntersector
virtual ~RungeKuttaIntersector()=default
Trk::RungeKuttaIntersector::intersectPlaneSurface
virtual std::optional< TrackSurfaceIntersection > intersectPlaneSurface(const PlaneSurface &surface, const TrackSurfaceIntersection &trackIntersection, const double qOverP) const override
IIntersector interface method for specific Surface type : PlaneSurface.
Definition: RungeKuttaIntersector.cxx:339
Trk::RungeKuttaIntersector::m_inDetZ2
const double m_inDetZ2
Definition: RungeKuttaIntersector.h:162
Trk::RungeKuttaIntersector::m_toroidZ8
const double m_toroidZ8
Definition: RungeKuttaIntersector.h:209
python.SystemOfUnits.MeV
int MeV
Definition: SystemOfUnits.py:154
Trk::RungeKuttaIntersector::initialize
virtual StatusCode initialize() override
Definition: RungeKuttaIntersector.cxx:56
Trk::RungeKuttaIntersector::distanceToCylinder
double distanceToCylinder(const TrackSurfaceIntersection &isect, const double cylinderRadius, const double offsetRadius, const Amg::Vector3D &offset, double &stepLength) const
Definition: RungeKuttaIntersector.h:220
Trk::RungeKuttaIntersector::m_stepsUntilTrapped
int m_stepsUntilTrapped
Definition: RungeKuttaIntersector.h:182
Trk::RungeKuttaIntersector::m_stepMax4
double m_stepMax4
Definition: RungeKuttaIntersector.h:181
Trk::RungeKuttaIntersector::m_toroidR3
const double m_toroidR3
Definition: RungeKuttaIntersector.h:191
Trk::RungeKuttaIntersector::RungeKuttaIntersector
RungeKuttaIntersector(const std::string &type, const std::string &name, const IInterface *parent)
Definition: RungeKuttaIntersector.cxx:50
Trk::RungeKuttaIntersector::m_shortStepMax
double m_shortStepMax
Definition: RungeKuttaIntersector.h:171
Trk::TrackSurfaceIntersection
Definition: TrackSurfaceIntersection.h:32
Trk::RungeKuttaIntersector::field
Amg::Vector3D field(const Amg::Vector3D &point, MagField::AtlasFieldCache &fieldCache) const
Definition: RungeKuttaIntersector.h:291
Trk::DiscSurface
Definition: DiscSurface.h:54
SG::VarHandleKey::key
const std::string & key() const
Return the StoreGate ID for the referenced object.
Definition: AthToolSupport/AsgDataHandles/Root/VarHandleKey.cxx:141
python.CaloAddPedShiftConfig.type
type
Definition: CaloAddPedShiftConfig.py:42
Trk::RungeKuttaIntersector::m_caloR4
const double m_caloR4
Definition: RungeKuttaIntersector.h:142
ReadCondHandle.h
Trk::RungeKuttaIntersector::m_caloR0
const double m_caloR0
Definition: RungeKuttaIntersector.h:134
Trk::RungeKuttaIntersector::m_stepMax2
double m_stepMax2
Definition: RungeKuttaIntersector.h:179
Trk::RungeKuttaIntersector::m_momentumThreshold
const double m_momentumThreshold
Definition: RungeKuttaIntersector.h:164
Trk::RungeKuttaIntersector::m_caloR1
const double m_caloR1
Definition: RungeKuttaIntersector.h:136
GeoPrimitives.h
Trk::RungeKuttaIntersector::m_stepMax3
double m_stepMax3
Definition: RungeKuttaIntersector.h:180
Trk::RungeKuttaIntersector::m_toroidZ3
const double m_toroidZ3
Definition: RungeKuttaIntersector.h:199
Trk::RungeKuttaIntersector::m_toroidZ7
const double m_toroidZ7
Definition: RungeKuttaIntersector.h:207
Trk::RungeKuttaIntersector::m_solenoidZ
const double m_solenoidZ
Definition: RungeKuttaIntersector.h:176
Trk::RungeKuttaIntersector::m_caloR3
const double m_caloR3
Definition: RungeKuttaIntersector.h:140
Trk::RungeKuttaIntersector::distanceToDisc
double distanceToDisc(const TrackSurfaceIntersection &isect, const double discZ, double &stepLength) const
Definition: RungeKuttaIntersector.h:236
Trk::TrackSurfaceIntersection::position
const Amg::Vector3D & position() const
Method to retrieve the position of the Intersection.
Definition: TrackSurfaceIntersection.h:80
Trk::RungeKuttaIntersector::m_toroidR1
const double m_toroidR1
Definition: RungeKuttaIntersector.h:187
Trk::RungeKuttaIntersector::m_caloZ2
const double m_caloZ2
Definition: RungeKuttaIntersector.h:148
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
Trk::CylinderSurface
Definition: CylinderSurface.h:55
Trk::RungeKuttaIntersector::step
void step(TrackSurfaceIntersection &isect, Amg::Vector3D &fieldValue, double &stepLength, const double qOverP, MagField::AtlasFieldCache &fieldCache) const
Definition: RungeKuttaIntersector.cxx:711
Trk::RungeKuttaIntersector::shortStep
void shortStep(TrackSurfaceIntersection &isect, const Amg::Vector3D &fieldValue, const double stepLength, const double qOverP) const
Definition: RungeKuttaIntersector.cxx:684
AthAlgTool.h
test_pyathena.parent
parent
Definition: test_pyathena.py:15
Trk::RungeKuttaIntersector::distanceToLine
double distanceToLine(const TrackSurfaceIntersection &isect, const Amg::Vector3D &linePosition, const Amg::Vector3D &lineDirection, double &stepLength) const
Definition: RungeKuttaIntersector.h:247
Trk::RungeKuttaIntersector::initializeFieldCache
void initializeFieldCache(MagField::AtlasFieldCache &fieldCache) const
Definition: RungeKuttaIntersector.h:280
Trk::RungeKuttaIntersector::m_caloZ1
const double m_caloZ1
Definition: RungeKuttaIntersector.h:146
Trk::RungeKuttaIntersector::m_countShortStep
std::atomic< unsigned long long > m_countShortStep
Definition: RungeKuttaIntersector.h:213
Trk::Intersection::position
Amg::Vector3D position
Definition: Intersection.h:25
python.SystemOfUnits.nanometer
int nanometer
Definition: SystemOfUnits.py:72
dot.dot
def dot(G, fn, nodesToHighlight=[])
Definition: dot.py:5
Trk::RungeKuttaIntersector::m_toroidZ5
const double m_toroidZ5
Definition: RungeKuttaIntersector.h:203
beamspotman.dir
string dir
Definition: beamspotman.py:623
ReadCondHandleKey.h
EventPrimitives.h
Trk
Ensure that the ATLAS eigen extensions are properly loaded.
Definition: FakeTrackBuilder.h:9
Trk::RungeKuttaIntersector::m_inDetZ1
const double m_inDetZ1
Definition: RungeKuttaIntersector.h:160
Trk::RungeKuttaIntersector::isValid
virtual bool isValid(Amg::Vector3D, Amg::Vector3D) const override
IIntersector interface method for validity check over a particular extrapolation range.
Definition: RungeKuttaIntersector.h:77
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:228
Trk::RungeKuttaIntersector::m_muonR0
const double m_muonR0
Definition: RungeKuttaIntersector.h:168
columnar::final
CM final
Definition: ColumnAccessor.h:106
Trk::RungeKuttaIntersector::assignStepLength
void assignStepLength(const TrackSurfaceIntersection &isect, double &stepLength) const
Definition: RungeKuttaIntersector.cxx:382
Trk::RungeKuttaIntersector::m_fieldCacheCondObjInputKey
SG::ReadCondHandleKey< AtlasFieldCacheCondObj > m_fieldCacheCondObjInputKey
Definition: RungeKuttaIntersector.h:125
Trk::RungeKuttaIntersector::m_solenoidR
const double m_solenoidR
Definition: RungeKuttaIntersector.h:174
Amg::Vector3D
Eigen::Matrix< double, 3, 1 > Vector3D
Definition: GeoPrimitives.h:47
IIntersector.h
Trk::TrackSurfaceIntersection::direction
const Amg::Vector3D & direction() const
Method to retrieve the direction at the Intersection.
Definition: TrackSurfaceIntersection.h:88
Trk::Intersection::valid
bool valid
Definition: Intersection.h:28
python.LumiBlobConversion.pos
pos
Definition: LumiBlobConversion.py:18
SG::ReadCondHandleKey< AtlasFieldCacheCondObj >
python.SystemOfUnits.mm
int mm
Definition: SystemOfUnits.py:83
Trk::RungeKuttaIntersector::m_inDetZ0
const double m_inDetZ0
Definition: RungeKuttaIntersector.h:158
taskman.fieldValue
fieldValue
Definition: taskman.py:493
Trk::RungeKuttaIntersector::intersectCylinderSurface
virtual std::optional< TrackSurfaceIntersection > intersectCylinderSurface(const CylinderSurface &surface, const TrackSurfaceIntersection &trackIntersection, const double qOverP) const override
IIntersector interface method for specific Surface type : CylinderSurface.
Definition: RungeKuttaIntersector.cxx:226
Trk::PlaneSurface
Definition: PlaneSurface.h:64
Trk::RungeKuttaIntersector::m_toroidR0
const double m_toroidR0
Definition: RungeKuttaIntersector.h:185
Trk::RungeKuttaIntersector::m_caloZ3
const double m_caloZ3
Definition: RungeKuttaIntersector.h:150
DEBUG
#define DEBUG
Definition: page_access.h:11
MagField::AtlasFieldCache
Local cache for magnetic field (based on MagFieldServices/AtlasFieldSvcTLS.h)
Definition: AtlasFieldCache.h:43
Trk::qOverP
@ qOverP
perigee
Definition: ParamDefs.h:67
Trk::RungeKuttaIntersector::m_inDetR0
const double m_inDetR0
Definition: RungeKuttaIntersector.h:152
convertTimingResiduals.offset
offset
Definition: convertTimingResiduals.py:71
Trk::RungeKuttaIntersector::m_caloZ0
const double m_caloZ0
Definition: RungeKuttaIntersector.h:144
MagField::AtlasFieldCache::getField
void getField(const double *ATH_RESTRICT xyz, double *ATH_RESTRICT bxyz, double *ATH_RESTRICT deriv=nullptr)
get B field value at given position xyz[3] is in mm, bxyz[3] is in kT if deriv[9] is given,...
Definition: AtlasFieldCache.cxx:42
Trk::RungeKuttaIntersector::m_toroidZ4
const double m_toroidZ4
Definition: RungeKuttaIntersector.h:201
Trk::RungeKuttaIntersector::m_stepMax0
double m_stepMax0
Definition: RungeKuttaIntersector.h:177
Trk::RungeKuttaIntersector::intersectSurface
virtual std::optional< TrackSurfaceIntersection > intersectSurface(const Surface &surface, const TrackSurfaceIntersection &trackIntersection, const double qOverP) const override
IIntersector interface method for general Surface type.
Definition: RungeKuttaIntersector.cxx:101
Trk::RungeKuttaIntersector::approachPerigeeSurface
virtual std::optional< TrackSurfaceIntersection > approachPerigeeSurface(const PerigeeSurface &surface, const TrackSurfaceIntersection &trackIntersection, const double qOverP) const override
IIntersector interface method for specific Surface type : PerigeeSurface.
Definition: RungeKuttaIntersector.cxx:142
Trk::RungeKuttaIntersector::m_countStep
std::atomic< unsigned long long > m_countStep
Definition: RungeKuttaIntersector.h:214
Trk::RungeKuttaIntersector::approachStraightLineSurface
virtual std::optional< TrackSurfaceIntersection > approachStraightLineSurface(const StraightLineSurface &surface, const TrackSurfaceIntersection &trackIntersection, const double qOverP) const override
IIntersector interface method for specific Surface type : StraightLineSurface.
Definition: RungeKuttaIntersector.cxx:184
Trk::RungeKuttaIntersector::m_inDetR2
const double m_inDetR2
Definition: RungeKuttaIntersector.h:156
Trk::RungeKuttaIntersector::intersectDiscSurface
virtual std::optional< TrackSurfaceIntersection > intersectDiscSurface(const DiscSurface &surface, const TrackSurfaceIntersection &trackIntersection, const double qOverP) const override
IIntersector interface method for specific Surface type : DiscSurface.
Definition: RungeKuttaIntersector.cxx:298
TrackSurfaceIntersection.h
Trk::Surface
Definition: Tracking/TrkDetDescr/TrkSurfaces/TrkSurfaces/Surface.h:75
Trk::RungeKuttaIntersector::newIntersection
std::optional< TrackSurfaceIntersection > newIntersection(TrackSurfaceIntersection &&isect, const Surface &surface, const double qOverP, const double rStart, const double zStart) const
Definition: RungeKuttaIntersector.h:299
Amg::distance
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
Definition: GeoPrimitivesHelpers.h:54
Trk::RungeKuttaIntersector::m_toroidZ2
const double m_toroidZ2
Definition: RungeKuttaIntersector.h:197
Trk::RungeKuttaIntersector::m_shortStepMin
const double m_shortStepMin
Definition: RungeKuttaIntersector.h:172
Trk::RungeKuttaIntersector::m_toroidZ1
const double m_toroidZ1
Definition: RungeKuttaIntersector.h:195
Trk::StraightLineSurface
Definition: StraightLineSurface.h:51
Trk::RungeKuttaIntersector::debugFailure
void debugFailure(TrackSurfaceIntersection &&isect, const Surface &surface, const double qOverP, const double rStart, const double zStart, const bool trapped) const
Definition: RungeKuttaIntersector.cxx:585