ATLAS Offline Software
Loading...
Searching...
No Matches
EMTrackMatchBuilder Class Referencefinal

EMTrackMatch data object builder. More...

#include <EMTrackMatchBuilder.h>

Inheritance diagram for EMTrackMatchBuilder:
Collaboration diagram for EMTrackMatchBuilder:

Classes

struct  TrackMatch
 A structure for keeping track match information. More...
class  TrackMatchSorter
 function object to sort track matches based on quality More...

Public Member Functions

 EMTrackMatchBuilder (const std::string &type, const std::string &name, const IInterface *parent)
 Default constructor.
virtual ~EMTrackMatchBuilder ()=default
 Destructor.
StatusCode initialize () override final
 Gaudi algorithm hooks.
virtual StatusCode executeRec (const EventContext &ctx, EgammaRecContainer *egammas) const override final
 execute method
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysInitialize () override
 Perform system initialization for an algorithm.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Static Public Member Functions

static const InterfaceID & interfaceID ()
 AlgTool interface methods.

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

StatusCode trackExecute (const EventContext &ctx, egammaRec *eg, const xAOD::TrackParticleContainer *trackPC, const CaloDetDescrManager &caloDD) const
 execute method
bool inBroadWindow (const EventContext &ctx, std::vector< TrackMatch > &trackMatches, const xAOD::CaloCluster &cluster, int trackNumber, const xAOD::TrackParticle &trkPB, const CaloDetDescrManager &caloDD) const
 Compute for tracks passing the loose matching the distance between track extrapolated to 2nd sampling and cluster.
bool isCandidateMatch (const xAOD::CaloCluster *cluster, const xAOD::TrackParticle *track, bool flip) const
 Loose track-cluster matching.
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

SG::ReadHandleKey< xAOD::TrackParticleContainerm_TrackParticlesKey
 name of TrackParticle container in TDS
SG::ReadCondHandleKey< CaloDetDescrManagerm_caloDetDescrMgrKey
Gaudi::Property< double > m_broadDeltaEta
 broad cut on deltaEta
Gaudi::Property< double > m_broadDeltaPhi
 broad cut on deltaPhi
Gaudi::Property< double > m_narrowDeltaEta
 narrow cut on deltaEta
Gaudi::Property< double > m_narrowDeltaPhi
 narrow cut on deltaPhiRescale
Gaudi::Property< double > m_narrowDeltaPhiBrem
 narrow cut on deltaPhi for electrons
Gaudi::Property< double > m_narrowDeltaPhiRescale
 narrow cut on deltaPhiRescale
Gaudi::Property< double > m_narrowDeltaPhiRescaleBrem
 narrow cut on deltaPhiRescale for electrons
Gaudi::Property< double > m_MaxDeltaPhiRescale
 @Maximum deltaPhi (Res) allowed for a match
Gaudi::Property< bool > m_useCandidateMatch
 flag to turn on/off use of isCandidateMatch
Gaudi::Property< bool > m_useScoring
 Boolean to apply heuristic when tracks have close deltaR.
Gaudi::Property< bool > m_useRescaleMetric
 Boolean to use Rescale in the metric.
Gaudi::Property< bool > m_SecondPassRescale
 Boolean to do second pass with Rescale.
Gaudi::Property< float > m_deltaEtaResolution
 The resolutions: might be good to split in barrel/end-cap in the future.
Gaudi::Property< float > m_deltaPhiResolution
Gaudi::Property< float > m_deltaPhiRescaleResolution
Gaudi::Property< float > m_distanceForScore
 The distance from which one goes from using better deltaR to using score.
ToolHandle< IEMExtrapolationToolsm_extrapolationTool
Gaudi::Property< bool > m_isCosmics
double m_deltaEtaWeight {}
double m_deltaPhiWeight {}
double m_deltaPhiRescaleWeight {}
TrackMatchSorter m_sorter
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

EMTrackMatch data object builder.

Author
H. Ma
RD schaffer
Thomas Koffas
Christos Anastopoulos The matching of a track to a cluster is driven by the EMTrackMatchBuilder tool located in the Reconstruction/egamma/egammaTools package.

Definition at line 39 of file EMTrackMatchBuilder.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< AlgTool > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ EMTrackMatchBuilder()

EMTrackMatchBuilder::EMTrackMatchBuilder ( const std::string & type,
const std::string & name,
const IInterface * parent )

Default constructor.

Definition at line 29 of file EMTrackMatchBuilder.cxx.

32 : AthAlgTool(type, name, parent)
33{
34 // declare interface
35 declareInterface<IEMTrackMatchBuilder>(this);
36}
AthAlgTool()
Default constructor:

◆ ~EMTrackMatchBuilder()

virtual EMTrackMatchBuilder::~EMTrackMatchBuilder ( )
virtualdefault

Destructor.

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ executeRec()

StatusCode EMTrackMatchBuilder::executeRec ( const EventContext & ctx,
EgammaRecContainer * egammas ) const
finaloverridevirtual

execute method

Implements IEMTrackMatchBuilder.

