ATLAS Offline Software
Loading...
Searching...
No Matches
EleEleOverlapTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3*/
4
5// System includes
6#include <typeinfo>
7#include <exception>
8
9// Framework includes
12
13// EDM includes
15
16// Local includes
18
19namespace
20{
21
24 bool isDummyVal(float x)
25 {
26 const float epsilon = 1e-5;
27 const float dummyVal = -999;
28 return std::abs(x - dummyVal) < epsilon;
29 }
30
31 // Define simple exception for Dummy/missing values
32 class DummyValError : public std::exception
33 {};
34
35}
36
37namespace ORUtils
38{
39
40 //---------------------------------------------------------------------------
41 // Constructor
42 //---------------------------------------------------------------------------
43 EleEleOverlapTool::EleEleOverlapTool(const std::string& name)
44 : BaseOverlapTool(name),
45 m_useTrackMatch(true),
46 m_useClusterMatch(false),
47 m_clusterDeltaEta(3*0.025),
48 m_clusterDeltaPhi(5*0.025)
49 {
50 declareProperty("UseTrackMatch", m_useTrackMatch,
51 "Match electrons by shared track");
52 declareProperty("UseClusterMatch", m_useClusterMatch,
53 "Match electrons by cluster proximity");
54 declareProperty("ClusterDeltaEta", m_clusterDeltaEta,
55 "Cluster matching delta eta");
56 declareProperty("ClusterDeltaPhi", m_clusterDeltaPhi,
57 "Cluster matching delta phi");
58 }
59
60 //---------------------------------------------------------------------------
61 // Initialize
62 //---------------------------------------------------------------------------
64 {
65 ATH_MSG_DEBUG("UseTrackMatch " << m_useTrackMatch <<
66 " UseClusterMatch " << m_useClusterMatch <<
67 " ClusterDeltaEta " << m_clusterDeltaEta <<
68 " ClusterDeltaPhi " << m_clusterDeltaPhi);
69
70 // Sanity check
72 ATH_MSG_ERROR("You must enable at least one: UseTrackMatch or UseClusterMatch");
73 return StatusCode::FAILURE;
74 }
76 {
77 // Initialize the cluster accessors
78 resetAccessor (m_accessors->m_clusterContainerAcc, *m_accessors, "egammaClusters");
79 resetAccessor (m_accessors->m_caloClusterAcc, *m_accessors, "caloClusterLinks");
80 m_accessors->m_etaBEAcc.emplace (*m_accessors);
81 m_accessors->m_phiBEAcc.emplace (*m_accessors);
82 }
84 {
85 // Initialize the track accessor
86 resetAccessor (m_accessors->m_track0Acc, *m_accessors, "GSFTrackParticles");
87 resetAccessor (m_accessors->m_trackAcc, *m_accessors, "trackParticleLinks");
88 }
89
90 return StatusCode::SUCCESS;
91 }
92
93 //---------------------------------------------------------------------------
94 // Identify overlaps
95 //---------------------------------------------------------------------------
99 columnar::EventContextId /*eventContext*/) const
100 {
102 {
103 // I require that the two containers are the same so that I can
104 // arbitrarily pick one of them to use.
105 if(&cont1.getXAODObject() != &cont2.getXAODObject()) {
106 ATH_MSG_ERROR("This tool expects both electron containers to be the " <<
107 "same for now");
108 return StatusCode::FAILURE;
109 }
110 }
111 // Check the container type
112 ATH_CHECK( checkForXAODContainer<xAOD::ElectronContainer>(cont1, "Container is not of type ElectronContainer!") );
113 // Call the type-specific method
115 return StatusCode::SUCCESS;
116 }
117
118 //---------------------------------------------------------------------------
119 // Identify overlaps
120 //---------------------------------------------------------------------------
123 {
124 ATH_MSG_DEBUG("Removing overlapping electrons");
125
126 // Initialize output decorations if necessary
127 initializeDecorations(electrons);
128
129 // TODO: consider adding cluster-based matching also
130
131 // Loop over surviving electron pairs
132 for(const auto el1 : electrons) {
133 if(!isSurvivingObject(el1)) continue;
134 for(const auto el2 : electrons) {
135 if(el1 == el2) continue;
136 if(!isSurvivingObject(el2)) continue;
137
138 // Perform the match and decide whether to reject el1
139 try {
140 if(electronsMatch(el1, el2) && rejectFirst(el1, el2)) {
141 ATH_CHECK( handleOverlap(el1, el2) );
142 }
143 }
144 catch(const DummyValError& e) {
145 ATH_MSG_ERROR("Cluster 2nd sampling eta/phi values are -999. " <<
146 "It seems you are missing the neccesary variables to "
147 "do the requested electron-electron cluster matching");
148 return StatusCode::FAILURE;
149 }
150 }
151 }
152
153 return StatusCode::SUCCESS;
154 }
155
156 //---------------------------------------------------------------------------
157 // Apply the ele-ele matching criteria
158 //---------------------------------------------------------------------------
161 {
162 auto& acc = *m_accessors;
163 // Look for a shared track
165 {
166 auto trk1 = el1(acc.m_trackAcc);
167 auto trk2 = el2(acc.m_trackAcc);
168 if(trk1.size() > 0 && trk2.size() > 0 && trk1[0] == trk2[0])
169 return true;
170 }
171
172 // Look for overlapping clusters
174 auto clus1 = el1(acc.m_caloClusterAcc)[0].value();
175 auto clus2 = el2(acc.m_caloClusterAcc)[0].value();
177
178 // We use coordinates from 2nd sampling
179 const unsigned layer = 2;
180 const float eta1 = clus1(acc.m_etaBEAcc.value(), layer);
181 const float eta2 = clus2(acc.m_etaBEAcc.value(), layer);
182 const float phi1 = clus1(acc.m_phiBEAcc.value(), layer);
183 const float phi2 = clus2(acc.m_phiBEAcc.value(), layer);
184
185 // Check validity of the eta/phi (no dummy -999 values)
186 if(isDummyVal(eta1) || isDummyVal(eta2) ||
187 isDummyVal(phi1) || isDummyVal(phi2)) {
188 throw DummyValError();
189 }
190
191 const float dEta = eta1 - eta2;
192 const float dPhi = deltaPhi(phi1, phi2);
193
194 if( std::abs(dEta) < m_clusterDeltaEta &&
195 std::abs(dPhi) < m_clusterDeltaPhi )
196 {
197 return true;
198 }
199 }
200
201 return false;
202 }
203
204 //---------------------------------------------------------------------------
205 // Decide whether to reject the first electron compared to the 2nd.
206 // This function assumes a matching criteria has already been applied.
207 //---------------------------------------------------------------------------
210 {
211 auto& acc = *m_accessors;
212 // TODO: consider incorporating track-match information in the priority
213 // selection.
214
215 // Always reject author "Ambiguous" when compared to author "Electron"
216 if(el1(acc.m_authorAcc) == xAOD::EgammaParameters::AuthorAmbiguous &&
217 el2(acc.m_authorAcc) == xAOD::EgammaParameters::AuthorElectron)
218 return true;
219 if(el1(acc.m_authorAcc) == xAOD::EgammaParameters::AuthorElectron &&
220 el2(acc.m_authorAcc) == xAOD::EgammaParameters::AuthorAmbiguous)
221 return false;
222
223 // Reject the softer electron
224 if(el1(acc.m_ptAcc) < el2(acc.m_ptAcc)) return true;
225 return false;
226 }
227
228} // namespace ORUtils
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
DataVector adapter that acts like it holds const pointers.
#define x
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
void initializeDecorations(columnar::Particle1Range container) const
BaseOverlapTool(const std::string &name)
Create proper constructor for Athena.
StatusCode handleOverlap(const columnar::ObjectId< CI1, CM > &testParticle, const columnar::ObjectId< CI2, CM > &refParticle) const
Common helper method to handle an overlap result.
bool isSurvivingObject(columnar::Particle1Id obj) const
StatusCode checkForXAODContainer(columnar::ObjectRange< CI, CM > cont, std::string_view message) const
check whether the container is of the right type
double m_clusterDeltaPhi
Cluster-matching dPhi.
virtual StatusCode findOverlaps(columnar::Particle1Range cont1, columnar::Particle2Range cont2, columnar::EventContextId eventContext) const override
Identify overlapping electrons.
EleEleOverlapTool(const std::string &name)
Create proper constructor for Athena.
virtual StatusCode initializeDerived() override
Initialize the tool.
bool rejectFirst(columnar::Particle1Id el1, columnar::Particle1Id el2) const
Helper method to decide which electron to reject.
bool electronsMatch(columnar::Particle1Id el1, columnar::Particle1Id el2) const
Helper method for matching electrons.
double m_clusterDeltaEta
Cluster-matching dEta.
std::unique_ptr< Accessors > m_accessors
virtual StatusCode internalFindOverlaps(columnar::Particle1Range electrons) const
Identify overlapping electrons and jets.
bool m_useTrackMatch
Match electrons by shared track (on by default)
bool m_useClusterMatch
Match electrons by cluster distance (off by default)
ObjectId< ContainerId::particle1 > Particle1Id
Definition ParticleDef.h:48
ObjectId< ContainerId::eventContext > EventContextId
ObjectRange< ContainerId::particle2 > Particle2Range
Definition ParticleDef.h:53
ObjectRange< ContainerId::particle1 > Particle1Range
Definition ParticleDef.h:47
const uint16_t AuthorElectron
Object Reconstructed by standard cluster-based algorithm.
Definition EgammaDefs.h:24
const uint16_t AuthorAmbiguous
Object Reconstructed by standard cluster-based algorithm.
Definition EgammaDefs.h:32
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
static constexpr bool isXAOD
Whether this is the xAOD mode.
Definition ColumnarDef.h:20