ATLAS Offline Software
Loading...
Searching...
No Matches
ExtrapolatorTest.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6// ExtrapolatorTest.cxx, (c) ATLAS Detector software
8
9// Tracking
16// std
17#include <cmath>
18
19//================ Constructor =================================================
20
21Trk::ExtrapolatorTest::ExtrapolatorTest(const std::string& name, ISvcLocator* pSvcLocator) :
22 AthAlgorithm(name,pSvcLocator) {}
23
24//================ Destructor =================================================
25
27{
28 delete m_gaussDist;
29 delete m_flatDist;
30 // cleanup of the surfaces
31 for (const auto& surfaceTriple : m_referenceSurfaceTriples) {
32 for (const auto* surface : surfaceTriple) {
33 delete surface;
34 }
35 }
36}
37
38
39//================ Initialisation =================================================
40
42{
43 // Code entered here will be executed once at program start.
44
45 ATH_MSG_INFO(" initialize()");
46
47 ATH_CHECK( m_extrapolator.retrieve() );
48 ATH_CHECK( m_propagator.retrieve() );
49
51
53 // assign the size
55 // loop over it and create the
56 std::vector<double>::iterator radiusIter = m_referenceSurfaceRadius.begin();
57 std::vector<double>::iterator radiusIterEnd = m_referenceSurfaceRadius.end();
58 std::vector<double>::iterator halfZIter = m_referenceSurfaceHalflength.begin();
59
60 for ( ; radiusIter != radiusIterEnd; ++radiusIter, ++halfZIter){
61 // create the Surface triplet
62 std::vector< const Trk::Surface*> surfaceTriplet;
63 surfaceTriplet.push_back(new Trk::DiscSurface(Amg::Transform3D(Amg::Translation3D(0.,0.,*halfZIter)),0.,*radiusIter));
64 surfaceTriplet.push_back(new Trk::CylinderSurface(Amg::Transform3D(),*radiusIter, *halfZIter));
65 surfaceTriplet.push_back(new Trk::DiscSurface(Amg::Transform3D(Amg::Translation3D(0.,0.,-(*halfZIter))),0.,*radiusIter));
66
67 ATH_MSG_INFO("Creating surfaces: R " << *radiusIter << " Z " << *halfZIter);
68 m_referenceSurfaceTriples.push_back(surfaceTriplet);
69
70 m_referenceSurfaceNegativeBoundary.push_back(atan2(*radiusIter,-(*halfZIter)));
71 m_referenceSurfacePositiveBoundary.push_back(atan2(*radiusIter,(*halfZIter)));
72
73 }
74 }
75
76 m_gaussDist = new Rndm::Numbers(randSvc(), Rndm::Gauss(0.,1.));
77 m_flatDist = new Rndm::Numbers(randSvc(), Rndm::Flat(0.,1.));
78
79 msg(MSG::INFO) << "initialize() successful in " << endmsg;
80
81 if( m_eventsPerExecute > 0 ){
83 for( int i=0;i<m_eventsPerExecute;++i ) m_perigees.push_back(generatePerigee());
84 }
85
86 return StatusCode::SUCCESS;
87}
88
89//================ Finalisation =================================================
90
92{
93 // Code entered here will be executed once at the end of the program run.
94 return StatusCode::SUCCESS;
95}
96
98 // generate with random number generator
99 double d0 = m_gaussDist->shoot() * m_sigmaD0;
100 double z0 = m_gaussDist->shoot() * m_sigmaZ0;
101 double phi = m_minPhi + (m_maxPhi-m_minPhi)* m_flatDist->shoot();
102 double eta = m_minEta + m_flatDist->shoot()*(m_maxEta-m_minEta);
103 double theta = 2.*atan(exp(-eta));
104 double p = m_minP + m_flatDist->shoot()*(m_maxP-m_minP);
105 double charge = (m_flatDist->shoot() > 0.5 ) ? -1. : 1.;
106 double qOverP = charge/(p);
107 const Trk::PerigeeSurface perSurface;
108 return {d0, z0, phi, theta, qOverP,perSurface};
109}
110
111//================ Execution ====================================================
112
114{
115
117 else for( int i=0;i<m_eventsPerExecute;++i ) runTest(m_perigees[i]);
118 return StatusCode::SUCCESS;
119}
120
121
122void Trk::ExtrapolatorTest::runTest( const Trk::Perigee& initialPerigee ) {
124
125 ATH_MSG_VERBOSE("Starting from : " << initialPerigee );
126 ATH_MSG_VERBOSE(" ---> with direction: " << propagationDirection );
127
128
129
130 std::vector< std::vector< const Trk::Surface* > >::const_iterator surfaceTripleIter = m_referenceSurfaceTriples.begin();
131 std::vector< std::vector< const Trk::Surface* > >::const_iterator surfaceTripleIterEnd = m_referenceSurfaceTriples.end();
132
133 std::vector<double>::iterator negRefIter = m_referenceSurfaceNegativeBoundary.begin();
134 std::vector<double>::iterator posRefIter = m_referenceSurfacePositiveBoundary.begin();
135
136 double theta = initialPerigee.parameters()[Trk::theta];
137
138 const EventContext& ctx = Gaudi::Hive::currentContext();
139 for (int refSurface = 0 ; surfaceTripleIter != surfaceTripleIterEnd; ++surfaceTripleIter, ++negRefIter, ++posRefIter ){
140 // decide which reference surface to take
141 refSurface = theta < (*posRefIter) ? 2 : 1;
142 refSurface = theta > (*negRefIter) ? 0 : 1;
143
144 const Trk::Surface* destinationSurface = (*surfaceTripleIter)[refSurface];
145
146 const auto* destParameters =
148 ? m_extrapolator->extrapolate(
149 ctx,
150 initialPerigee,
151 *destinationSurface,
152 propagationDirection,
153 false,
154 static_cast<Trk::ParticleHypothesis>(m_particleType.value())).release()
155 :
156
158 ->propagate(ctx,
159 initialPerigee,
160 *destinationSurface,
161 propagationDirection,
162 false,
164 static_cast<Trk::ParticleHypothesis>(m_particleType.value()))
165 .release();
166
167 if (destParameters) {
168 // intersection output
169 ATH_MSG_VERBOSE(" [ intersection ] with surface at (x,y,z) = " <<
170 destParameters->position().x() << ", " <<
171 destParameters->position().y() << ", " <<
172 destParameters->position().z() );
173 } else if (!destParameters)
174 ATH_MSG_DEBUG(" Extrapolation not successful! " );
175
176 delete destParameters;
177 }
178}
179
180//============================================================================================
181
Scalar eta() const
pseudorapidity method
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
double charge(const T &p)
Definition AtlasPID.h:997
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
Class for a CylinderSurface in the ATLAS detector.
Class for a DiscSurface in the ATLAS detector.
Definition DiscSurface.h:54
Rndm::Numbers * m_gaussDist
Random Number setup.
ToolHandle< IExtrapolator > m_extrapolator
The Extrapolator to be retrieved.
StatusCode initialize()
standard Athena-Algorithm method
ExtrapolatorTest(const std::string &name, ISvcLocator *pSvcLocator)
Standard Athena-Algorithm Constructor.
DoubleArrayProperty m_referenceSurfaceRadius
~ExtrapolatorTest()
Default Destructor.
IntegerProperty m_particleType
MagneticFieldProperties * m_magFieldProperties
magnetic field properties
std::vector< double > m_referenceSurfaceNegativeBoundary
unsigned int m_referenceSurfaces
member variables for algorithm properties:
PublicToolHandle< IPropagator > m_propagator
Trk::Perigee generatePerigee()
IntegerProperty m_eventsPerExecute
IntegerProperty m_direction
std::vector< Trk::Perigee > m_perigees
BooleanProperty m_useExtrapolator
std::vector< std::vector< const Surface * > > m_referenceSurfaceTriples
Rndm::Numbers * m_flatDist
StatusCode finalize()
standard Athena-Algorithm method
DoubleArrayProperty m_referenceSurfaceHalflength
void runTest(const Trk::Perigee &perigee)
std::vector< double > m_referenceSurfacePositiveBoundary
StatusCode execute()
standard Athena-Algorithm method
magnetic field properties to steer the behavior of the extrapolation
Class describing the Line to which the Perigee refers to.
Abstract Base Class for tracking surfaces.
Eigen::Affine3d Transform3D
Eigen::Translation< double, 3 > Translation3D
PropDirection
PropDirection, enum for direction of the propagation.
@ oppositeMomentum
@ alongMomentum
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ phi
Definition ParamDefs.h:75
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
ParticleHypothesis
Enumeration for Particle hypothesis respecting the interaction with material.
MsgStream & msg
Definition testRead.cxx:32