Definition at line 56 of file EMTrackMatchBuilder.cxx.

58{
59 // protection against bad pointers
60 if (egammas == nullptr) {
61 return StatusCode::SUCCESS;
62 }
63 // retrieve the trackparticle container
64 SG::ReadHandle<xAOD::TrackParticleContainer> trackPC(m_TrackParticlesKey,
65 ctx);
66
67 SG::ReadCondHandle<CaloDetDescrManager> caloDetDescrMgrHandle{
69 };
70 ATH_CHECK(caloDetDescrMgrHandle.isValid());
71
72 const CaloDetDescrManager* caloDD = *caloDetDescrMgrHandle;
73
74 // check is only used for serial running; remove when MT scheduler used
75 ATH_CHECK(trackPC.isValid());
76 // Loop over calling the trackExecute method
77 for (egammaRec* eg : *egammas) {
78 // retrieve the cluster
79 ATH_CHECK(trackExecute(ctx, eg, trackPC.cptr(), *caloDD));
80 }
81 return StatusCode::SUCCESS;
82}
#define ATH_CHECK
Evaluate an expression and check for errors.
SG::ReadCondHandleKey< CaloDetDescrManager > m_caloDetDescrMgrKey
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_TrackParticlesKey
name of TrackParticle container in TDS
StatusCode trackExecute(const EventContext &ctx, egammaRec *eg, const xAOD::TrackParticleContainer *trackPC, const CaloDetDescrManager &caloDD) const
execute method

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ inBroadWindow()

bool EMTrackMatchBuilder::inBroadWindow ( const EventContext & ctx,
std::vector< TrackMatch > & trackMatches,
const xAOD::CaloCluster & cluster,
int trackNumber,
const xAOD::TrackParticle & trkPB,
const CaloDetDescrManager & caloDD ) const
private

Compute for tracks passing the loose matching the distance between track extrapolated to 2nd sampling and cluster.

Definition at line 152 of file EMTrackMatchBuilder.cxx.

