ATLAS Offline Software
Loading...
Searching...
No Matches
JetForwardJvtTool.h
Go to the documentation of this file.
1
2
3/*
4 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
5*/
6
7// JetForwardJvtTool.h
8// Header file for class JetForwardJvtTool
9// Author: Matt Klein<matthew.henry.klein@cern.ch>
11#ifndef FORWARDJVTTOOL_JVT_FORWARDJVTTOOL_H
12#define FORWARDJVTTOOL_JVT_FORWARDJVTTOOL_H 1
13
14// STL includes
15#include <string>
16
17
18// FrameWork includes
19#include "AsgTools/ToolHandle.h"
20#include "AsgTools/AsgTool.h"
24
25// EDM includes
31#include "AsgTools/IAsgTool.h"
32
33//ASG_TOOL_INTERFACE(IJetUpdateJvt)
34
36 : public asg::AsgTool,
37 virtual public IJetDecorator{
39
40 // This macro defines the constructor with the interface declaration
41 //ASG_TOOL_CLASS(JetForwardJvtTool, IJetForwardJvtTool)
42
43
44 // Public methods:
46 public:
47
48 // Copy constructor:
49
51 JetForwardJvtTool(const std::string& name);
52
55
56 // Athena algtool's Hooks
57 virtual StatusCode initialize() override;
58
59 virtual StatusCode decorate(const xAOD::JetContainer& jetCont) const override;
60
61 float getFJVT(const xAOD::Jet *jet, const std::vector<TVector2>& pileupMomenta, std::size_t pvind) const;
62 bool forwardJet(const xAOD::Jet *jet) const;
63 bool centralJet(const xAOD::Jet *jet) const;
64 float getDrpt(const xAOD::Jet *jet) const;
65 int getJetVertex(const xAOD::Jet *jet) const;
66
67 StatusCode tagTruth(const xAOD::JetContainer *jets,const xAOD::JetContainer *truthJets);
68 std::vector<TVector2> calculateVertexMomenta(const xAOD::JetContainer *jets, std::size_t pvind) const;
69 float getCombinedWidth(const xAOD::Jet *jet) const;
70
71 private:
72
73 Gaudi::Property<double> m_etaThresh{this, "EtaThresh", 2.5, "Eta threshold"};
74 Gaudi::Property<double> m_forwardMinPt{this, "ForwardMinPt", 20e3, "Forward minimum Pt"};
75 Gaudi::Property<double> m_forwardMaxPt{this, "ForwardMaxPt", 50e3, "Forward maximum Pt"};
76 Gaudi::Property<double> m_centerMinPt{this, "CentralMinPt", 20e3, "Central minimum Pt"};
77 Gaudi::Property<double> m_centerMaxPt{this, "CentralMaxPt", -1, "Central maximum Pt"};
78 Gaudi::Property<double> m_centerJvtThresh{this, "CentralJvtThresh", 0.14, "Central JVT threshold"};
79 Gaudi::Property<double> m_centerDrptThresh{this, "CentralDrptThresh", 0.2, "Central drpt threshold"};
80 Gaudi::Property<double> m_maxStochPt{this, "CentralMaxStochPt", 35e3, "Central maximum StochPt"};
81 Gaudi::Property<double> m_jetScaleFactor{this, "JetScaleFactor", 0.4, "Jet scale factor"};
82 Gaudi::Property<double> m_fjvtThresh{this, "FjvtThresh", 15e3, "FJVT threshold"}; //15GeV->92%,11GeV->85%
83 Gaudi::Property<bool> m_tightOP{this, "UseTightOP", false, "Use tight (true) or loose (false)"};
84 Gaudi::Property<bool> m_recalculateFjvt{this, "RecalculateFjvt", true, "Recalculate Fjvt or use stored value"};
85 std::size_t getPV() const;
86
89
90 Gaudi::Property<bool> m_renounceOutputs{this, "RenounceOutputs", false, "If true, then we're not interested in output dependencies. In that case, JetContainer is irrelevant."};
91 Gaudi::Property<std::string> m_jetContainerName{this, "JetContainer", "", "SG key for the input jet container"};
92 SG::ReadDecorHandleKey<xAOD::JetContainer> m_orKey{this, "OverlapDec", "", "Overlap decoration key"};
93 SG::WriteDecorHandleKey<xAOD::JetContainer> m_outKey{this, "OutputDec", "passFJVT", "Output decoration key"};
94 SG::WriteDecorHandleKey<xAOD::JetContainer> m_isHSKey{this, "IsJvtHSName", "isJvtHS", "Decoration key for isJvtHS"};
95 SG::WriteDecorHandleKey<xAOD::JetContainer> m_isPUKey{this, "IsJvtPUName", "isJvtPU", "Decoration key for isJvtPU"};
96 SG::WriteDecorHandleKey<xAOD::JetContainer> m_fjvtDecKey{this, "FJVTName", "fJvt", "Decoration key for fJvt"};
97
98 SG::ReadHandleKey<xAOD::VertexContainer> m_vertexContainerName{this, "VertexContainerName", "PrimaryVertices", "SG key for vertex container"};
99 SG::ReadHandleKey<xAOD::MissingETContainer> m_trkMETName{this, "Met_TrackName", "MET_Track", "SG key for MET track container"};
100
101 SG::ReadDecorHandleKey<xAOD::JetContainer> m_widthKey{this, "WidthName", "Width", "SG key for jet width"};
102 SG::ReadDecorHandleKey<xAOD::JetContainer> m_jvtMomentKey{this, "JvtMomentName", "Jvt", "JVT moment name"};
103 SG::ReadDecorHandleKey<xAOD::JetContainer> m_sumPtsKey{this, "SumPtsName", "SumPtTrkPt500", "SG key for SumPt vector"};
104
105 };
106#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.
Interface for adding a decoration to a jet container.
Gaudi::Property< bool > m_recalculateFjvt
Gaudi::Property< double > m_maxStochPt
Gaudi::Property< double > m_centerJvtThresh
virtual ~JetForwardJvtTool()
Destructor:
SG::ReadDecorHandleKey< xAOD::JetContainer > m_orKey
SG::ReadDecorHandleKey< xAOD::JetContainer > m_widthKey
SG::ReadDecorHandleKey< xAOD::JetContainer > m_jvtMomentKey
Gaudi::Property< double > m_etaThresh
SG::ReadHandleKey< xAOD::MissingETContainer > m_trkMETName
JetForwardJvtTool()
Default constructor:
float getCombinedWidth(const xAOD::Jet *jet) const
Gaudi::Property< bool > m_renounceOutputs
Gaudi::Property< bool > m_tightOP
StatusCode tagTruth(const xAOD::JetContainer *jets, const xAOD::JetContainer *truthJets)
Gaudi::Property< double > m_fjvtThresh
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Gaudi::Property< double > m_centerMinPt
SG::WriteDecorHandleKey< xAOD::JetContainer > m_isPUKey
bool forwardJet(const xAOD::Jet *jet) const
bool centralJet(const xAOD::Jet *jet) const
Gaudi::Property< double > m_centerDrptThresh
SG::WriteDecorHandleKey< xAOD::JetContainer > m_fjvtDecKey
Gaudi::Property< double > m_forwardMaxPt
std::vector< TVector2 > calculateVertexMomenta(const xAOD::JetContainer *jets, std::size_t pvind) const
SG::WriteDecorHandleKey< xAOD::JetContainer > m_outKey
SG::ReadDecorHandleKey< xAOD::JetContainer > m_sumPtsKey
virtual StatusCode decorate(const xAOD::JetContainer &jetCont) const override
Decorate a jet collection without otherwise modifying it.
JetForwardJvtTool(const std::string &name)
Constructor with parameters:
int getJetVertex(const xAOD::Jet *jet) const
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexContainerName
Gaudi::Property< double > m_jetScaleFactor
Gaudi::Property< double > m_forwardMinPt
float getDrpt(const xAOD::Jet *jet) const
float getFJVT(const xAOD::Jet *jet, const std::vector< TVector2 > &pileupMomenta, std::size_t pvind) const
Gaudi::Property< std::string > m_jetContainerName
std::size_t getPV() const
Gaudi::Property< double > m_centerMaxPt
SG::WriteDecorHandleKey< xAOD::JetContainer > m_isHSKey
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".
JetContainer_v1 JetContainer
Definition of the current "jet container version".