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

Applies ML corrections to PFO ML corrections are stored in linked CaloClusters as decorations. More...

#include <NeutralPFOClusterMLCorrectionTool.h>

Inheritance diagram for NeutralPFOClusterMLCorrectionTool:
Collaboration diagram for NeutralPFOClusterMLCorrectionTool:

Public Member Functions

 NeutralPFOClusterMLCorrectionTool (const std::string &type, const std::string &name, const IInterface *parent)
virtual ~NeutralPFOClusterMLCorrectionTool ()=default
virtual StatusCode initialize () override
virtual void correctContainer (xAOD::FlowElementContainer &pfos) const override

Private Member Functions

void scaleEnergyToAlternativeSignalState (xAOD::FlowElement &pfo, const xAOD::CaloCluster &cls) const
const xAOD::CaloClustergetLinkedCluster (const xAOD::FlowElement &pfo) const

Private Attributes

Gaudi::Property< std::string > m_clusterMLCorrectedEnergyKey

Detailed Description

Applies ML corrections to PFO ML corrections are stored in linked CaloClusters as decorations.

The correction is applied as a scale factor to the PFO.

Definition at line 26 of file NeutralPFOClusterMLCorrectionTool.h.

Constructor & Destructor Documentation

◆ NeutralPFOClusterMLCorrectionTool()

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

Definition at line 7 of file NeutralPFOClusterMLCorrectionTool.cxx.

8 : base_class(type, name, parent) {}

◆ ~NeutralPFOClusterMLCorrectionTool()

virtual NeutralPFOClusterMLCorrectionTool::~NeutralPFOClusterMLCorrectionTool ( )
virtualdefault

Member Function Documentation

◆ correctContainer()

void NeutralPFOClusterMLCorrectionTool::correctContainer ( xAOD::FlowElementContainer & pfos) const
overridevirtual

Definition at line 15 of file NeutralPFOClusterMLCorrectionTool.cxx.

16{
17
18 for (xAOD::FlowElement *pfo : container)
19 {
20 if (pfo->isCharged())
21 {
22 ATH_MSG_WARNING("NeutralPFOClusterMLCorrectionTool: Charged FlowElement found in neutral FlowElementContainer with index " + std::to_string(pfo->index()));
23 continue;
24 }
25
27 if (cls != nullptr)
29 }
30}
#define ATH_MSG_WARNING(x)
const xAOD::CaloCluster * getLinkedCluster(const xAOD::FlowElement &pfo) const
void scaleEnergyToAlternativeSignalState(xAOD::FlowElement &pfo, const xAOD::CaloCluster &cls) const
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
FlowElement_v1 FlowElement
Definition of the current "pfo version".
Definition FlowElement.h:16

◆ getLinkedCluster()

const xAOD::CaloCluster * NeutralPFOClusterMLCorrectionTool::getLinkedCluster ( const xAOD::FlowElement & pfo) const
private

Definition at line 58 of file NeutralPFOClusterMLCorrectionTool.cxx.

59{
60 // Returns xAOD::Type::CaloCluster type link. There should be at most one such link.
61 // It can happen that no such link exists. This can happen due to negative cells subtraction.
62 // Empty link is returned if no valid cluster link is found.
63 const std::vector<ElementLink<xAOD::IParticleContainer>> &otherObjectLinks = pfo.otherObjectLinks();
64 if (otherObjectLinks.size() > 1)
65 ATH_MSG_ERROR("NeutralPFOClusterMLCorrectionTool: Multiple links found for neutral FlowElement with index " + std::to_string(pfo.index()));
66
67 bool hasValidLink = !otherObjectLinks.empty() && otherObjectLinks[0].isValid();
68 if (hasValidLink)
69 {
70 if ((**otherObjectLinks[0]).type() != xAOD::Type::CaloCluster)
71 ATH_MSG_ERROR("NeutralPFOClusterMLCorrectionTool: Link for neutral FlowElement with index " + std::to_string(pfo.index()) + " is not of type CaloCluster");
72
73 return static_cast<const xAOD::CaloCluster *>(*otherObjectLinks[0]);
74 }
75 else
76 return nullptr;
77}
#define ATH_MSG_ERROR(x)
size_t index() const
Return the index of this element within its container.
const std::vector< ElementLink< IParticleContainer > > & otherObjectLinks() const
@ CaloCluster
The object is a calorimeter cluster.
Definition ObjectType.h:39