158{
159
162
163 // Now get the delta eta/phi and eta correction at the calorimeter
164 // final arrays that we will write
165 // Save the value of deltaPhiRescale. If we do not use rescaled
166 // perigee, we recalculate deltaPhi using rescaled momentum. This
167 // will be saved in EMTrackMatch
168 std::array<double, 4> eta = { -999.0, -999.0, -999.0, -999.0 };
169 std::array<double, 4> phi = { -999.0, -999.0, -999.0, -999.0 };
170 std::array<double, 4> deltaEta = { -999.0, -999.0, -999.0, -999.0 };
171 std::array<double, 4> deltaPhi = { -999.0, -999.0, -999.0, -999.0 };
172
173 /*
174 * Try both from perigee
175 * and from perigee Rescale.
176 *
177 * We need anyhow both to be there at the end.
178 */
179 std::pair<std::vector<CaloSampling::CaloSample>,
180 std::vector<std::unique_ptr<Trk::Surface>>>
181 layersAndSurfaces =
182 m_extrapolationTool->getClusterLayerSurfaces(cluster, caloDD);
184 ->getMatchAtCalo(ctx,
185 cluster,
186 trkPB,
187 layersAndSurfaces.first,
188 layersAndSurfaces.second,
189 eta,
190 phi,
191 deltaEta,
192 deltaPhi,
193 extrapFrom)
194 .isFailure()) {
195 return false;
196 }
197
200 std::array<double, 4> etaRes = { -999.0, -999.0, -999.0, -999.0 };
201 std::array<double, 4> phiRes = { -999.0, -999.0, -999.0, -999.0 };
202 std::array<double, 4> deltaEtaRes = { -999.0, -999.0, -999.0, -999.0 };
203 std::array<double, 4> deltaPhiRes = { -999.0, -999.0, -999.0, -999.0 };
204
206 ->getMatchAtCalo(ctx,
207 cluster,
208 trkPB,
209 layersAndSurfaces.first,
210 layersAndSurfaces.second,
211 etaRes,
212 phiRes,
213 deltaEtaRes,
214 deltaPhiRes,
215 extrapFromRes)
216 .isFailure()) {
217 return false;
218 }
219
220 double deltaPhiRescale = deltaPhiRes[2];
221 /*
222 * Sanity check for very far away matches
223 * The assumption is when we rescale we should be in the
224 * correct neighborhood for a valid track-cluster pair.
225 */
226 if (std::abs(deltaPhiRes[2]) > m_MaxDeltaPhiRescale) {
227 ATH_MSG_DEBUG("DeltaPhiRescaled above maximum: "
228 << deltaPhiRes[2] << " (max: " << m_MaxDeltaPhiRescale
229 << ")");
230 return false;
231 }
232 /*
233 * Try to match : First standard way.
234 * If this fails and the cluster Et is larger than the track Pt
235 * it might get matched only under the rescaled assumption that
236 * should be less sensitive to radiative losses.
237 */
238 if (std::abs(deltaEta[2]) < m_narrowDeltaEta && deltaPhi[2] < m_narrowDeltaPhi &&
240 ATH_MSG_DEBUG("Matched with Perigee");
241 } else if (m_SecondPassRescale && cluster.et() > trkPB.pt() &&
242 std::abs(deltaEtaRes[2]) < m_narrowDeltaEta &&
243 deltaPhiRes[2] < m_narrowDeltaPhiRescale &&
244 deltaPhiRes[2] > -m_narrowDeltaPhiRescaleBrem) {
245 ATH_MSG_DEBUG("Not Perigee but matched with Rescale");
246 } else {
247 ATH_MSG_DEBUG("Normal matched Failed deltaPhi/deltaEta "
248 << deltaPhi[2] << " / " << deltaEta[2]);
249 ATH_MSG_DEBUG("Rescaled matched Failed deltaPhi/deltaEta "
250 << deltaPhiRes[2] << " / " << deltaEtaRes[2]);
251 return false;
252 }
253
254 // Always the deltaPhiLast will be from the last measurement
257 std::array<double, 4> eta1 = { -999.0, -999.0, -999.0, -999.0 };
258 std::array<double, 4> phi1 = { -999.0, -999.0, -999.0, -999.0 };
259 std::array<double, 4> deltaEta1 = { -999.0, -999.0, -999.0, -999.0 };
260 std::array<double, 4> deltaPhi1 = { -999.0, -999.0, -999.0, -999.0 };
261
263 ->getMatchAtCalo(ctx,
264 cluster,
265 trkPB,
266 layersAndSurfaces.first,
267 layersAndSurfaces.second,
268 eta1,
269 phi1,
270 deltaEta1,
271 deltaPhi1,
272 extrapFrom1)
273 .isFailure()) {
274 ATH_MSG_DEBUG("Extrapolation from last measurement failed");
275 return false;
276 }
277 double deltaPhiLast = deltaPhi1[2];
278 ATH_MSG_DEBUG("Rescale dPhi " << deltaPhiRescale);
279 ATH_MSG_DEBUG("dPhi Last measurement " << deltaPhiLast);
280 /*
281 * Done with extrapolation
282 * Lets do the matching logic
283 */
284 TrackMatch trkmatch{};
285 // Add the matching variable to the TrackMAtch
286 trkmatch.deltaEta = deltaEta;
287 trkmatch.deltaPhi = deltaPhi;
288 trkmatch.deltaPhiRescaled = deltaPhiRes;
289 trkmatch.deltaPhiLast = deltaPhiLast;
290
291 // Variables used for the sorting. Note both dPhi's will be used.
292 trkmatch.trackNumber = trackNumber;
293 if (m_useRescaleMetric) {
294 trkmatch.dR = sqrt(std::pow(m_deltaEtaWeight * deltaEta[2], 2) +
295 std::pow(m_deltaPhiRescaleWeight * deltaPhiRescale, 2));
296 trkmatch.seconddR = sqrt(std::pow(m_deltaEtaWeight * deltaEta[2], 2) +
297 std::pow(m_deltaPhiWeight * deltaPhi[2], 2));
298 } else {
299 trkmatch.dR = sqrt(std::pow(m_deltaEtaWeight * deltaEta[2], 2) +
300 std::pow(m_deltaPhiWeight * deltaPhi[2], 2));
301 trkmatch.seconddR =
302 sqrt(std::pow(m_deltaEtaWeight * deltaEta[2], 2) +
303 std::pow(m_deltaPhiRescaleWeight * deltaPhiRescale, 2));
304 }
305 ATH_MSG_DEBUG(" DR " << trkmatch.dR << " deltaPhi " << deltaPhi[2]
306 << " deltaEta " << deltaEta[2]);
307 /*
308 * The first thing is
309 * Prefer pixel over SCT only
310 */
311 // Check number of pixel hits
312 int nPixel = summaryValueInt(trkPB, xAOD::numberOfPixelDeadSensors, 0);
313 nPixel += summaryValueInt(trkPB, xAOD::numberOfPixelHits, 0);
314 trkmatch.hasPix = (nPixel > 0);
315
316 /*
317 * Seconday hitsScore score based on hits to be used
318 * for track that are very close
319 * to each other at the calo i.e similar dR with cluster,
320 * pick the longest possible one
321 */
322 trkmatch.hitsScore = 0;
323 if (m_useScoring) {
324 // Check the 2 innermost layers
325 int nInnerMost = summaryValueInt(trkPB, xAOD::numberOfInnermostPixelLayerHits, 0);
327 int nNextToInnerMost = summaryValueInt(trkPB, xAOD::numberOfNextToInnermostPixelLayerHits, 0);
329
330 // Secondary score , find the longest track possible,
331 // i.e the one with the most inner hists in the pixel
332 // npixel*5
333 trkmatch.hitsScore += (nPixel * 5);
334 // Extra points for NextToInnermost
335 if (!expectNextToInnermostPixelLayerHit || nNextToInnerMost > 0) {
336 trkmatch.hitsScore += 5;
337 }
338 // Extra points for Innermost
339 if (!expectInnermostPixelLayerHit || nInnerMost > 0) {
340 trkmatch.hitsScore += 10;
341 }
342 }
343 ATH_MSG_DEBUG("hasPix : " << trkmatch.hasPix
344 << " hitsScore : " << trkmatch.hitsScore);
345
346 trackMatches.push_back(trkmatch);
347 return true;
348}
int summaryValueInt(const xAOD::TrackParticle &tp, const xAOD::SummaryType &info, int deflt=-999)
return the summary value for a TrackParticle or default value (-999) (to be used mostly in python whe...
Scalar eta() const
pseudorapidity method
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
Scalar phi() const
phi method
#define ATH_MSG_DEBUG(x)
Gaudi::Property< bool > m_useScoring
Boolean to apply heuristic when tracks have close deltaR.
Gaudi::Property< double > m_narrowDeltaPhi
narrow cut on deltaPhiRescale
ToolHandle< IEMExtrapolationTools > m_extrapolationTool
Gaudi::Property< double > m_MaxDeltaPhiRescale
@Maximum deltaPhi (Res) allowed for a match
Gaudi::Property< double > m_narrowDeltaPhiRescale
narrow cut on deltaPhiRescale
Gaudi::Property< double > m_narrowDeltaPhiRescaleBrem
narrow cut on deltaPhiRescale for electrons
Gaudi::Property< bool > m_useRescaleMetric
Boolean to use Rescale in the metric.
Gaudi::Property< double > m_narrowDeltaEta
narrow cut on deltaEta
Gaudi::Property< bool > m_SecondPassRescale
Boolean to do second pass with Rescale.
Gaudi::Property< double > m_narrowDeltaPhiBrem
narrow cut on deltaPhi for electrons
TrkExtrapDef
Enum for track extrapolation to calo.
@ fromPerigee
from the perigee of TrackParticle
@ fromLastMeasurement
from the last measurement of TrackParticle
@ fromPerigeeRescaled
from the perigee of TrackParticle recaled by Ecluster
virtual double pt() const override final
The transverse momentum ( ) of the particle.
double deltaEta(const I4Momentum &p1, const I4Momentum &p2)
Computes efficiently .
Definition P4Helpers.h:66
@ deltaPhiLast
deltaPhi from the las point
@ deltaPhi1
difference between the cluster eta (1st sampling) and the eta of the track extrapolated to the 1st sa...
@ deltaEta1
difference between the cluster eta (first sampling) and the eta of the track extrapolated to the firs...
setEt setPhi setE277 setWeta2 eta1
@ expectInnermostPixelLayerHit
Do we expect a 0th-layer barrel hit for this track?
@ numberOfNextToInnermostPixelLayerHits
these are the hits in the 1st pixel barrel layer
@ expectNextToInnermostPixelLayerHit
Do we expect a 1st-layer barrel hit for this track?
@ numberOfInnermostPixelLayerHits
these are the hits in the 0th pixel barrel layer
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
@ numberOfPixelDeadSensors
number of dead pixel sensors crossed [unit8_t].
A structure for keeping track match information.

◆ initialize()

StatusCode EMTrackMatchBuilder::initialize ( )
finaloverridevirtual

Gaudi algorithm hooks.

Implements IEMTrackMatchBuilder.

Definition at line 39 of file EMTrackMatchBuilder.cxx.

40{
41 ATH_CHECK(m_TrackParticlesKey.initialize());
42 ATH_CHECK(m_caloDetDescrMgrKey.initialize());
43 // the extrapolation tool
45
46 // set things up for the sorting
51
52 return StatusCode::SUCCESS;
53}
function object to sort track matches based on quality
Gaudi::Property< float > m_distanceForScore
The distance from which one goes from using better deltaR to using score.
Gaudi::Property< float > m_deltaPhiRescaleResolution
Gaudi::Property< float > m_deltaPhiResolution
Gaudi::Property< float > m_deltaEtaResolution
The resolutions: might be good to split in barrel/end-cap in the future.
TrackMatchSorter m_sorter

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ interfaceID()

const InterfaceID & IEMTrackMatchBuilder::interfaceID ( )
inlinestaticinherited

AlgTool interface methods.

Definition at line 50 of file IEMTrackMatchBuilder.h.

51{
53}
static const InterfaceID IID_IEMTrackMatchBuilder("IEMTrackMatchBuilder", 1, 0)

◆ isCandidateMatch()

bool EMTrackMatchBuilder::isCandidateMatch ( const xAOD::CaloCluster * cluster,
const xAOD::TrackParticle * track,
bool flip ) const
private

Loose track-cluster matching.

Definition at line 351 of file EMTrackMatchBuilder.cxx.

354{
355 // loose cluster-track matching
356 if (!m_useCandidateMatch) {
357 return true;
358 }
359
360 // Tracking
361 const Trk::Perigee& candidatePerigee = track->perigeeParameters();
362 // Decide whether to try the opposite direction (cosmics)
363 const double trkPhi = (!flip) ? candidatePerigee.parameters()[Trk::phi]
364 : -candidatePerigee.parameters()[Trk::phi];
365 const double trkEta =
366 (!flip) ? candidatePerigee.eta() : -candidatePerigee.eta();
367 const double z_perigee = candidatePerigee.position().z();
368 const double r_perigee = candidatePerigee.position().perp();
369 const Amg::Vector3D PerigeeXYZPosition(candidatePerigee.position().x(),
370 candidatePerigee.position().y(),
371 z_perigee);
372 // Cluster variables
373 const double clusterEta = cluster->eta();
374 const bool isEndCap = !xAOD::EgammaHelpers::isBarrel(cluster);
375 const double Et = cluster->e() / cosh(trkEta);
376 const double clusterPhi = cluster->phi();
377
378 // Avoid clusters with |eta| > 10 or Et less than 10 MeV
379 if (std::abs(clusterEta) > 10.0 || Et < 10) {
380 return false;
381 }
382 // Calculate the eta/phi of the cluster as would be seen from the perigee
383 // position of the Track
384 const Amg::Vector3D XYZClusterWrtTrackPerigee =
386 *cluster, PerigeeXYZPosition, isEndCap);
387
388 const double clusterEtaCorrected = XYZClusterWrtTrackPerigee.eta();
389 // check eta match . Both metrics need to fail in order to disgard the track
390 if ((std::abs(clusterEta - trkEta) > 2. * m_broadDeltaEta) &&
391 (std::abs(clusterEtaCorrected - trkEta) > 2. * m_broadDeltaEta)) {
392 ATH_MSG_DEBUG(" Fails broad window eta match (track eta, cluster eta, "
393 "cluster eta corrected): ( "
394 << trkEta << ", " << clusterEta << ", " << clusterEtaCorrected
395 << ")");
396 return false;
397 }
398 // Calculate the possible rotation of the track
399 // Once assuming the cluster Et being the better estimate (e.g big brem)
400 const double phiRotRescaled = CandidateMatchHelpers::PhiROT(
401 Et, trkEta, track->charge(), r_perigee, isEndCap);
402 // And also assuming the track Pt being correct
403 const double phiRotTrack = CandidateMatchHelpers::PhiROT(
404 track->pt(), trkEta, track->charge(), r_perigee, isEndCap);
405 //
406 const double clusterPhiCorrected = XYZClusterWrtTrackPerigee.phi();
407 // deltaPhi between the track and the cluster
408 const double deltaPhiStd = P4Helpers::deltaPhi(clusterPhiCorrected, trkPhi);
409 // deltaPhi between the track and the cluster accounting for rotation assuming
410 // cluster Et is a better estimator
411 const double trkPhiRescaled = P4Helpers::deltaPhi(trkPhi, phiRotRescaled);
412 const double deltaPhiRescaled =
413 P4Helpers::deltaPhi(clusterPhiCorrected, trkPhiRescaled);
414 // deltaPhi between the track and the cluster accounting for rotation
415 const double trkPhiCorrTrack = P4Helpers::deltaPhi(trkPhi, phiRotTrack);
416 const double deltaPhiTrack =
417 P4Helpers::deltaPhi(clusterPhiCorrected, trkPhiCorrTrack);
418
419 // It has to fail all phi metrics in order to be disgarded
420 if ((std::abs(deltaPhiRescaled) > 2. * m_broadDeltaPhi) &&
421 (std::abs(deltaPhiTrack) > 2. * m_broadDeltaPhi) &&
422 (std::abs(deltaPhiStd) > 2. * m_broadDeltaPhi)) {
423
425 "FAILS broad window phi match (track phi, phirotCluster , phiRotTrack , "
426 << "cluster phi corrected, cluster phi): ( " << trkPhi << ", "
427 << phiRotRescaled << ", " << phiRotTrack << ", " << clusterPhiCorrected
428 << ", " << clusterPhi << ")");
429
430 return false;
431 }
432 // if not false returned we end up here
433 return true;
434}
Gaudi::Property< bool > m_useCandidateMatch
flag to turn on/off use of isCandidateMatch
Gaudi::Property< double > m_broadDeltaEta
broad cut on deltaEta
Gaudi::Property< double > m_broadDeltaPhi
broad cut on deltaPhi
double eta() const
Access method for pseudorapidity - from momentum.
const Amg::Vector3D & position() const
Access method for the position.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double e() const
The total energy of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
Eigen::Matrix< double, 3, 1 > Vector3D
double PhiROT(const double pt, const double eta, const int charge, const double r_start, const bool isEndCap)
Function to calculate the approximate rotation in phi/bending of a track until it reaches the calo.
Amg::Vector3D approxXYZwrtPoint(const xAOD::CaloCluster &cluster, const Amg::Vector3D &point, const bool isEndCap)
Function to get the (x,y,z) of the cluster wrt to a point (x0,y0,z0)
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
Definition P4Helpers.h:34
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ phi
Definition ParamDefs.h:75
@ deltaPhiRescaled
difference between the cluster phi (sampling 2) and the phi of the track extrapolated from the perige...
bool isBarrel(const xAOD::Egamma *eg)
return true if the cluster is in the barrel

