ATLAS Offline Software
Loading...
Searching...
No Matches
GSCCalibStep.cxx
Go to the documentation of this file.
1
2
3/*
4 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
5*/
6
7// GSCCalibStep.cxx
8// Implementation file for class GSCCalibStep
10
12
13GSCCalibStep::GSCCalibStep(const std::string& name)
14 : asg::AsgTool( name ){ }
15
17// Public methods:
19
21 ATH_MSG_DEBUG ("Initializing " << name() );
22
23 ATH_MSG_DEBUG("Reading from " << m_jetInScale << " and writing to " << m_jetOutScale);
24
25 ATH_CHECK( m_histTool_EM3.retrieve());
27 ATH_CHECK( m_histTool_Tile0.retrieve());
29 ATH_CHECK( m_histTool_nTrk.retrieve());
31
32 ATH_CHECK(m_vertexContainer_key.initialize());
33
34 return StatusCode::SUCCESS;
35}
36
37
39
40 ATH_MSG_DEBUG("calibrating jet collection.");
41
42 // Retrieve the primary vertex location:
44 const xAOD::VertexContainer& vertices = *vertexHandle;
45 const xAOD::Vertex *HSvertex = findHSVertex(vertices);
46 if(!HSvertex) {
47 ATH_MSG_WARNING("Invalid primary vertex found, will not continue applying the GSC.");
48 return StatusCode::FAILURE;
49 }
50 int PVindex = HSvertex->index();
51 ATH_MSG_DEBUG("PV index:" << PVindex);
52
53 // Calibrate the jets
54 for (xAOD::Jet* jet : jets){
55
57
58 xAOD::JetFourMom_t jetconstitP4 = jet->getAttribute<xAOD::JetFourMom_t>("JetConstitScaleMomentum");
59 std::vector<float> samplingFrac = jet->getAttribute<std::vector<float> >("EnergyPerSampling");
60 // get detector Eta
61 float detectorEta = jet->getAttribute<float>("DetectorEta");
62 // get trackWIDTHPVX
63 float trackWIDTHPVX = 0;
64 static const SG::ConstAccessor<std::vector<float> > TrackWidthPt1000Acc ("TrackWidthPt1000");
65 if(TrackWidthPt1000Acc.isAvailable(*jet))
66 {
67 trackWIDTHPVX = TrackWidthPt1000Acc(*jet).at(PVindex);
68 ATH_MSG_DEBUG("trackWIDTHPVX found set to: " << trackWIDTHPVX);
69 }
70 jc.setValue("trackWIDTH", trackWIDTHPVX);
71 // get nTrkPVX
72 int nTrkPVX = 0;
73 static const SG::ConstAccessor<std::vector<int>> NumTrkPt1000Acc ("NumTrkPt1000");
74 if(NumTrkPt1000Acc.isAvailable(*jet))
75 {
76 nTrkPVX = NumTrkPt1000Acc(*jet).at(PVindex);
77 ATH_MSG_DEBUG("nTrkPVX found set to: " << nTrkPVX);
78 }
79 jc.setValue("nTrk", nTrkPVX);
80 // get Charged Fraction
81 float ChargedFraction = 0;
82 static const SG::ConstAccessor<std::vector<float>> SumPtChargedPFOPt500Acc ("SumPtChargedPFOPt500");
83 if(SumPtChargedPFOPt500Acc.isAvailable(*jet))
84 {
85 ChargedFraction = SumPtChargedPFOPt500Acc(*jet).at(PVindex)/jetconstitP4.Pt();
86 ATH_MSG_DEBUG("ChargedFraction found set to: " << ChargedFraction);
87 }
88 jc.setValue("ChargedFraction", ChargedFraction);
89 // get EM3
90 float EM3 = (samplingFrac[3]+samplingFrac[7])/jetconstitP4.e();
91 ATH_MSG_DEBUG("EM3 found set to: " << EM3);
92 jc.setValue("EM3", EM3);
93 // get Tile0
94 float Tile0 = (samplingFrac[12]+samplingFrac[18])/jetconstitP4.e();
95 ATH_MSG_DEBUG("Tile0 found set to: " << Tile0);
96 jc.setValue("Tile0", Tile0);
97 // get N90Constituents
98 double N90Constituents = 0;
99 static const SG::ConstAccessor<float> N90ConstituentsAcc ("N90Constituents");
100 if(N90ConstituentsAcc.isAvailable(*jet))
101 {
102 N90Constituents = N90ConstituentsAcc(*jet);
103 ATH_MSG_DEBUG("N90Constituents found set to: " << N90Constituents);
104 }
105 jc.setValue("N90Constituents", N90Constituents);
106 // get caloWIDTH
107 double caloWIDTH = 0;
108 static const SG::ConstAccessor<float> WidthAcc ("Width");
109 if(WidthAcc.isAvailable(*jet))
110 {
111 caloWIDTH = WidthAcc(*jet);
112 ATH_MSG_DEBUG("caloWIDTH found set to: " << caloWIDTH);
113 }
114 jc.setValue("caloWIDTH", caloWIDTH);
115 // get TG3
116 float TG3 = (samplingFrac[17])/jetconstitP4.e();
117 ATH_MSG_DEBUG("TG3 found set to: " << TG3);
118 jc.setValue("TG3", TG3);
119 // get Muon segments
120 int Nsegments = 0;
121 static const SG::ConstAccessor<int> GhostMuonSegmentCountAcc ("GhostMuonSegmentCount");
122 if(GhostMuonSegmentCountAcc.isAvailable(*jet))
123 {
124 Nsegments = GhostMuonSegmentCountAcc(*jet);
125 ATH_MSG_DEBUG("Nsegments found set to: " << Nsegments);
126 }
127 jc.setValue("Nsegments", Nsegments);
128
129 float getGSCCorrection = 1.0;
130 int etabin = std::abs(detectorEta)/0.1;// m_binSize in old version
131
132 const xAOD::JetFourMom_t startingP4 = jet->getAttribute<xAOD::JetFourMom_t>(m_jetInScale);
133 jet->setJetP4(startingP4);
134
135 ATH_MSG_DEBUG("Jet pt original ("<<m_jetInScale<<"): " << jet->pt()*1e-3);
136
137 ATH_MSG_DEBUG("ChargedFraction Response: " << getChargedFractionResponse(*jet, jc, etabin));
138 ATH_MSG_DEBUG("Tile0 Response: " <<getTile0Response(*jet, jc, etabin));
139 ATH_MSG_DEBUG("EM3 Response: " << getEM3Response(*jet, jc, etabin));
140 ATH_MSG_DEBUG("NTrk Response: " <<getNTrkResponse(*jet, jc, etabin));
141 ATH_MSG_DEBUG("TrkWidth Response: " <<getTrackWIDTHResponse(*jet, jc, etabin));
142
143 getGSCCorrection*=1./getChargedFractionResponse(*jet, jc, etabin);
144 jet->setJetP4( startingP4*getGSCCorrection );
145 getGSCCorrection*=1./getTile0Response(*jet, jc, etabin);
146 jet->setJetP4( startingP4*getGSCCorrection );
147 getGSCCorrection*=1./getEM3Response(*jet, jc, etabin);
148 jet->setJetP4( startingP4*getGSCCorrection );
149 getGSCCorrection*=1./getNTrkResponse(*jet, jc, etabin);
150 jet->setJetP4( startingP4*getGSCCorrection );
151 getGSCCorrection*=1./getTrackWIDTHResponse(*jet, jc, etabin);
152
153 if(m_applyPunchThrough && startingP4.Pt() >= m_punchThroughMinPt){
154 jet->setJetP4( startingP4*getGSCCorrection );
155 getGSCCorrection*=1./getPunchThroughResponse(*jet, jc, std::abs(detectorEta));
156 }
157
158 ATH_MSG_DEBUG("GSC full correction: " << getGSCCorrection);
159
160 jet->setAttribute<xAOD::JetFourMom_t>(m_jetOutScale,startingP4*getGSCCorrection);
161 jet->setJetP4( startingP4*getGSCCorrection );
162
163 ATH_MSG_DEBUG("Jet pt calibrated:" << jet->pt()*1e-3);
164
165
166 }// loop jets
167
168 return StatusCode::SUCCESS;
169}
170
172 if (jc.getValue<float>("ChargedFraction")<=0) return 1; //ChargedFraction < 0 is unphysical, ChargedFraction = 0 is a special case, so we return 1 for ChargedFraction <= 0
173 if ( etabin >= m_histTool_ChargedFraction.size() ) return 1.;
174 double ChargedFractionResponse = m_histTool_ChargedFraction[etabin]->getValue(jet, jc);
175 return ChargedFractionResponse;
176}
177
179 if (jc.getValue<float>("Tile0")<0) return 1; //Tile0 < 0 is unphysical, so we return 1
180 if ( etabin >= m_histTool_Tile0.size() ) return 1.;
181 double Tile0Response = m_histTool_Tile0[etabin]->getValue(jet, jc);
182 return Tile0Response;
183}
184
186 if (jc.getValue<float>("EM3")<=0) return 1; //EM3 < 0 is unphysical, EM3 = 0 is a special case, so we return 1 for EM3 <= 0
187 if ( etabin >= m_histTool_EM3.size() ) return 1.;
188 float EM3Response = m_histTool_EM3[etabin]->getValue(jet, jc);
189 return EM3Response;
190}
191
192float GSCCalibStep::getPunchThroughResponse(const xAOD::Jet& jet, const JetHelper::JetContext& jc, double eta_det) const {
193 int etabin=-99;
195 ATH_MSG_WARNING("Please check that the punch through eta binning is properly set in your config file");
196 if ( eta_det >= m_punchThroughEtaBins[m_punchThroughEtaBins.size()-1] || jc.getValue<int>("Nsegments") < 20 ) return 1;
197 for (uint i=0; i<m_punchThroughEtaBins.size()-1; ++i) {
198 if(eta_det >= m_punchThroughEtaBins[i] && eta_det < m_punchThroughEtaBins[i+1]) etabin = i;
199 }
200 if(etabin<0) {
201 ATH_MSG_WARNING("There was a problem determining the eta bin to use for the punch through correction.");
202 //this could probably be improved, but to avoid a seg fault...
203 return 1;
204 }
205 double PunchThroughResponse = m_histTool_PunchThrough[etabin]->getValue(jet, jc);
206 if(PunchThroughResponse>1) return 1;
207 return PunchThroughResponse;
208}
209
211 if (jc.getValue<int>("nTrk")<=0) return 1; //nTrk < 0 is unphysical, nTrk = 0 is a special case, so return 1 for nTrk <= 0
212 if ( etabin >= m_histTool_nTrk.size() ) return 1.;
213 double nTrkResponse = m_histTool_nTrk[etabin]->getValue(jet, jc);
214 return nTrkResponse;
215}
216
218 if (jc.getValue<float>("trackWIDTH")<=0) return 1;
219 if ( etabin >= m_histTool_trackWIDTH.size() ) return 1.;
220 //jets with no tracks are assigned a trackWIDTH of -1, we use the trackWIDTH=0 correction in those cases
221 double trackWIDTHResponse = m_histTool_trackWIDTH[etabin]->getValue(jet, jc);
222 return trackWIDTHResponse;
223}
224
226 for ( const xAOD::Vertex* vertex : vertices ) {
227 if(vertex->vertexType() == xAOD::VxType::PriVtx) {
228 ATH_MSG_VERBOSE("GSCCalibStep " << name() << " Found HS vertex at index: "<< vertex->index());
229 return vertex;
230 }
231 }
232 if (vertices.size()==1) {
233 ATH_MSG_VERBOSE("GSCCalibStep " << name() << " Found no HS vertex, return dummy");
234 if (vertices.back()->vertexType() == xAOD::VxType::NoVtx)
235 return vertices.back();
236 }
237 ATH_MSG_VERBOSE("No vertex found in container.");
238 return nullptr;
239}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
unsigned int uint
const T * back() const
Access the last element in the collection as an rvalue.
size_type size() const noexcept
Returns the number of elements in the collection.
float getPunchThroughResponse(const xAOD::Jet &jet, const JetHelper::JetContext &jc, double eta_det) const
Gaudi::Property< std::vector< double > > m_punchThroughEtaBins
ToolHandleArray< JetHelper::IVarTool > m_histTool_ChargedFraction
virtual StatusCode calibrate(xAOD::JetContainer &) const override
Apply calibration to a jet container.
ToolHandleArray< JetHelper::IVarTool > m_histTool_trackWIDTH
float getNTrkResponse(const xAOD::Jet &jet, const JetHelper::JetContext &jc, uint etabin) const
Gaudi::Property< std::string > m_jetInScale
ToolHandleArray< JetHelper::IVarTool > m_histTool_EM3
float getTile0Response(const xAOD::Jet &jet, const JetHelper::JetContext &jc, uint etabin) const
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexContainer_key
Gaudi::Property< std::string > m_jetOutScale
const xAOD::Vertex * findHSVertex(const xAOD::VertexContainer &vertices) const
Retrieve hard scatter vertex for its index. Return nullptr if one cannot be found.
float getEM3Response(const xAOD::Jet &jet, const JetHelper::JetContext &jc, uint etabin) const
ToolHandleArray< JetHelper::IVarTool > m_histTool_Tile0
float getTrackWIDTHResponse(const xAOD::Jet &jet, const JetHelper::JetContext &jc, uint etabin) const
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
ToolHandleArray< JetHelper::IVarTool > m_histTool_PunchThrough
ToolHandleArray< JetHelper::IVarTool > m_histTool_nTrk
float getChargedFractionResponse(const xAOD::Jet &jet, const JetHelper::JetContext &jc, uint etabin) const
Functions for retrieving the correction factors.
Gaudi::Property< float > m_punchThroughMinPt
Gaudi::Property< bool > m_applyPunchThrough
Properties for the punch-through correction:
GSCCalibStep(const std::string &name="GSCCalibStep")
Constructor with parameters:
Class JetContext Designed to read AOD information related to the event, N vertices,...
Definition JetContext.h:24
bool setValue(const std::string &name, const T value, bool allowOverwrite=false)
Definition JetContext.h:52
void getValue(const std::string &name, T &value) const
Definition JetContext.h:39
size_t index() const
Return the index of this element within its container.
Helper class to provide constant type-safe access to aux data.
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
VxType::VertexType vertexType() const
The type of the vertex.
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
@ PriVtx
Primary vertex.
@ NoVtx
Dummy vertex. TrackParticle was not used in vertex fit.
Jet_v1 Jet
Definition of the current "jet version".
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
JetContainer_v1 JetContainer
Definition of the current "jet container version".
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition JetTypes.h:17