◆ initialize()

StatusCode NeutralPFOClusterMLCorrectionTool::initialize ( )
overridevirtual

Definition at line 10 of file NeutralPFOClusterMLCorrectionTool.cxx.

10 {
11 ATH_MSG_DEBUG("Initializing with ClusterMLCorrectedEnergyDecorationKey: " << m_clusterMLCorrectedEnergyKey.value());
12 return StatusCode::SUCCESS;
13}
#define ATH_MSG_DEBUG(x)
Gaudi::Property< std::string > m_clusterMLCorrectedEnergyKey

◆ scaleEnergyToAlternativeSignalState()

void NeutralPFOClusterMLCorrectionTool::scaleEnergyToAlternativeSignalState ( xAOD::FlowElement & pfo,
const xAOD::CaloCluster & cls ) const
private

Definition at line 32 of file NeutralPFOClusterMLCorrectionTool.cxx.

33{
34 // Scale factor is defined as the ratio of the cluster energy stored in decoration
35 // to the energy in the EM calibration state (UNCALIBRATED). This scale factor is then applied to the PFO energy.
36 // If the EM energy is zero or if decoration is not found, no scaling is applied.
37 const SG::AuxElement::Accessor<double> clusterMLCorrectedEnergyAccessor(m_clusterMLCorrectedEnergyKey.value());
38 const double clusterEMEnergy = cls.rawE();
39 if (!clusterMLCorrectedEnergyAccessor.isAvailable(cls))
40 {
41 ATH_MSG_WARNING("NeutralPFOClusterMLCorrectionTool: No ML energy decoration '" << m_clusterMLCorrectedEnergyKey.value()
42 << "' found for cluster with index " << cls.index()
43 << ". PFO energy will not be scaled.");
44 }
45 const double clusterDecorEnergy = clusterMLCorrectedEnergyAccessor.isAvailable(cls) ? clusterMLCorrectedEnergyAccessor(cls) : clusterEMEnergy;
46 const double scaleFactor = clusterEMEnergy > FLT_MIN ? clusterDecorEnergy / clusterEMEnergy : 1.0;
47
48 ATH_MSG_DEBUG("NeutralPFOClusterMLCorrectionTool: Scaling PFO with index " << pfo.index()
49 << " energy from " << pfo.e() << " to " << (pfo.e() * scaleFactor)
50 << " using cluster index " << cls.index()
51 << " EM energy " << clusterEMEnergy
52 << " Decor energy " << clusterDecorEnergy
53 << " scale factor " << scaleFactor);
54
55 pfo.setP4(pfo.pt() * scaleFactor, pfo.eta(), pfo.phi(), pfo.m() * scaleFactor);
56}
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
virtual double pt() const override
void setP4(float pt, float eta, float phi, float m)
virtual double m() const override
The invariant mass of the particle.
virtual double phi() const override
The azimuthal angle ( ) of the particle.
virtual double eta() const override
The pseudorapidity ( ) of the particle.
virtual double e() const override
The total energy of the particle.

Member Data Documentation

◆ m_clusterMLCorrectedEnergyKey

Gaudi::Property<std::string> NeutralPFOClusterMLCorrectionTool::m_clusterMLCorrectedEnergyKey
private
Initial value:
{this, "ClusterMLCorrectedEnergyDecorationKey", "clusterE_ML",
"Name of the decoration storing the ML-corrected cluster energy"}

Definition at line 37 of file NeutralPFOClusterMLCorrectionTool.h.

37 {this, "ClusterMLCorrectedEnergyDecorationKey", "clusterE_ML",
38 "Name of the decoration storing the ML-corrected cluster energy"};

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