◆ msg()

MsgStream & AthCommonMsg< AlgTool >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< AlgTool >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< AlgTool > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysInitialize ( )
overridevirtualinherited

Perform system initialization for an algorithm.

We override this to declare all the elements of handle key arrays at the end of initialization. See comments on updateVHKA.

Reimplemented in asg::AsgMetadataTool, AthCheckedComponent< AthAlgTool >, AthCheckedComponent<::AthAlgTool >, and DerivationFramework::CfAthAlgTool.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ trackExecute()

StatusCode EMTrackMatchBuilder::trackExecute ( const EventContext & ctx,
egammaRec * eg,
const xAOD::TrackParticleContainer * trackPC,
const CaloDetDescrManager & caloDD ) const
private

execute method

Definition at line 85 of file EMTrackMatchBuilder.cxx.

89{
90 if (!eg || !trackPC) {
92 "trackExecute: NULL pointer to egammaRec or TrackParticleContainer");
93 return StatusCode::SUCCESS;
94 }
95 // retrieve corresponding cluster
96 const xAOD::CaloCluster* cluster = eg->caloCluster();
97 // check if the cluster is sane
98 if (cluster && cluster->e() == 0.0) {
99 ATH_MSG_WARNING("trackExecute: cluster energy is 0.0! Ignoring cluster.");
100 return StatusCode::SUCCESS;
101 }
102
103 // Loop over tracks and fill TrackMatch vector
104 std::vector<TrackMatch> trkMatches;
106 for (unsigned int trackNumber = 0; trkIt != trackPC->end();
107 ++trkIt, ++trackNumber) {
108 // Avoid TRT alone
109 if (xAOD::EgammaHelpers::numberOfSiHits(*trkIt) < 4) {
110 continue;
111 }
112 /*
113 * Try with normal directions.
114 * For cosmics allow a retry with inverted direction.
115 */
116 if (isCandidateMatch(cluster, (*trkIt), false)) {
117 inBroadWindow(ctx, trkMatches, *cluster, trackNumber, (**trkIt), caloDD);
118 }
119 }
120
121 if (!trkMatches.empty()) {
122 // sort the track matches
123 std::sort(trkMatches.begin(), trkMatches.end(), m_sorter);
124 // set the matching values
125 TrackMatch bestTrkMatch = trkMatches.at(0);
126 eg->setDeltaEta(bestTrkMatch.deltaEta);
127 eg->setDeltaPhi(bestTrkMatch.deltaPhi);
128 eg->setDeltaPhiRescaled(bestTrkMatch.deltaPhiRescaled);
129 eg->setDeltaPhiLast(bestTrkMatch.deltaPhiLast);
130
131 // set the element Links
132 using EL = ElementLink<xAOD::TrackParticleContainer>;
133 std::vector<EL> trackParticleLinks;
134 trackParticleLinks.reserve(trkMatches.size());
135 const std::string key = EL(*trackPC, 0, ctx).dataID();
136 for (const TrackMatch& m : trkMatches) {
137 ATH_MSG_DEBUG("Match dR: " << m.dR << " second dR: " << m.seconddR
138 << " hasPix: " << m.hasPix
139 << " hitsScore: " << m.hitsScore);
140 if (key.empty()) {
141 trackParticleLinks.emplace_back(*trackPC, m.trackNumber, ctx);
142 } else {
143 trackParticleLinks.emplace_back(key, m.trackNumber, ctx);
144 }
145 }
146 eg->setTrackParticles(trackParticleLinks);
147 }
148 return StatusCode::SUCCESS;
149}
#define ATH_MSG_WARNING(x)
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
bool isCandidateMatch(const xAOD::CaloCluster *cluster, const xAOD::TrackParticle *track, bool flip) const
Loose track-cluster matching.
bool inBroadWindow(const EventContext &ctx, std::vector< TrackMatch > &trackMatches, const xAOD::CaloCluster &cluster, int trackNumber, const xAOD::TrackParticle &trkPB, const CaloDetDescrManager &caloDD) const
Compute for tracks passing the loose matching the distance between track extrapolated to 2nd sampling...
ElementLink< CVec > EL
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
std::size_t numberOfSiHits(const xAOD::TrackParticle *tp)
return the number of Si hits in the track particle
std::vector< ElementLink< xAOD::TrackParticleContainer > > trackParticleLinks(const xAOD::TauJet *tau, xAOD::TauJetParameters::TauTrackFlag flag=xAOD::TauJetParameters::TauTrackFlag::classifiedCharged)
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_broadDeltaEta

