ATLAS Offline Software
Loading...
Searching...
No Matches
MuidTrackIsolation.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3*/
4
6// MuidTrackIsolation
7// AlgTool for estimating the number, total charged momentum and most
8// energetic inner detector tracks in a cone surrounding a muon
9//
11
12#include "MuidTrackIsolation.h"
13
14#include <cmath>
15#include <iomanip>
16
18#include "GaudiKernel/SystemOfUnits.h"
23#include "TrkSurfaces/Surface.h"
24#include "TrkTrack/Track.h"
26namespace Rec {
27
28 MuidTrackIsolation::MuidTrackIsolation(const std::string& type, const std::string& name, const IInterface* parent) :
29 AthAlgTool(type, name, parent), m_etaSafetyFactor(0.1) {
30 declareInterface<IMuidTrackIsolation>(this);
31 }
32
33 //<<<<<< PUBLIC MEMBER FUNCTION DEFINITIONS >>>>>>
34
36 ATH_MSG_INFO("MuidTrackIsolation::initialize()");
37
38 // get the Tools
39 ATH_CHECK(m_intersector.retrieve());
41 // create the calo barrel surfaces (cylinder) and 2 endcap discs)
42 double radius = 2.0 * Gaudi::Units::meter;
43 double halfLength = 4.0 * Gaudi::Units::meter;
44 Amg::Transform3D transform;
45 transform.setIdentity();
46 m_caloCylinder = std::make_unique<Trk::CylinderSurface>(transform, radius, halfLength);
47
48 // the corresponding max barrel cotTheta
49 m_barrelCotTheta = halfLength / radius;
50
51 // and the forward/backward endcap disks
52 Amg::Transform3D discRotation;
53 discRotation.setIdentity();
54 Amg::Vector3D forwardDiscPosition(0., 0., halfLength);
55 auto transform1 = std::make_unique<Amg::Transform3D>(discRotation * forwardDiscPosition);
56 m_caloForwardDisc = std::make_unique<Trk::DiscSurface>(*transform1, 0., radius);
57 Amg::Vector3D backwardDiscPosition(0., 0., -halfLength);
58 auto transform2 = std::make_unique<Amg::Transform3D>(discRotation * backwardDiscPosition);
59 m_caloBackwardDisc = std::make_unique<Trk::DiscSurface>(*transform2, 0., radius);
60
61 ATH_CHECK(m_inDetTracksLocation.initialize());
62
63 return StatusCode::SUCCESS;
64 }
65 std::pair<int, double> MuidTrackIsolation::trackIsolation(const EventContext& ctx, double eta, double phi) const {
66 // debug input quantities
67 ATH_MSG_DEBUG(" MuidTrackIsolation:: " << std::setiosflags(std::ios::fixed)
68 << (m_trackExtrapolation ? "applied after extrapolation to calo, " : "applied at perigee, ")
69 << " for muon at calo with eta,phi " << std::setw(8) << std::setprecision(3) << eta
70 << std::setw(8) << std::setprecision(3) << phi);
71
72 // set initial state
73 std::pair<int, double> isolation{0, 0.};
74
75 // retrieve track collection
77 if (!inDetTracks.isPresent()) {
78 ATH_MSG_DEBUG(" no ID Track container at location " << m_inDetTracksLocation.key());
79 return isolation;
80 }
81
82 if (!inDetTracks.isValid()) {
83 ATH_MSG_WARNING(" ID Track container " << m_inDetTracksLocation.key() << " not valid!");
84 return isolation;
85 }
86
87 // evaluate isolation according to configuration
89 isolation = trackExtrapolated(inDetTracks.cptr(), eta, phi);
90 } else {
91 isolation = trackVertex(inDetTracks.cptr(), eta, phi);
92 }
93
94 // debug result
95 ATH_MSG_DEBUG("Found " << isolation.first << std::setiosflags(std::ios::fixed) << " InDet tracks with total momentum "
96 << std::setw(8) << std::setprecision(1) << isolation.second / Gaudi::Units::GeV << " GeV");
97
98 return isolation;
99 }
100
101 std::pair<int, double> MuidTrackIsolation::trackVertex(const TrackCollection* inDetTracks, double eta, double phi) const {
102 // set initial state
103 double sumP = 0.;
104 int numberTracks = 0;
105
106 // choose tracks in cone
107 for (const Trk::Track* id : *inDetTracks) {
108 const Trk::Perigee& perigee = *id->perigeeParameters();
109 if (id->info().trackProperties(Trk::TrackInfo::StraightTrack) || perigee.pT() < m_minPt) continue;
110
111 double inDetPhi = perigee.parameters()[Trk::phi];
112 double inDetEta = perigee.eta();
113
114 double diffEta = std::abs(eta - inDetEta);
115 double diffPhi = xAOD::P4Helpers::deltaPhi(phi, inDetPhi);
116
117 ATH_MSG_DEBUG(std::endl
118 << std::setiosflags(std::ios::fixed) << " Id track: momentum " << std::setw(8) << std::setprecision(1)
119 << perigee.momentum().mag() / Gaudi::Units::GeV << " with perigee eta and difference " << std::setw(8)
120 << std::setprecision(3) << perigee.eta() << std::setw(8) << std::setprecision(3) << diffEta
121 << " and same for phi " << std::setw(8) << std::setprecision(3) << perigee.parameters()[Trk::phi] << std::setw(8)
122 << std::setprecision(3) << diffPhi);
123
124 if ((diffPhi * diffPhi + diffEta * diffEta) > m_trackCone2) continue;
125 ++numberTracks;
126 const double p = perigee.momentum().mag();
127 sumP += p;
128
129 ATH_MSG_VERBOSE("inside cone, track#" << std::setw(3) << numberTracks);
130 }
131
132 return std::make_pair(numberTracks, sumP);
133 }
134
135 std::pair<int, double> MuidTrackIsolation::trackExtrapolated(const TrackCollection* inDetTracks, double eta, double phi) const {
136 // set initial state
137 double sumP = 0.;
138 int numberTracks = 0;
139
140 // extrapolate close in eta tracks to calorimeter surface
141 for (const Trk::Track* id : *inDetTracks) {
142 const Trk::Perigee& perigee = *id->perigeeParameters();
143 if (id->info().trackProperties(Trk::TrackInfo::StraightTrack) || perigee.pT() < m_minPt) continue;
144
145 double inDetEta = perigee.eta();
146 if (std::abs(eta - inDetEta) > m_trackCone + m_etaSafetyFactor) continue;
147
148 // track has sufficient momentum and is close in eta:
149 // find intersection at calo surface
150 double qOverP = perigee.parameters()[Trk::qOverP];
152 double cotTheta = std::tan(M_PI_2 - perigee.parameters()[Trk::theta]);
153 Amg::Vector3D direction(std::cos(perigee.parameters()[Trk::phi]), std::sin(perigee.parameters()[Trk::phi]), cotTheta);
154 direction /= direction.mag();
155
156 const Trk::TrackSurfaceIntersection idIntersection(perigee.position(), direction, 0.);
157 const Trk::Surface* surface = m_caloCylinder.get();
158 if (cotTheta > m_barrelCotTheta) {
159 surface = m_caloForwardDisc.get();
160 } else if (cotTheta < -m_barrelCotTheta) {
161 surface = m_caloBackwardDisc.get();
162 }
163 std::optional<Trk::TrackSurfaceIntersection> caloIntersection(
164 m_intersector->intersectSurface(*surface, idIntersection, qOverP));
165
166 // no intersection - should never happen !
167 if (!caloIntersection) {
168 ATH_MSG_DEBUG(" track didn't find intersection !!! "
169 << std::setiosflags(std::ios::fixed) << " Id track: momentum " << std::setw(8) << std::setprecision(1)
170 << perigee.momentum().mag() / Gaudi::Units::GeV << " with initial eta " << std::setw(8)
171 << std::setprecision(3) << perigee.eta() << " and phi " << std::setw(8) << std::setprecision(3)
172 << perigee.parameters()[Trk::phi]);
173
174 continue;
175 }
176
177 double diffEta = eta - caloIntersection->position().eta();
178 double diffPhi = xAOD::P4Helpers::deltaPhi(phi, caloIntersection->position().phi());
179 ATH_MSG_VERBOSE(std::endl
180 << std::setiosflags(std::ios::fixed) << " Id track: momentum " << std::setw(8) << std::setprecision(1)
181 << perigee.momentum().mag() / Gaudi::Units::GeV << " with initial,extrapolated and calo difference for eta "
182 << std::setw(8) << std::setprecision(3) << perigee.eta() << std::setw(8) << std::setprecision(3)
183 << caloIntersection->position().eta() << std::setw(8) << std::setprecision(3) << diffEta << " and phi "
184 << std::setw(8) << std::setprecision(3) << perigee.parameters()[Trk::phi] << std::setw(8)
185 << std::setprecision(3) << caloIntersection->position().phi() << std::setw(8) << std::setprecision(3)
186 << diffPhi);
187
188 // check if inside cone
189 if ((diffPhi * diffPhi + diffEta * diffEta) < m_trackCone2) {
190 ++numberTracks;
191 const double p = perigee.momentum().mag();
192 sumP += p;
193
194 ATH_MSG_VERBOSE(" inside cone, track#" << std::setw(3) << numberTracks);
195 }
196 }
197
198 return std::make_pair(numberTracks, sumP);
199 }
200
201} // namespace Rec
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
DataVector< Trk::Track > TrackCollection
This typedef represents a collection of Trk::Track objects.
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
std::pair< int, double > trackIsolation(const EventContext &ctx, double eta, double phi) const override
IMuidTrackIsolation interface: get the number of tracks and summed momentum in a cone at the producti...
Gaudi::Property< double > m_minPt
std::unique_ptr< const Trk::Surface > m_caloForwardDisc
std::unique_ptr< const Trk::Surface > m_caloBackwardDisc
Gaudi::Property< double > m_trackCone
ToolHandle< Trk::IIntersector > m_intersector
Gaudi::Property< bool > m_trackExtrapolation
MuidTrackIsolation(const std::string &type, const std::string &name, const IInterface *parent)
std::unique_ptr< const Trk::Surface > m_caloCylinder
SG::ReadHandleKey< TrackCollection > m_inDetTracksLocation
StatusCode initialize() override
std::pair< int, double > trackExtrapolated(const TrackCollection *indetTracks, double eta, double phi) const
std::pair< int, double > trackVertex(const TrackCollection *indetTracks, double eta, double phi) const
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
bool isPresent() const
Is the referenced object present in SG?
double eta() const
Access method for pseudorapidity - from momentum.
const Amg::Vector3D & momentum() const
Access method for the momentum.
const Amg::Vector3D & position() const
Access method for the position.
double pT() const
Access method for transverse momentum.
Abstract Base Class for tracking surfaces.
An intersection with a Surface is given by.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
Gaudi Tools.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ phi
Definition ParamDefs.h:75
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[