ATLAS Offline Software
Loading...
Searching...
No Matches
JetForwardPFlowJvtTool.h
Go to the documentation of this file.
1
2
3/*
4 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
5*/
6
7// JetForwardPFlowJvtTool.h
8// Header file for class JetForwardPFlowJvtTool
9// Author: Anastasia Kotsokechagia <anastasia.kotsokechagia@cern.ch>
10
11// Tool for calculating fjvt values for pflow jets.
12// Short describtion of the tool;
13// First central PU jets are built per vertex.
14// Reconstructed calibrated jets are then used to calculate the per vertex missing momentum (miss-mom).
15// The per vertex missing momentum is defined as: The vector some of the calibrated jet momenta (for jets with pt>20GeV && Rpt>0.1 wrt to the vertex) + tracks assosiated to the vertex (otherwise).
16// PU Jets closeby (dR<0.3) to a HS jet are not considered.
17// The fJVT value for every forward jet (fj) of the event is then calculated choosing the vertex with the largest negative miss-mom projection on the fj.
18// User action: After initializing the tool the user has to call the modify(xAOD::JetContainer& jetCont) function. Argument in this fuction is the PFlow jet container of the event.
19// The fjvt value for every forward jet of the container is then calculated and can be retrieved.
21 //Parameters
22 // m_orLabel: ""
23 // m_jetsName : "Container name for the output reconstructed PU jets "
24 // m_tightOP: "If true a tight fjvt threshold value is applied"
25 // m_outLabelFjvt: "Decorator for passing fJVT threshold (tight or loose)"
26 // m_jetchargedp4: "Name of the jet charged momentum 4-vector"
27 // m_etaThresh: "Maximum eta value for considering a jet as central"
28 // m_forwardMinPt: "Minimum forward jet pt"
29 // m_forwardMaxPt: "Maximum forward jet pt. If -1 no threshold is applied"
30 // m_centerMinPt: "Minimum central jet pt"
31 // m_centerMaxPt: "Maximum central jet pt. If -1 no threshold is applied"
32 // m_pvind: "Hard-Scatter primary vertex index of the event. If -1 it's automatically retrieved from the event"
33 // m_rptCut: "Rpt cut value for central PU jets contributing in the missing momentum calculation"
34 // m_jvtCut: "JVT threshold value for considering a central PU jet as HS"
35 // m_dzCut: "Dz=z-z0 cut value for pfo objects participating in the HS vertex jet reco"
36 // m_vertices: "Number of vertices for which the missing momentum is calculated"
37 // m_maxRap: "Maximum rapidity value in fastjet::AreaDefinition"
38 // m_neutMaxRap: "Maximum rapidity value for neutral pfos participating in jet reco"
39 // m_weight: "PFO weight value"
40 // m_pfoToolName: "Name of PFO retriever tool"
41 // m_wpfoToolName: "Name of PFO weighting tool"
42 // m_pfoJESName: "Name of jet calibration tool"
43 // m_jetAlgo: "Jet calibration collection name"
44 // m_calibconfig: "Calibration config for PFlow jets, need to be updated with latest one"
45 // m_calibSeq: "Calibration sequence to be applied"
46 // m_calibArea: "Calibration area"
47 // m_isdata: "True if data"
48
49
50#ifndef FORWARDPFLOWJVTTOOL_JVT_FORWARDPFLOWJVTTOOL_H
51#define FORWARDPFLOWJVTTOOL_JVT_FORWARDPFLOWJVTTOOL_H 1
52
53// STL includes
54#include <string>
55
56// FrameWork includes
57#include "AsgTools/ToolHandle.h"
58#include "AsgTools/AsgTool.h"
66
67
72
73// Pflow / FE tools
78
79#include "AsgTools/ToolHandle.h"
81
82namespace pflow {
83 struct puJets {
84 std::shared_ptr<xAOD::JetContainer> jetCont;
85 std::shared_ptr<xAOD::JetAuxContainer> jetAuxCont;
86 };
87}
88
90 : public asg::AsgTool,
91 virtual public IJetDecorator{
93
94
95 // Public methods:
97 public:
98
100 JetForwardPFlowJvtTool(const std::string& name);
101
104
105 virtual StatusCode initialize() override;
106
107
108 virtual StatusCode decorate(const xAOD::JetContainer& jetCont) const override;
109
110 float getFJVT(const xAOD::Jet *jet,const std::vector<TVector2>& pileupMomenta) const;
111 bool isForwardJet(const xAOD::Jet *jet) const;
112 bool isCentralJet(const xAOD::Jet *jet) const;
113
114 StatusCode tagTruth(const xAOD::JetContainer *jets,const xAOD::JetContainer *truthJets);
115 virtual std::vector<TVector2> calculateVertexMomenta(const xAOD::JetContainer *jets,int pvind, int vertices) const;
117 bool hasCloseByHSjet(const xAOD::Jet *jet, const xAOD::JetContainer *pjets ) const;
118 double getRpt(const xAOD::Jet *jet) const;
119 fastjet::PseudoJet pfoToPseudoJet(const xAOD::PFO* pfo, const CP::PFO_JetMETConfig_charge& theCharge, const xAOD::Vertex *vx) const;
120 fastjet::PseudoJet feToPseudoJet(const xAOD::FlowElement* fe, const CP::PFO_JetMETConfig_charge& theCharge, const xAOD::Vertex *vx) const;
121 std::size_t getPV() const;
122
123 protected:
124
125 SG::ReadHandleKey<jet::TrackVertexAssociation> m_tvaKey{this, "TrackVertexAssociation", "", "Input track-vertex association"};
126 Gaudi::Property<std::string> m_jetContainerName{this, "JetContainer", "", "SG key for the input jet container"};
127 Gaudi::Property<std::string> m_jetsName{this, "jetsName", "AntiKt4PUPFlowJets", "Container name for the output reconstructed PU jets"};
128 Gaudi::Property<std::string> m_jetchargedp4{this, "jetchargedp4", "JetChargedScaleMomentum", "Name of the jet charged momentum 4-vector"};
129
130 Gaudi::Property<bool> m_isdata{this, "isdata", false, "True if data"};
131 Gaudi::Property<int> m_pvind{this, "pvind", -1, "Hard-Scatter primary vertex index of the event. If -1 it will be automatically retrieved from the event"};
132 Gaudi::Property<int> m_vertices{this, "vertices", 10, "Number of vertices for which the missing momentum is calculated"};
133 Gaudi::Property<bool> m_includePV{this, "includePV", false, "Flag to include jets and tracks associated to PV in the calculation"};
134 Gaudi::Property<double> m_etaThresh{this, "etaThresh", 2.5, "Maximum eta value for considering a jet as central"};
135 Gaudi::Property<double> m_forwardMinPt{this, "forwardMinPt", 18e3, "Minimum forward jet pt"};
136 Gaudi::Property<double> m_forwardMaxPt{this, "forwardMaxPt", -1, "Maximum forward jet pt. If -1 no threshold is applied"};
137 Gaudi::Property<double> m_centerMinPt{this, "centralMinPt", 20e3, "Minimum central jet pt"};
138 Gaudi::Property<double> m_centerMaxPt{this, "centralMaxPt", -1, "Maximum central jet pt. If -1 no threshold is applied"};
139 Gaudi::Property<double> m_fjvtThresh{this, "fjvtThresh", 15e3, "fjvt threshold value"};
140 Gaudi::Property<double> m_rptCut{this, "rptCut", 0.1, "Rpt cut value for central PU jets contributing in the missing momentum calculation"};
141 Gaudi::Property<double> m_jvtCut{this, "jvtCut", 0.2, "JVT threshold value for considering a central PU jet as HS"};
142 Gaudi::Property<double> m_dzCut{this, "dzCut", 2.0, "Dz=z=-z0 cut for pfo objects participating in the HS vertex jet reco"};
143 Gaudi::Property<double> m_maxRap{this, "maxRap", 2.5, "Maximum rapidity value in fastjet::AreaDefinition"};
144 Gaudi::Property<double> m_neutMaxRap{this, "neutMaxRap", 2.5, "Maximum rapidity value for neutral pfos participating in jet reco"};
145 Gaudi::Property<float> m_weight{this, "weight", 0, "PFO weight value"};
146 Gaudi::Property<bool> m_tightOP{this, "tightOP", false, "If true a tight fjvt threshold value is applied"};
147
148 // not used?
149 //Gaudi::Property<std::string> m_jvtMomentName{"jvtMomentName", "", ""};
150 //Gaudi::Property<double> m_centerJvtThresh{"", 0, ""};
151
152 SG::ReadHandleKey<xAOD::VertexContainer> m_vxContKey{this, "verticesName", "PrimaryVertices", "Container name of vertices to be retrieved"};
153 SG::ReadHandleKey<xAOD::PFOContainer> m_PFOKey{this, "PFOName", "CHSParticleFlowObjects", "SG Key for CHS PFO Container"};
154 SG::ReadHandleKey<xAOD::FlowElementContainer> m_FEKey{this, "FEName", "", "SG Key for CHS FlowElement Container (overrides PFO if not empty)"};
156
157 SG::ReadDecorHandleKey<xAOD::JetContainer> m_passJvtKey{this, "passJvtName", "NNJvtPass", "SG key for output pass-JVT decoration"};
158 SG::ReadDecorHandleKey<xAOD::PFO> m_orKey{this, "ORName", "", "OR label"};
159 SG::ReadDecorHandleKey<xAOD::FlowElement> m_orFEKey{this, "ORNameFE", "", "OR label"};
160
161 SG::WriteDecorHandleKey<xAOD::JetContainer> m_fjvtKey{this, "FjvtName", "passOnlyFJVT", "Decorator for passing fJVT threshold (tight or loose)"};
162 SG::WriteDecorHandleKey<xAOD::JetContainer> m_fjvtRawKey{this, "FjvtRawName", "fJvt", "Decorator for raw fJVT variable"};
163 SG::WriteDecorHandleKey<xAOD::JetContainer> m_isHSKey{this, "isHSName", "isJVTHS", "SG key for output isJVTHS decoration"};
164 SG::WriteDecorHandleKey<xAOD::JetContainer> m_isPUKey{this, "isPUName", "isJvtPU", "SG key for output isJVTPU decoration"};
165
166 ToolHandle<CP::WeightPFOTool> m_wpfotool{this,"WeightPFOTool", "", "Weight PFO tool name"};
167 ToolHandle<IJetCalibrationTool> m_pfoJES{this,"JetCalibrationTool", "", "Jet calibration tool name"};
168
169 Gaudi::Property<bool> m_suppressInputDependence{this, "SuppressInputDependence", false};
170
171 };
172#endif //> !FORWARDJVTTOOL_JVT_FORWARDJVTTOOL_H
#define ASG_TOOL_CLASS(CLASSNAME, INT1)
Property holding a SG store/key/clid/attr name from which a ReadDecorHandle is made.
Handle class for reading a decoration on an object.
Handle class for adding a decoration to an object.
Interface for adding a decoration to a jet container.
SG::ReadDecorHandleKey< xAOD::PFO > m_orKey
Gaudi::Property< bool > m_suppressInputDependence
SG::WriteDecorHandleKey< xAOD::JetContainer > m_fjvtKey
Gaudi::Property< double > m_etaThresh
bool hasCloseByHSjet(const xAOD::Jet *jet, const xAOD::JetContainer *pjets) const
Gaudi::Property< double > m_rptCut
SG::ReadDecorHandleKey< xAOD::JetContainer > m_passJvtKey
Gaudi::Property< double > m_forwardMinPt
Gaudi::Property< double > m_jvtCut
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Gaudi::Property< std::string > m_jetContainerName
Gaudi::Property< double > m_centerMaxPt
fastjet::PseudoJet feToPseudoJet(const xAOD::FlowElement *fe, const CP::PFO_JetMETConfig_charge &theCharge, const xAOD::Vertex *vx) const
SG::ReadHandleKey< xAOD::VertexContainer > m_vxContKey
Gaudi::Property< bool > m_isdata
SG::WriteDecorHandleKey< xAOD::JetContainer > m_isHSKey
SG::ReadHandleKey< xAOD::FlowElementContainer > m_FEKey
virtual std::vector< TVector2 > calculateVertexMomenta(const xAOD::JetContainer *jets, int pvind, int vertices) const
Gaudi::Property< std::string > m_jetchargedp4
Gaudi::Property< float > m_weight
pflow::puJets buildPFlowPUjets(const xAOD::Vertex &vx) const
Gaudi::Property< double > m_centerMinPt
Gaudi::Property< bool > m_includePV
Gaudi::Property< double > m_fjvtThresh
ToolHandle< CP::WeightPFOTool > m_wpfotool
float getFJVT(const xAOD::Jet *jet, const std::vector< TVector2 > &pileupMomenta) const
SG::ReadDecorHandleKey< xAOD::FlowElement > m_orFEKey
Gaudi::Property< double > m_dzCut
ToolHandle< IJetCalibrationTool > m_pfoJES
Gaudi::Property< double > m_maxRap
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
Gaudi::Property< int > m_pvind
bool isForwardJet(const xAOD::Jet *jet) const
virtual ~JetForwardPFlowJvtTool()
Destructor:
Gaudi::Property< bool > m_tightOP
bool isCentralJet(const xAOD::Jet *jet) const
StatusCode tagTruth(const xAOD::JetContainer *jets, const xAOD::JetContainer *truthJets)
SG::ReadHandleKey< xAOD::PFOContainer > m_PFOKey
Gaudi::Property< std::string > m_jetsName
SG::WriteDecorHandleKey< xAOD::JetContainer > m_isPUKey
virtual StatusCode decorate(const xAOD::JetContainer &jetCont) const override
Decorate a jet collection without otherwise modifying it.
Gaudi::Property< double > m_forwardMaxPt
Gaudi::Property< int > m_vertices
fastjet::PseudoJet pfoToPseudoJet(const xAOD::PFO *pfo, const CP::PFO_JetMETConfig_charge &theCharge, const xAOD::Vertex *vx) const
double getRpt(const xAOD::Jet *jet) const
JetForwardPFlowJvtTool(const std::string &name)
Constructor with parameters:
Gaudi::Property< double > m_neutMaxRap
SG::WriteDecorHandleKey< xAOD::JetContainer > m_fjvtRawKey
SG::ReadHandleKey< jet::TrackVertexAssociation > m_tvaKey
Property holding a SG store/key/clid/attr name from which a ReadDecorHandle is made.
Property holding a SG store/key/clid from which a ReadHandle is made.
Property holding a SG store/key/clid/attr name from which a WriteDecorHandle is made.
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
Jet_v1 Jet
Definition of the current "jet version".
PFO_v1 PFO
Definition of the current "pfo version".
Definition PFO.h:17
FlowElement_v1 FlowElement
Definition of the current "pfo version".
Definition FlowElement.h:16
Vertex_v1 Vertex
Define the latest version of the vertex class.
JetContainer_v1 JetContainer
Definition of the current "jet container version".
std::shared_ptr< xAOD::JetAuxContainer > jetAuxCont
std::shared_ptr< xAOD::JetContainer > jetCont