ATLAS Offline Software
Loading...
Searching...
No Matches
HIJetClusterSubtractorTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
13#include <TH3F.h>
14#include <TFile.h>
15#include <iomanip>
16
18 m_init(false),
19 m_h3W(nullptr),
20 m_h3Eta(nullptr),
21 m_h3Phi(nullptr)
22{
23 setUseCells(false);
24}
25
26void HIJetClusterSubtractorTool::subtract(xAOD::IParticle::FourMom_t& subtr_mom, const xAOD::IParticle* cl_in, const xAOD::HIEventShapeContainer* shape, const HIEventShapeIndex* index, const ToolHandle<IHIUEModulatorTool>& modulator, const xAOD::HIEventShape* eshape) const
27{
28
29 const xAOD::CaloCluster* cl=static_cast<const xAOD::CaloCluster*>(cl_in);
30 float E_unsubtr=cl->e(HIJetRec::unsubtractedClusterState());
31 float eta_unsubtr=cl->eta(HIJetRec::unsubtractedClusterState());
32 float phi_unsubtr=cl->phi(HIJetRec::unsubtractedClusterState());
33
34 float E_subtr=0;
35 float E_eta_subtr=0;
36 float E_phi_subtr=0;
37 const float eta0=cl->eta0();
38 const float phi0=cl->phi0();
39 float mod=1;
40 if(modulator) mod=modulator->getModulation(phi0, eshape);
41
42 float energy=E_unsubtr;
43 float eta=eta0;
44 float phi=phi0;
45
46 static const SG::ConstAccessor<float> HIEtaPhiWeightAcc("HIEtaPhiWeight");
47 if(HIEtaPhiWeightAcc.isAvailable(*cl))
48 {
49 float DF_weight=HIEtaPhiWeightAcc(*cl);
50 if(DF_weight!=0.) mod/=DF_weight;
51 }
52
53 if(shape==nullptr)
54 {
55 ATH_MSG_ERROR("No HIEventShape supplied, cannot do subtraction");
56 return;
57 }
58
60 {
61 if(index==nullptr)
62 {
63 ATH_MSG_ERROR("No HIEventShapeIndex supplied, cannot do subtraction");
64 return;
65 }
66
67 for(unsigned int sample=0; sample<24; sample++)
68 {
70 if(!cl->hasSampling(s)) continue;
71
72 float nCells=index->getShape(eta0,sample,shape)->nCells();
73 float rho=0;
74 if(nCells!=0.) rho=index->getShape(eta0,sample,shape)->rho()/nCells;
75 rho*=mod*std::cosh(eta0);
76
77 E_subtr+=rho*getWeight(eta0,phi0,sample);
78 E_eta_subtr+=rho*getWeightEta(eta0,phi0,sample);
79 E_phi_subtr+=rho*getWeightPhi(eta0,phi0,sample);
80 }
81
82 energy=E_unsubtr-E_subtr;
83 if(energy!=0.)
84 {
85 eta=eta_unsubtr+(eta_unsubtr*E_subtr-E_eta_subtr)/energy;
86 phi=phi_unsubtr+(phi_unsubtr*E_subtr-E_phi_subtr)/energy;
87 }
88 setSubtractedEtaPhi(energy,eta,phi,eta0,phi0,E_subtr==0 ? energy : energy/E_subtr);
89 }
90 else
91 {
92 unsigned int es_bin=HI::TowerBins::findBinEta(eta0);
93
94 float area=shape->at(es_bin)->area();
95 float rho=(area==0.) ? 0. : shape->at(es_bin)->et()/area;
96
97 rho*=mod*std::cosh(eta0);
98 E_subtr=rho*HI::TowerBins::getBinArea();
99 energy=E_unsubtr-E_subtr;
100 eta=eta_unsubtr;
101 phi=phi_unsubtr;
102
103 ATH_MSG_VERBOSE("Subtracting"
104 << std::setw(8) << eta0
105 << std::setw(8) << phi0
106 << std::setw(10) << std::setprecision(3) << modulator->getModulation(phi0, eshape)
107 << std::setw(10) << std::setprecision(3) << mod
108 << std::setw(10) << std::setprecision(3) << modulator->getModulation(phi0, eshape)/mod
109 << std::setw(10) << std::setprecision(3) << area
110 << std::setw(10) << std::setprecision(3) << rho*area
111 << std::setw(10) << std::setprecision(3) << E_unsubtr*1e-3
112 << std::setw(10) << std::setprecision(3) << E_subtr*1e-3
113 << std::setw(10) << std::setprecision(3) << E_unsubtr*1e-3/std::cosh(eta0)
114 << std::setw(10) << std::setprecision(3) << E_subtr*1e-3/std::cosh(eta0));
115
116
117
118
119 }
120 float ET=energy/std::cosh(eta);
121 subtr_mom.SetPxPyPzE(ET*std::cos(phi),ET*std::sin(phi),ET*std::sinh(eta),energy);
122}
123
125{
126
127 float eta0=cl->eta0();
128 float phi0=cl->phi0();
129
131 {
132 for(unsigned int sample=0; sample<24; sample++)
133 {
135 if(!cl->hasSampling(s)) continue;
136 float esamp=cl->eSample(s);
137 float ET=esamp/std::cosh(eta0);
138 xAOD::HIEventShape* slice=shape->at(index->getIndex(eta0,sample));
139 float area_cluster=getWeight(eta0,phi0,sample);
140 //update members
141 updateSlice(slice,ET,phi0,area_cluster);
142 }
143 }
144 else
145 {
146 unsigned int es_bin=HI::TowerBins::findBinEta(eta0);
147 xAOD::HIEventShape* slice=shape->at(es_bin);
148 constexpr float area_cluster=HI::TowerBins::getBinArea();
149 float ET=cl->e(HIJetRec::unsubtractedClusterState())/std::cosh(eta0);
150 static const SG::ConstAccessor<float> HIEtaPhiWeightAcc("HIEtaPhiWeight");
151 if(HIEtaPhiWeightAcc.isAvailable(*cl))
152 {
153 float HI_weight=HIEtaPhiWeightAcc(*cl);
154 ET*=HI_weight;
155 }
156 //update members
157 updateSlice(slice,ET,phi0,area_cluster);
158 }
159}
160
161
162void HIJetClusterSubtractorTool::subtractWithMoments(xAOD::CaloCluster* cl, const xAOD::HIEventShapeContainer* shape, const HIEventShapeIndex* index, const ToolHandle<IHIUEModulatorTool>& modulator, const xAOD::HIEventShape* eshape) const
163{
165 subtract(subtr_mom, cl, shape, index, modulator, eshape);
167}
168
170{
172 {
173 std::string local_path=static_cast<std::string>(m_configDir)+m_inputFile;
174 std::string full_path=PathResolverFindCalibFile(local_path);
175 ATH_MSG_INFO("Reading input file "<< m_inputFile << " from " << full_path);
176 TFile* f=TFile::Open(full_path.c_str());
177 if(f==nullptr)
178 {
179 ATH_MSG_FATAL("Cannot read config file " << full_path );
180 return StatusCode::FAILURE;
181 }
182
183 m_h3W=(TH3F*)f->GetObjectChecked("h3_w","TH3F");
184 if(m_h3W==nullptr)
185 {
186 ATH_MSG_FATAL("Cannot find TH3F m_h3W in config file " << full_path );
187 return StatusCode::FAILURE;
188 }
189
190 m_h3Eta=(TH3F*)f->GetObjectChecked("h3_eta","TH3F");
191 if(m_h3Eta==nullptr)
192 {
193 ATH_MSG_FATAL("Cannot find TH3F m_h3Eta in config file " << full_path );
194 return StatusCode::FAILURE;
195 }
196
197 m_h3Phi=(TH3F*)f->GetObjectChecked("h3_phi","TH3F");
198 if(m_h3Phi==nullptr)
199 {
200 ATH_MSG_FATAL("Cannot find TH3F m_h3Phi in config file " << full_path );
201 return StatusCode::FAILURE;
202 }
203 m_h3W->SetDirectory(0);
204 m_h3Eta->SetDirectory(0);
205 m_h3Phi->SetDirectory(0);
206
207 f->Close();
208 }
209 m_init=true;
210 return StatusCode::SUCCESS;
211
212}
213
214
216{
217 return initializeTool();
218
219}
220
221
222float HIJetClusterSubtractorTool::getWeight(float eta, float phi, int sample) const
223{
224 return m_h3W->GetBinContent(m_h3W->FindFixBin(eta,phi,sample));
225}
226float HIJetClusterSubtractorTool::getWeightEta(float eta, float phi, int sample) const
227{
228 return m_h3Eta->GetBinContent(m_h3Eta->FindFixBin(eta,phi,sample));
229}
230float HIJetClusterSubtractorTool::getWeightPhi(float eta, float phi, int sample) const
231{
232 return m_h3Phi->GetBinContent(m_h3Phi->FindFixBin(eta,phi,sample));
233}
234
235void HIJetClusterSubtractorTool::updateSlice(xAOD::HIEventShape* slice, float ET, float phi0, float area_cluster) const
236{
237 if(area_cluster==0.)
238 {
239 area_cluster=1.;
240
241 if(ET==0.)
242 {
243 ATH_MSG_INFO("Provided cluster area was 0, check if the input file " << m_inputFile << " from " << PathResolverFindCalibFile(static_cast<std::string>(m_configDir)+m_inputFile) << " was generated with the currently used geometry. See ATLHI-472");
244 }
245 else
246 {
247 ATH_MSG_ERROR("Provided cluster area was 0, the HIEventShape will be incorrect! Check if the input file " << m_inputFile << " from " << PathResolverFindCalibFile(static_cast<std::string>(m_configDir)+m_inputFile) << " was generated with the currently used geometry. See ATLHI-472");
248 }
249 }
250
251 float area_slice=slice->area();
252 float area_sf=1.;
253
254 if(area_slice!=0.) area_sf-=area_cluster/area_slice;
255 slice->setNCells( std::floor(area_sf*slice->nCells()) );
256 slice->setEt(slice->et()-ET);
257 slice->setArea(area_slice-area_cluster);
258
259 slice->setRho(slice->rho() - ET/area_cluster);
260
261
262 for(unsigned int ih=0; ih<slice->etCos().size(); ih++)
263 {
264 float ih_f=ih+1;
265 float tmp_cos = slice->etCos().at(ih);
266 slice->etCos()[ih] = tmp_cos + cos(ih_f*phi0)*ET;
267
268 float tmp_sin = slice->etSin().at(ih);
269 slice->etSin()[ih] = tmp_sin + sin(ih_f*phi0)*ET;
270 }
271}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
Helper class to provide constant type-safe access to aux data.
double area(double R)
std::string PathResolverFindCalibFile(const std::string &logical_file_name)
const T * at(size_type n) const
Access an element, as an rvalue.
virtual void subtractWithMoments(xAOD::CaloCluster *, const xAOD::HIEventShapeContainer *, const HIEventShapeIndex *index, const ToolHandle< IHIUEModulatorTool > &, const xAOD::HIEventShape *eshape) const override
Gaudi::Property< std::string > m_configDir
HIJetClusterSubtractorTool(const std::string &myname)
Implements method defined in base First argument is reference to four vector that is updated to refle...
void updateSlice(xAOD::HIEventShape *slice, float ET, float phi0, float area_cluster) const
Gaudi::Property< std::string > m_inputFile
virtual void updateUsingCluster(xAOD::HIEventShapeContainer *shape, const HIEventShapeIndex *index, const xAOD::CaloCluster *cl) const override
Method to update the shape based on a given cluster two sets of indices are passed by reference and u...
float getWeightPhi(float eta, float phi, int sample) const
virtual void subtract(xAOD::IParticle::FourMom_t &, const xAOD::IParticle *, const xAOD::HIEventShapeContainer *, const HIEventShapeIndex *, const ToolHandle< IHIUEModulatorTool > &, const xAOD::HIEventShape *eshape) const override
Abstract method where particle itself is not modified IParticle::FourMom_t containing kinematics afte...
float getWeight(float eta, float phi, int sample) const
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Gaudi::Property< bool > m_useSamplings
float getWeightEta(float eta, float phi, int sample) const
void setSubtractedEtaPhi(float E, float &eta, float &phi, float eta0, float phi0, float sig) const
HIJetSubtractorToolBase(const std::string &myname)
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.
Class providing the definition of the 4-vector interface.
TLorentzVector FourMom_t
Definition of the 4-momentum type.
constexpr xAOD::CaloCluster::State unsubtractedClusterState()
constexpr xAOD::CaloCluster::State subtractedClusterState()
void setClusterP4(const xAOD::CaloCluster::FourMom_t &p, xAOD::CaloCluster *cl, xAOD::CaloCluster::State s)
constexpr float getBinArea()
Definition HIEventDefs.h:34
unsigned int findBinEta(float eta)
Definition HIEventDefs.h:46
Definition index.py:1
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
HIEventShapeContainer_v2 HIEventShapeContainer
Define the latest version of the container.
HIEventShape_v2 HIEventShape
Definition of the latest event info version.