Gaudi::Property<double> EMTrackMatchBuilder::m_broadDeltaEta
private
Initial value:
{ this,
"broadDeltaEta",
0.1,
"Value of broad cut for delta eta" }

broad cut on deltaEta

Definition at line 125 of file EMTrackMatchBuilder.h.

125 { this,
126 "broadDeltaEta",
127 0.1,
128 "Value of broad cut for delta eta" };

◆ m_broadDeltaPhi

Gaudi::Property<double> EMTrackMatchBuilder::m_broadDeltaPhi
private
Initial value:
{ this,
"broadDeltaPhi",
0.1,
"Value of broad cut for delta phi" }

broad cut on deltaPhi

Definition at line 131 of file EMTrackMatchBuilder.h.

131 { this,
132 "broadDeltaPhi",
133 0.1,
134 "Value of broad cut for delta phi" };

◆ m_caloDetDescrMgrKey

SG::ReadCondHandleKey<CaloDetDescrManager> EMTrackMatchBuilder::m_caloDetDescrMgrKey
private
Initial value:
{
this,
"CaloDetDescrManager",
"CaloDetDescrManager",
"SG Key for CaloDetDescrManager in the Condition Store"
}

Definition at line 117 of file EMTrackMatchBuilder.h.

117 {
118 this,
119 "CaloDetDescrManager",
120 "CaloDetDescrManager",
121 "SG Key for CaloDetDescrManager in the Condition Store"
122 };

