ATLAS Offline Software
Loading...
Searching...
No Matches
CalibrationNtupleMakerTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5
7// CalibrationNtupleMakerTool.cxx, (c) ATLAS Detector software
9
11
14
17
19
20// Gaudi
21#include "GaudiKernel/ITHistSvc.h"
22#include "GaudiKernel/SystemOfUnits.h"
23
24// Root
25#include "TTree.h"
26#include "TString.h"
27#include "TH1F.h"
28
29#include <vector>
30
31using Gaudi::Units::GeV;
32
33CalibrationNtupleMakerTool::CalibrationNtupleMakerTool(const std::string& name, ISvcLocator* pSvcLocator) :
34 AthAlgorithm(name, pSvcLocator),
35 m_treeFolder("/calibration/"),
36 m_treeDescription("Calibration Ntuple"),
37 m_truthJetContainerName("MyAntiKt10TruthJets"),
38 m_vertexContainerName("PrimaryVertices"),
39 m_recoIsoDR(1.5),
40 m_recoIsoPtCut(100.*GeV),
41 m_trueIsoDR(2.5),
42 m_trueIsoPtCut(100.*GeV),
43 m_matchingCut(0.6),
44 m_h_events(nullptr),
45 m_index(nullptr),
46 m_etaCalo(nullptr),
47 m_etaDetCalo(nullptr),
48 m_phiCalo(nullptr),
49 m_eCalo(nullptr),
50 m_mCalo(nullptr),
51 m_etaCorr(nullptr),
52 m_etaDetCorr(nullptr),
53 m_phiCorr(nullptr),
54 m_eCorr(nullptr),
55 m_mCorr(nullptr),
56 m_etaTrue(nullptr),
57 m_phiTrue(nullptr),
58 m_eTrue(nullptr),
59 m_mTrue(nullptr)
60 {
61 declareProperty("FolderName" , m_treeFolder);
62 declareProperty("Description", m_treeDescription);
63
64 // the jets collections to calibrate
65 declareProperty("JetCollections" , m_collectionNames);
66 declareProperty("TruthJetContainerName", m_truthJetContainerName);
67 declareProperty("VertexContainerName" , m_vertexContainerName);
68
69 declareProperty("RecoIsolationDR" , m_recoIsoDR);
70 declareProperty("RecoIsolationPtCut" , m_recoIsoPtCut);
71 declareProperty("TruthIsolationDR" , m_trueIsoDR);
72 declareProperty("TruthIsolationPtCut", m_trueIsoPtCut);
73 declareProperty("MatchingCut" , m_matchingCut);
74 }
75
77= default;
78
80{
81 ATH_MSG_INFO( "initialize()" );
82
83 if (bookTree().isFailure()){
84 ATH_MSG_FATAL( "Could not book the TTree object" );
85 return StatusCode::FAILURE;
86 }
87
88 ATH_CHECK( m_evt.initialize() );
89
90 return StatusCode::SUCCESS;
91}
92
94{
95 ATH_MSG_INFO( "bookTree()" );
96
97
98 for (auto& name : m_collectionNames) {
99
100 // creating the tree for thw jet collection
101 TTree * tree = new TTree (name.c_str(), m_treeDescription.c_str());
102
103 // add the branches
104 tree->Branch("EventWeight" , &m_eventWeight );
105
106 tree->Branch("eta_calo" , &m_etaCalo );
107 tree->Branch("eta_det_calo" , &m_etaDetCalo );
108 tree->Branch("phi_calo" , &m_phiCalo );
109 tree->Branch("E_calo" , &m_eCalo );
110 tree->Branch("m_calo" , &m_mCalo );
111
112 tree->Branch("eta_corr1" , &m_etaCorr );
113 tree->Branch("eta_det_corr1" , &m_etaDetCorr );
114 tree->Branch("phi_corr1" , &m_phiCorr );
115 tree->Branch("E_corr1" , &m_eCorr );
116 tree->Branch("m_corr1" , &m_mCorr );
117
118 tree->Branch("eta_true" , &m_etaTrue );
119 tree->Branch("phi_true" , &m_phiTrue );
120 tree->Branch("E_true" , &m_eTrue );
121 tree->Branch("m_true" , &m_mTrue );
122
123 tree->Branch("index" , &m_index );
124 tree->Branch("mu" , &m_mu );
125 tree->Branch("NPV" , &m_npv );
126
127
128 m_trees.insert(std::pair<std::string, TTree*>(name, tree));
129
130 }
131
132 // now register the Tree
133 SmartIF<ITHistSvc> tHistSvc{service("THistSvc")};
134 ATH_CHECK( tHistSvc.isValid() );
135
136 if (tHistSvc) {
137 for (const auto& name : m_collectionNames) {
138 if((tHistSvc->regTree(m_treeFolder+name, m_trees.at(name))).isFailure()) {
139 ATH_MSG_ERROR( "initialize() Could not register the validation Tree!" );
140 return StatusCode::FAILURE;
141 }
142 }
143 }
144
145 // now register the Event histogram
146 m_h_events = new TH1F("h_events","total events", 10, 0, 10);
147
148 if (tHistSvc and (tHistSvc->regHist(m_treeFolder+m_h_events->GetName(), m_h_events)).isFailure()) {
149 ATH_MSG_ERROR( "Can not register histogram" << m_h_events->GetName() );
150 return StatusCode::FAILURE;
151 }
152
153 ATH_MSG_INFO("Calibration Tree booked and registered successfully!");
154
155 return StatusCode::SUCCESS;
156
157}
158
160{
161 m_h_events->Fill(0);
162
164 if(!evt.isValid()) {
165 ATH_MSG_FATAL( "Unable to retrieve Event Info" );
166 }
167 float ev_weight = evt->mcEventWeight();
168
170 if (not truths) return StatusCode::FAILURE;
171
173 if (not vertices) return StatusCode::FAILURE;
174
175 // get mu
176 float mu= evt->averageInteractionsPerCrossing();
177
178 //get NPV
179 float npv = 0.;
180
181 for (const auto vertex : *vertices) {
182 if (vertex->nTrackParticles()>=2)
183 npv++;
184 }
185
186 for (auto& name : m_collectionNames) {
187 const auto *const jets = getContainer<xAOD::JetContainer>(name);
188
189 m_etaCalo ->clear();
190 m_etaDetCalo ->clear();
191 m_phiCalo ->clear();
192 m_eCalo ->clear();
193 m_mCalo ->clear();
194 m_etaCorr ->clear();
195 m_etaDetCorr ->clear();
196 m_phiCorr ->clear();
197 m_eCorr ->clear();
198 m_mCorr ->clear();
199 m_etaTrue ->clear();
200 m_phiTrue ->clear();
201 m_eTrue ->clear();
202 m_mTrue ->clear();
203 m_index ->clear();
204
205 for (const auto truth: *truths) {
206
207 // here we match to the reco
208 int index = 0;
209 std::vector<const xAOD::Jet*> matched = {};
210
211 int Nmatches = Matched(truth, jets, matched, index);
212
213 // skip truth jets that don't match any reco jets
214 if (Nmatches==0) continue;
215
216 // skip the jets that are not isolated
217 if ( m_recoIsoDR > 0 ) {
218 double DRminReco = DRmin(matched.at(0),jets,m_recoIsoPtCut);
219 if ( DRminReco < m_recoIsoDR ) continue;
220 }
221
222 if ( m_trueIsoDR > 0 ) {
223 double DRminTruth = DRmin(truth,truths,m_trueIsoPtCut);
224 if ( DRminTruth < m_trueIsoDR ) continue;
225 }
226
227 //Storing variables
228 m_etaTrue->push_back(truth->eta());
229 m_phiTrue->push_back(truth->phi());
230 m_mTrue->push_back(truth->m()/GeV);
231 m_eTrue->push_back(truth->e()/GeV);
232
233 m_etaCalo->push_back(jets->at(index)->eta());
234 m_phiCalo->push_back(jets->at(index)->phi());
235 m_mCalo->push_back(jets->at(index)->m()/GeV);
236 m_eCalo->push_back(jets->at(index)->e()/GeV);
237
238 float detectorEta = DetectorEta(jets->at(index));
239 m_etaDetCalo->push_back(detectorEta);
240
241 m_etaCorr->push_back(jets->at(index)->eta());
242 m_phiCorr->push_back(jets->at(index)->phi());
243 m_mCorr->push_back(jets->at(index)->m()/GeV);
244 m_eCorr->push_back(jets->at(index)->e()/GeV);
245 m_etaDetCorr->push_back(detectorEta);
246
247 m_index->push_back(index);
248
249 }
250
251 m_eventWeight = ev_weight;
252 m_mu = mu;
253 m_npv = npv;
254
255 m_trees.at(name)->Fill();
256
257 }
258
259 return StatusCode::SUCCESS;
260
261}
262
263int CalibrationNtupleMakerTool::Matched(const xAOD::Jet* truth, const xAOD::JetContainer* jets, std::vector<const xAOD::Jet*>& matched, int& index) const {
264
265 int Nmatches = 0;
266 double drmin = 999.;
267 int Min_index=-1;
268
269 for (unsigned int ind = 0; ind < jets->size(); ind++) {
270 double dr = truth->p4().DeltaR(jets->at(ind)->p4());
271 if (dr < m_matchingCut) ++Nmatches;
272 //find minimum:
273 if (dr < drmin) {
274 drmin = dr;
275 Min_index = ind;
276 }
277 }
278
279 if (drmin<m_matchingCut) {
280 matched.push_back(jets->at(Min_index));
281 index = Min_index;
282 }
283
284 return Nmatches;
285}
286
287double CalibrationNtupleMakerTool::DRmin(const xAOD::Jet* myjet, const xAOD::JetContainer* jets, double PtMin) {
288
289 double DRmin=9999;
290 for (const auto jet : *jets) {
291 if (PtMin>0. and jet->pt()<PtMin) continue;
292 double Dr = myjet->p4().DeltaR(jet->p4());
293 if (Dr>0.0001 and Dr<DRmin)
294 DRmin=Dr;
295 }
296 return DRmin;
297}
298
300
301 xAOD::IParticle::FourMom_t corrP4(0,0,0,0);
302
303 const auto& partLinks = jet->constituentLinks();
304 for (const xAOD::IParticleLink& link : partLinks) {
305
306 if (not link.isValid()) {
307 ATH_MSG_WARNING("Got an invalid element link. Returning jet eta...");
308 return jet->eta();
309 }
310
311 const xAOD::TrackCaloCluster* tcc = dynamic_cast<const xAOD::TrackCaloCluster*>(*link);
312
313 static const SG::AuxElement::Accessor< float > acc_detEta( "DetectorEta" );
314 float det_eta = tcc->eta();
315
316 if (acc_detEta.isAvailable(*tcc)) {
317 det_eta = acc_detEta(*tcc);
318 } else
319 ATH_MSG_WARNING("DetectorEta decoration not found for TCCs! Using eta...");
320
321 double pt = tcc->p4().P()/cosh(det_eta);
322
324 p4CorrCl.SetPtEtaPhiE(pt, det_eta, tcc->p4().Phi(), tcc->p4().E());
325 if(tcc->p4().E() < 0.) p4CorrCl*=-1.;
326 corrP4 += p4CorrCl;
327 }
328
329 return corrP4.Eta();
330}
331
333{
334 ATH_MSG_INFO( "finalize()" );
335
336 delete m_index;
337
338 delete m_etaCalo;
339 delete m_etaDetCalo;
340 delete m_phiCalo;
341 delete m_eCalo;
342 delete m_mCalo;
343
344 delete m_etaCorr;
345 delete m_etaDetCorr;
346 delete m_phiCorr;
347 delete m_eCorr;
348 delete m_mCorr;
349
350 delete m_etaTrue;
351 delete m_phiTrue;
352 delete m_eTrue;
353 delete m_mTrue;
354
355 return StatusCode::SUCCESS;
356
357}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
const T * getContainer(const std::string &containerName)
static double DRmin(const xAOD::Jet *myjet, const xAOD::JetContainer *jets, double PtMin)
CalibrationNtupleMakerTool(const std::string &name, ISvcLocator *pSvcLocator)
Standard Athena-Algorithm Constructor.
std::vector< std::string > m_collectionNames
virtual ~CalibrationNtupleMakerTool()
Default Destructor.
float DetectorEta(const xAOD::Jet *jet)
int Matched(const xAOD::Jet *truth, const xAOD::JetContainer *jets, std::vector< const xAOD::Jet * > &matched, int &index) const
StatusCode finalize()
standard Athena-Algorithm method
StatusCode execute()
standard Athena-Algorithm method
StatusCode initialize()
standard Athena-Algorithm method
SG::ReadHandleKey< xAOD::EventInfo > m_evt
std::map< std::string, TTree * > m_trees
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
TLorentzVector FourMom_t
Definition of the 4-momentum type.
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition Jet_v1.cxx:71
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition index.py:1
Jet_v1 Jet
Definition of the current "jet version".
IParticleLink_v1 IParticleLink
Define the latest version of the IParticleLink class.
TrackCaloCluster_v1 TrackCaloCluster
Reference the current persistent version:
JetContainer_v1 JetContainer
Definition of the current "jet container version".
TChain * tree