ATLAS Offline Software
Loading...
Searching...
No Matches
JetBadChanCorrTool.h
Go to the documentation of this file.
1// this file is -*- C++ -*-
2
3/*
4 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
5*/
6
7#ifndef XAOD_ANALYSIS
8
9#ifndef JETMOMENTTOOLS_JETBADCHANCORRTOOL_H
10#define JETMOMENTTOOLS_JETBADCHANCORRTOOL_H
11
33
34#include "AsgTools/AsgTool.h"
37
39
41
42#include "TH1D.h"
43#include <utility>
44
45
46class Identifier;
47
49class ITHistSvc;
50
52 virtual public IJetDecorator
53{
55public:
56 JetBadChanCorrTool(const std::string& name);
57
59
60 virtual StatusCode initialize() override;
61
62 virtual StatusCode decorate(const xAOD::JetContainer& jets) const override;
63
64 virtual StatusCode setupEvent();
65
66
67protected:
68 // These two apply the moments they compute as decorations to the jet collection
69 StatusCode correctionFromClustersBadCells(const xAOD::JetContainer& jets) const;
70 StatusCode correctionFromCellsInJet(const xAOD::JetContainer& jets, const jet::CaloCellFastMap * badCellMap) const;
71 // This one computes the moment without applying it
72 float correctionFromCellsInCone(const xAOD::Jet* jet, const jet::CaloCellFastMap * badCellMap) const;
73
74 private:
75
77
78 Gaudi::Property<std::string> m_jetContainerName{this, "JetContainer", "", "SG key for the input jet container"};
79 Gaudi::Property<int> m_nBadCellLimit{this, "NBadCellLimit", 10000, "Limit to calculate moments"};
80
81 // for jet-level correction
82 Gaudi::Property<std::string> m_streamName{this, "StreamName", "/JetBadChanCorrTool/", "Stream name"};
83 Gaudi::Property<std::string> m_profileName{this, "ProfileName", "JetBadChanCorrTool.root", "Profile name"};
84 Gaudi::Property<std::string> m_profileTag{this, "ProfileTag", "", "Profile tag"};
85
86 Gaudi::Property<bool> m_useCone{this, "UseCone", true, "Use cone?"};
87 Gaudi::Property<bool> m_useClusters{this, "UseClusters", false, "Use clusters?"};
88
89 SG::ReadHandleKey<jet::CaloCellFastMap> m_badCellMap_key{this, "MissingCellMap", "MissingCaloCellsMap", "SG key for missing cell map"};
90 SG::WriteDecorHandleKey<xAOD::JetContainer> m_corrCellKey{this, "CorrCellDecorKey", "BchCorrCell", "SG key for cell level decoration"};
91 SG::WriteDecorHandleKey<xAOD::JetContainer> m_corrDotxKey{this, "CorrDotxDecorKey", "BchCorrDotx", "SG key for DOTX decoration"};
92 SG::WriteDecorHandleKey<xAOD::JetContainer> m_corrJetKey{this, "CorrJetDecorKey", "BchCorrJet", "SG key for jet level decoration"};
93 SG::WriteDecorHandleKey<xAOD::JetContainer> m_corrJetForCellKey{this, "CorrJetForCellDecorKey", "BchCorrJetForCell", "SG key for JetForCell decoration"};
94
95 // jet profiles
97 public:
98 ProfileData(TH1* th, int sample,
99 double ptMin=0, double ptMax=9999,
100 double etaMin=0, double etaMax=5.0,
101 double phiMin=-M_PI, double phiMax=M_PI):
102 m_th(th), m_sample(sample),
103 m_ptMin(ptMin), m_ptMax(ptMax),
104 m_etaMin(etaMin), m_etaMax(etaMax),
105 m_phiMin(phiMin), m_phiMax(phiMax) {}
106
107 virtual ~ProfileData(){}
108
109 bool match(double pt, int sample, double eta, double phi) const {
110 return ( pt>=m_ptMin && pt<m_ptMax
111 && sample==m_sample
112 && fabs(eta)>=m_etaMin && fabs(eta)<m_etaMax
113 && phi>=m_phiMin && phi<m_phiMax);
114 }
115
116 double frac(double dr) const {
117 int idr = std::as_const(m_th)->FindBin(dr);
118 return m_th->GetBinContent(idr);
119 }
120 private:
121 TH1* m_th;
123 double m_ptMin;
124 double m_ptMax;
125 double m_etaMin;
126 double m_etaMax;
127 double m_phiMin;
128 double m_phiMax;
129
130 };
131 std::vector<ProfileData> m_profileDatas[CaloCell_ID::Unknown];//24
132
133 double getProfile(double pt, double dr, int sample, double eta, double phi) const;
134};
135#endif
136
137#endif
138
139// DoxygenDocumentation
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ASG_TOOL_CLASS0(CLASSNAME)
Property holding a SG store/key/clid/attr name from which a WriteDecorHandle is made.
Interface for adding a decoration to a jet container.
An abstract inteface to get Tile channel and ADC status.
ProfileData(TH1 *th, int sample, double ptMin=0, double ptMax=9999, double etaMin=0, double etaMax=5.0, double phiMin=-M_PI, double phiMax=M_PI)
bool match(double pt, int sample, double eta, double phi) const
SG::WriteDecorHandleKey< xAOD::JetContainer > m_corrJetKey
virtual StatusCode setupEvent()
Gaudi::Property< int > m_nBadCellLimit
Gaudi::Property< std::string > m_jetContainerName
virtual ~JetBadChanCorrTool()
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Gaudi::Property< bool > m_useCone
std::vector< ProfileData > m_profileDatas[CaloCell_ID::Unknown]
float correctionFromCellsInCone(const xAOD::Jet *jet, const jet::CaloCellFastMap *badCellMap) const
SG::WriteDecorHandleKey< xAOD::JetContainer > m_corrCellKey
double getProfile(double pt, double dr, int sample, double eta, double phi) const
SG::WriteDecorHandleKey< xAOD::JetContainer > m_corrDotxKey
StatusCode correctionFromClustersBadCells(const xAOD::JetContainer &jets) const
SG::WriteDecorHandleKey< xAOD::JetContainer > m_corrJetForCellKey
SG::ReadHandleKey< jet::CaloCellFastMap > m_badCellMap_key
Gaudi::Property< std::string > m_streamName
Gaudi::Property< std::string > m_profileTag
JetBadChanCorrTool(const std::string &name)
StatusCode correctionFromCellsInJet(const xAOD::JetContainer &jets, const jet::CaloCellFastMap *badCellMap) const
Gaudi::Property< std::string > m_profileName
Gaudi::Property< bool > m_useClusters
ServiceHandle< ITHistSvc > m_thistSvc
virtual StatusCode decorate(const xAOD::JetContainer &jets) const override
Decorate a jet collection without otherwise modifying it.
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".