◆ m_deltaEtaResolution

Gaudi::Property<float> EMTrackMatchBuilder::m_deltaEtaResolution
private
Initial value:
{ this,
"DeltaEtaResolution",
1.0,
"The deltaEta resolution" }

The resolutions: might be good to split in barrel/end-cap in the future.

Definition at line 213 of file EMTrackMatchBuilder.h.

213 { this,
214 "DeltaEtaResolution",
215 1.0,
216 "The deltaEta resolution" };

◆ m_deltaEtaWeight

double EMTrackMatchBuilder::m_deltaEtaWeight {}
private

Definition at line 256 of file EMTrackMatchBuilder.h.

256{};

◆ m_deltaPhiRescaleResolution

Gaudi::Property<float> EMTrackMatchBuilder::m_deltaPhiRescaleResolution
private
Initial value:
{
this,
"DeltaPhiRescaleResolution",
1.0,
"The deltaPhiRescale resolution"
}

Definition at line 223 of file EMTrackMatchBuilder.h.

223 {
224 this,
225 "DeltaPhiRescaleResolution",
226 1.0,
227 "The deltaPhiRescale resolution"
228 };

◆ m_deltaPhiRescaleWeight

double EMTrackMatchBuilder::m_deltaPhiRescaleWeight {}
private

