ATLAS Offline Software
Loading...
Searching...
No Matches
CorrectPFOTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5// CorrectPFOTool.cxx
6
10
11#include <cmath>
12
13CorrectPFOTool::CorrectPFOTool(const std::string &name):
15 m_weightPFOTool("",this) {
16
17 // Configuration
18 declareProperty("WeightPFOTool", m_weightPFOTool, "Name of tool that extracts the cPFO weights.");
19 declareProperty("InputIsEM", m_inputIsEM = true, "True if neutral PFOs are EM scale clusters.");
20 declareProperty("CalibratePFO", m_calibrate = false, "True if LC calibration should be applied to EM PFOs.");
21 declareProperty("CorrectNeutral", m_correctneutral = true, "True to use the neutral component of PFlow.");
22 declareProperty("CorrectCharged", m_correctcharged = true, "True if use the charged component of PFlow.");
23 declareProperty("UseChargedWeights",m_useChargedWeights = true, "True if we make use of weighting scheme for charged PFO");
24 declareProperty("DoByVertex", m_doByVertex = false, "True to add vertex-by-vertex corrections for neutral PFOs");
25
26 // Input properties
27 declareProperty("VertexContainerKey",
28 m_vertexContainer_key="PrimaryVertices",
29 "Datahandle key for the primary vertex container");
30}
31
34 ATH_MSG_ERROR("This tool is configured to do nothing!");
35 return StatusCode::FAILURE;
36 }
38 ATH_MSG_ERROR("CorrectPFOTool requires PFO inputs. It cannot operate on objects of type "
39 << m_inputType);
40 return StatusCode::FAILURE;
41 }
43 ATH_CHECK( m_weightPFOTool.retrieve() );
44 }
45 ATH_CHECK( m_vertexContainer_key.initialize() );
46
47 if(m_doByVertex){
48 ATH_MSG_INFO("Running CorrectPFOTool by vertex:" << m_doByVertex);
49 }
50
51 return StatusCode::SUCCESS;
52}
53
55 // Type-checking happens in the JetConstituentModifierBase class
56 // so it is safe just to static_cast
58 xAOD::FlowElementContainer* feCont = static_cast<xAOD::FlowElementContainer*>(cont);
59 if(!feCont->empty() && !(feCont->front()->signalType() & xAOD::FlowElement::PFlow)){
60 ATH_MSG_ERROR("CorrectPFOTool received FlowElements that aren't PFOs");
61 return StatusCode::FAILURE;
62 }
63 return m_doByVertex ? correctPFOByVertex(*feCont) : correctPFO(*feCont);
64 }
65 xAOD::PFOContainer* pfoCont = static_cast<xAOD::PFOContainer*> (cont);
66 return m_doByVertex ? correctPFOByVertex(*pfoCont) : correctPFO(*pfoCont);
67}
68
70 // Retrieve Primary Vertices
72 if (!handle.isValid()){
73 ATH_MSG_WARNING(" This event has no primary vertex container" );
74 return nullptr;
75 }
76
77 const xAOD::VertexContainer* pvtxs = handle.cptr();
78 if(pvtxs->empty()){
79 ATH_MSG_WARNING(" Failed to retrieve valid primary vertex container" );
80 return nullptr;
81 }
82
83 //Usually the 0th vertex is the primary one, but this is not always
84 // the case. So we will choose the first vertex of type PriVtx
85 for (const auto *theVertex : *pvtxs) {
86 if (theVertex->vertexType()==xAOD::VxType::PriVtx) {
87 return theVertex;
88 }//If we have a vertex of type primary vertex
89 }//iterate over the vertices and check their type
90
91 // If we failed to find an appropriate vertex, return the dummy vertex
92 ATH_MSG_DEBUG("Could not find a primary vertex in this event " );
93 for (const auto *theVertex : *pvtxs) {
94 if (theVertex->vertexType()==xAOD::VxType::NoVtx) {
95 return theVertex;
96 }
97 }
98
99 // If there is no primary vertex, then we cannot do PV matching.
100 ATH_MSG_WARNING("Primary vertex container was empty");
101 return nullptr;
102}
103
105
106 const xAOD::Vertex* vtx = nullptr;
107 if(m_correctneutral) {
108 vtx = getPrimaryVertex();
109 if(vtx==nullptr) {
110 ATH_MSG_ERROR("Primary vertex container was empty or no valid vertex found!");
111 return StatusCode::FAILURE;
112 } else if (vtx->vertexType()==xAOD::VxType::NoVtx) {
113 ATH_MSG_VERBOSE("No genuine primary vertex found. Will not apply origin correction");
114 }
115 }
116
117 for ( xAOD::PFO* ppfo : cont ) {
118
119 if ( std::abs(ppfo->charge())<FLT_MIN) { // Neutral PFOs
120 if(m_correctneutral) {
121 ATH_CHECK( applyNeutralCorrection(*ppfo, *vtx) );
122 }
123 } else { // Charged PFOs
124 if(m_correctcharged) {
126 }
127 }
128 } // PFO loop
129
130 return StatusCode::SUCCESS;
131}
132
134
135 const xAOD::Vertex* vtx = nullptr;
136 if(m_correctneutral) {
137 vtx = getPrimaryVertex();
138 if(vtx==nullptr) {
139 ATH_MSG_ERROR("Primary vertex container was empty or no valid vertex found!");
140 return StatusCode::FAILURE;
141 } else if (vtx->vertexType()==xAOD::VxType::NoVtx) {
142 ATH_MSG_VERBOSE("No genuine primary vertex found. Will not apply origin correction");
143 }
144 }
145
146 for ( xAOD::FlowElement* ppfo : cont ) {
147
148 if ( !ppfo->isCharged()) { // Neutral PFOs
149 if(m_correctneutral) {
150 ATH_CHECK( applyNeutralCorrection(*ppfo, *vtx) );
151 }
152 } else { // Charged PFOs
153 if(m_correctcharged) {
155 }
156 }
157 } // PFO loop
158
159 return StatusCode::SUCCESS;
160}
161
163 static const SG::AuxElement::Accessor<unsigned> copyIndex("ConstituentCopyIndex");
164 // Retrieve Primary Vertices
166 if (!handle.isValid()){
167 ATH_MSG_WARNING(" This event has no primary vertex container" );
168 return StatusCode::FAILURE;
169 }
170
171 const xAOD::VertexContainer* pvtxs = handle.cptr();
172 if(pvtxs->empty()){
173 ATH_MSG_WARNING(" Failed to retrieve valid primary vertex container" );
174 return StatusCode::FAILURE;
175 }
176
177 for ( xAOD::PFO* ppfo : cont ) {
178
179 if ( std::abs(ppfo->charge())<FLT_MIN) { // Neutral PFOs
180 if(m_correctneutral) {
181 // Neutral PFOs - there are copies, one per vertex, already created
182 // We need to now need to correct each copy to point to the corresponding vertex
183 if (!copyIndex.isAvailable(*ppfo))
184 {
185 ATH_MSG_WARNING("Encountered a neutral per-vertex PFO object without the corresponding vertex index attribute");
186 continue;
187 }
188
189 const unsigned iVtx = copyIndex(*ppfo);
190 if (iVtx >= pvtxs->size())
191 {
192 ATH_MSG_WARNING("Encountered a neutral per-vertex PFO object with an index beyond the size of the vertex container");
193 continue;
194 }
195 const xAOD::Vertex* vtx = pvtxs->at(iVtx);
196 // Determine the correct four-vector interpretation for a neutral PFO associated with this vertex
197 // Cannot use applyNeutralCorrection as that only used VxType::PriVtx, but we want to apply to all primary vertices
198 // As such, extract the relevant parts and modify appropriately here
199 // This is necessary to avoid changing sign of pT for pT<0 PFO
200 if (ppfo->e() < FLT_MIN) {
201 ppfo->setP4(0, 0, 0, 0);
202 }
203 if (!m_inputIsEM || m_calibrate) { // Use LC four-vector
204 ppfo->setP4(ppfo->GetVertexCorrectedFourVec(*vtx));
205 } else { // Use EM four-vector
206 ppfo->setP4(ppfo->GetVertexCorrectedEMFourVec(*vtx));
207 }
208 }
209 } else { // Charged PFOs
210 if(m_correctcharged) {
212 }
213 }
214 } // PFO loop
215
216
217 return StatusCode::SUCCESS;
218}
219
221 static const SG::AuxElement::Accessor<unsigned> copyIndex("ConstituentCopyIndex");
222
223 // Retrieve Primary Vertices
225 if (!handle.isValid()){
226 ATH_MSG_WARNING(" This event has no primary vertex container" );
227 return StatusCode::FAILURE;
228 }
229
230 const xAOD::VertexContainer* pvtxs = handle.cptr();
231 if(pvtxs->empty()){
232 ATH_MSG_WARNING(" Failed to retrieve valid primary vertex container" );
233 return StatusCode::FAILURE;
234 }
235
236 for ( xAOD::FlowElement* ppfo : cont ) {
237
238 if ( !ppfo->isCharged()) { // Neutral PFOs
239 if(m_correctneutral) {
240 // Neutral PFOs - there are copies, one per vertex, already created
241 // We need to now need to correct each copy to point to the corresponding vertex
242 if (!copyIndex.isAvailable(*ppfo))
243 {
244 ATH_MSG_WARNING("Encountered a neutral per-vertex PFO object without the corresponding vertex index attribute");
245 continue;
246 }
247
248 const unsigned iVtx = copyIndex(*ppfo);
249 if (iVtx >= pvtxs->size())
250 {
251 ATH_MSG_WARNING("Encountered a neutral per-vertex PFO object with an index beyond the size of the vertex container");
252 continue;
253 }
254 const xAOD::Vertex* vtx = pvtxs->at(iVtx);
255 // Determine the correct four-vector interpretation for a neutral PFO associated with this vertex
256 // Cannot use applyNeutralCorrection as that only used VxType::PriVtx, but we want to apply to all primary vertices
257 // As such, extract the relevant parts and modify appropriately here
258 // This is necessary to avoid changing sign of pT for pT<0 PFO
259 if (ppfo->e() < FLT_MIN) {
260 ppfo->setP4(0, 0, 0, 0);
261 }
262 ppfo->setP4(FEHelpers::getVertexCorrectedFourVec(*ppfo, *vtx));
263 }
264 } else { // Charged PFOs
265 if(m_correctcharged) {
267 }
268 }
269 } // PFO loop
270
271
272 return StatusCode::SUCCESS;
273}
274
275
277 if (pfo.e() < FLT_MIN) { //This is necessary to avoid changing sign of pT for pT<0 PFO
278 pfo.setP4(0,0,0,0);
279 } else {
280 if ( !m_inputIsEM || m_calibrate ) { // Use LC four-vector
281 // Only correct if we really reconstructed a vertex
283 pfo.setP4(pfo.GetVertexCorrectedFourVec(vtx));
284 }
285 } else { // Use EM four-vector
286 // Only apply origin correction if we really reconstructed a vertex
289 } else {pfo.setP4(pfo.p4EM());} // Just set EM 4-vec
290 }
291 }
292 return StatusCode::SUCCESS;
293}
294
296 if (pfo.e() < FLT_MIN) { //This is necessary to avoid changing sign of pT for pT<0 PFO
297 pfo.setP4(0,0,0,0);
298 }
299 // Only apply origin correction if we really reconstructed a vertex
300 else if(vtx.vertexType() == xAOD::VxType::PriVtx) {
302 }
303 return StatusCode::SUCCESS;
304}
305
308 float weight = 0.0;
309 ATH_CHECK( m_weightPFOTool->fillWeight( pfo, weight ) );
310 ATH_MSG_VERBOSE("Fill pseudojet for CPFO with weighted pt: " << pfo.pt()*weight);
311 pfo.setP4(pfo.p4()*weight);
312 }//if should use charged PFO weighting scheme
313 return StatusCode::SUCCESS;
314}
315
318 float weight = 0.0;
319 ATH_CHECK( m_weightPFOTool->fillWeight( pfo, weight ) );
320 ATH_MSG_VERBOSE("Fill pseudojet for CPFO with weighted pt: " << pfo.pt()*weight);
321 pfo.setP4(pfo.p4()*weight);
322 }//if should use charged PFO weighting scheme
323 return StatusCode::SUCCESS;
324}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Handle class for reading from StoreGate.
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
StatusCode correctPFOByVertex(xAOD::PFOContainer &cont) const
bool m_calibrate
If true EM clusters are used for neutral PFOs.
StatusCode correctPFO(xAOD::PFOContainer &cont) const
StatusCode process_impl(xAOD::IParticleContainer *cont) const override final
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexContainer_key
StatusCode applyChargedCorrection(xAOD::PFO &pfo) const
CorrectPFOTool(const std::string &name)
StatusCode applyNeutralCorrection(xAOD::PFO &pfo, const xAOD::Vertex &vtx) const
ToolHandle< CP::IWeightPFOTool > m_weightPFOTool
const xAOD::Vertex * getPrimaryVertex() const
bool m_correctneutral
If true, EM PFOs are calibrated to LC.
StatusCode initialize() override final
Dummy implementation of the initialisation function.
const T * at(size_type n) const
Access an element, as an rvalue.
const T * front() const
Access the first element in the collection as an rvalue.
size_type size() const noexcept
Returns the number of elements in the collection.
bool empty() const noexcept
Returns true if the collection is empty.
JetConstituentModifierBase(const std::string &name)
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:573
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
virtual double pt() const override
void setP4(float pt, float eta, float phi, float m)
signal_t signalType() const
virtual double e() const override
The total energy of the particle.
virtual FourMom_t p4() const override
The full 4-momentum of the particle.
TLorentzVector GetVertexCorrectedFourVec(const xAOD::Vertex &vertexToCorrectTo) const
Correct 4-vector to point at a vertex.
Definition PFO_v1.cxx:722
FourMom_t p4EM() const
get EM scale 4-vector
Definition PFO_v1.cxx:144
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition PFO_v1.cxx:95
void setP4(const FourMom_t &vec)
set the 4-vec
Definition PFO_v1.cxx:107
TLorentzVector GetVertexCorrectedEMFourVec(const xAOD::Vertex &vertexToCorrectTo) const
Correct EM scale 4-vector to point at a vertex.
Definition PFO_v1.cxx:737
virtual double e() const
The total energy of the particle.
Definition PFO_v1.cxx:81
virtual double pt() const
The transverse momentum ( ) of the particle.
Definition PFO_v1.cxx:52
VxType::VertexType vertexType() const
The type of the vertex.
TLorentzVector getVertexCorrectedFourVec(const xAOD::FlowElement &fe, const xAOD::Vertex &vertexToCorrectTo)
Definition FEHelpers.cxx:13
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
@ ParticleFlow
The object is a particle-flow object.
Definition ObjectType.h:41
@ FlowElement
The object is a track-calo-cluster.
Definition ObjectType.h:52
@ PriVtx
Primary vertex.
@ NoVtx
Dummy vertex. TrackParticle was not used in vertex fit.
FlowElementContainer_v1 FlowElementContainer
Definition of the current "pfo container version".
PFO_v1 PFO
Definition of the current "pfo version".
Definition PFO.h:17
PFOContainer_v1 PFOContainer
Definition of the current "pfo container version".
FlowElement_v1 FlowElement
Definition of the current "pfo version".
Definition FlowElement.h:16
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
DataVector< IParticle > IParticleContainer
Simple convenience declaration of IParticleContainer.