Definition at line 258 of file EMTrackMatchBuilder.h.

258{};

◆ m_deltaPhiResolution

Gaudi::Property<float> EMTrackMatchBuilder::m_deltaPhiResolution
private
Initial value:
{ this,
"DeltaPhiResolution",
1.0,
"The deltaPhi resolution" }

Definition at line 218 of file EMTrackMatchBuilder.h.

218 { this,
219 "DeltaPhiResolution",
220 1.0,
221 "The deltaPhi resolution" };

◆ m_deltaPhiWeight

double EMTrackMatchBuilder::m_deltaPhiWeight {}
private

Definition at line 257 of file EMTrackMatchBuilder.h.

257{};

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_distanceForScore

Gaudi::Property<float> EMTrackMatchBuilder::m_distanceForScore
private
Initial value:
{
this,
"DistanceForScore",
0.01,
"The distance from which one goes from using better deltaR to using score."
}

The distance from which one goes from using better deltaR to using score.

Note that this distance varies depending on the resolutions entered above. If you don't use resolutions (resolution = 1.0) this becomes deltaR distance.

Definition at line 235 of file EMTrackMatchBuilder.h.

235 {
236 this,
237 "DistanceForScore",
238 0.01,
239 "The distance from which one goes from using better deltaR to using score."
240 };

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extrapolationTool

ToolHandle<IEMExtrapolationTools> EMTrackMatchBuilder::m_extrapolationTool
private
Initial value:
{
this,
"ExtrapolationTool",
"EMExtrapolationTools",
"Name of the extrapolation tool"
}

Definition at line 242 of file EMTrackMatchBuilder.h.

242 {
243 this,
244 "ExtrapolationTool",
245 "EMExtrapolationTools",
246 "Name of the extrapolation tool"
247 };

◆ m_isCosmics

Gaudi::Property<bool> EMTrackMatchBuilder::m_isCosmics
private
Initial value:
{ this,
"isCosmics",
false,
"Boolean for use of cosmics" }

Definition at line 250 of file EMTrackMatchBuilder.h.

250 { this,
251 "isCosmics",
252 false,
253 "Boolean for use of cosmics" };

◆ m_MaxDeltaPhiRescale

Gaudi::Property<double> EMTrackMatchBuilder::m_MaxDeltaPhiRescale
private
Initial value:
{
this,
"MaxDeltaPhiRescale",
0.25,
"Maximum Value of the deltaPhi rescale"
}

@Maximum deltaPhi (Res) allowed for a match

Definition at line 177 of file EMTrackMatchBuilder.h.

177 {
178 this,
179 "MaxDeltaPhiRescale",
180 0.25,
181 "Maximum Value of the deltaPhi rescale"
182 };

◆ m_narrowDeltaEta

Gaudi::Property<double> EMTrackMatchBuilder::m_narrowDeltaEta
private
Initial value:
{
this,
"narrowDeltaEta",
0.05,
"Value of narrow cut for delta eta"
}

narrow cut on deltaEta

Definition at line 137 of file EMTrackMatchBuilder.h.

137 {
138 this,
139 "narrowDeltaEta",
140 0.05,
141 "Value of narrow cut for delta eta"
142 };

◆ m_narrowDeltaPhi

Gaudi::Property<double> EMTrackMatchBuilder::m_narrowDeltaPhi
private
Initial value:
{
this,
"narrowDeltaPhi",
0.05,
"Value of the narrowd cut for delta phi"
}

narrow cut on deltaPhiRescale

Definition at line 145 of file EMTrackMatchBuilder.h.

145 {
146 this,
147 "narrowDeltaPhi",
148 0.05,
149 "Value of the narrowd cut for delta phi"
150 };

◆ m_narrowDeltaPhiBrem

Gaudi::Property<double> EMTrackMatchBuilder::m_narrowDeltaPhiBrem
private
Initial value:
{
this,
"narrowDeltaPhiBrem",
0.1,
"Value of the narrow cut for delta phi Brem"
}

narrow cut on deltaPhi for electrons

Definition at line 153 of file EMTrackMatchBuilder.h.

153 {
154 this,
155 "narrowDeltaPhiBrem",
156 0.1,
157 "Value of the narrow cut for delta phi Brem"
158 };

◆ m_narrowDeltaPhiRescale

Gaudi::Property<double> EMTrackMatchBuilder::m_narrowDeltaPhiRescale
private
Initial value:
{
this,
"narrowDeltaPhiRescale",
0.05,
"Value of the narrow cut for delta phi Rescale"
}

narrow cut on deltaPhiRescale

Definition at line 161 of file EMTrackMatchBuilder.h.

161 {
162 this,
163 "narrowDeltaPhiRescale",
164 0.05,
165 "Value of the narrow cut for delta phi Rescale"
166 };

◆ m_narrowDeltaPhiRescaleBrem

Gaudi::Property<double> EMTrackMatchBuilder::m_narrowDeltaPhiRescaleBrem
private
Initial value:
{
this,
"narrowDeltaPhiRescaleBrem",
0.1,
"Value of the narrow cut for delta phi Rescale Brem"
}

narrow cut on deltaPhiRescale for electrons

Definition at line 169 of file EMTrackMatchBuilder.h.

169 {
170 this,
171 "narrowDeltaPhiRescaleBrem",
172 0.1,
173 "Value of the narrow cut for delta phi Rescale Brem"
174 };

◆ m_SecondPassRescale

Gaudi::Property<bool> EMTrackMatchBuilder::m_SecondPassRescale
private
Initial value:
{ this,
"SecondPassRescale",
true,
"Do second pass with rescale" }

Boolean to do second pass with Rescale.

Definition at line 207 of file EMTrackMatchBuilder.h.

207 { this,
208 "SecondPassRescale",
209 true,
210 "Do second pass with rescale" };

◆ m_sorter

TrackMatchSorter EMTrackMatchBuilder::m_sorter
private

Definition at line 260 of file EMTrackMatchBuilder.h.

◆ m_TrackParticlesKey

SG::ReadHandleKey<xAOD::TrackParticleContainer> EMTrackMatchBuilder::m_TrackParticlesKey
private
Initial value:
{
this,
"TrackParticlesName",
"",
"Name of the input track particle container"
}

name of TrackParticle container in TDS

Definition at line 110 of file EMTrackMatchBuilder.h.

110 {
111 this,
112 "TrackParticlesName",
113 "",
114 "Name of the input track particle container"
115 };

◆ m_useCandidateMatch

Gaudi::Property<bool> EMTrackMatchBuilder::m_useCandidateMatch
private
Initial value:
{
this,
"useCandidateMatch",
true,
"Boolean to use candidate matching"
}

flag to turn on/off use of isCandidateMatch

Definition at line 185 of file EMTrackMatchBuilder.h.

185 {
186 this,
187 "useCandidateMatch",
188 true,
189 "Boolean to use candidate matching"
190 };

◆ m_useRescaleMetric

Gaudi::Property<bool> EMTrackMatchBuilder::m_useRescaleMetric
private
Initial value:
{ this,
"UseRescaleMetric",
true,
"Use Rescale Metric" }

Boolean to use Rescale in the metric.

Definition at line 201 of file EMTrackMatchBuilder.h.

201 { this,
202 "UseRescaleMetric",
203 true,
204 "Use Rescale Metric" };

◆ m_useScoring

Gaudi::Property<bool> EMTrackMatchBuilder::m_useScoring
private
Initial value:
{
this,
"useScoring",
true,
"Boolean to apply heuristic when tracks have close deltaR"
}

Boolean to apply heuristic when tracks have close deltaR.

Definition at line 193 of file EMTrackMatchBuilder.h.

193 {
194 this,
195 "useScoring",
196 true,
197 "Boolean to apply heuristic when tracks have close deltaR"
198 };

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< AlgTool > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< AlgTool > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


The documentation for this class was generated from the following files: