ATLAS Offline Software
Loading...
Searching...
No Matches
FFJetSmearingTool_MyExample.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6// //
7// Name: FFJetSmearingTool_MyExample //
8// Purpose: example code for using the FFJetSmearingTool //
9// ---> For more information about the tool, check TWiki: https://twiki.cern.ch/twiki/bin/viewauth/AtlasProtected/FFJetSmearingTool //
10// //
11// * Authors: Alberto Prades and Davide Melini //
12// * Updates: Daniel Camarero Munoz //
13// * July 2024: Migration to R22+ for Phase-I large-R UFO pre-recommendations //
14// * November 2024: Available in Athena by default //
15// //
17
18// System
19#include <memory>
20#include <AsgTools/ToolHandle.h>
21#include <AsgTools/AsgTool.h>
23
24// xAOD EDM classes
27
34
35// Jet tools
39
40// ROOT
41#include <TFile.h>
42#include <TCanvas.h>
43#include <TH1F.h>
44#include <TTree.h>
45
46// Error checking macro
47#define CHECK( ARG ) \
48 do { \
49 const bool result = ARG; \
50 if (!result){ \
51 std::cout << "Failed to execute: " << ARG << std::endl; \
52 return 1; \
53 } \
54 } while( false )
55
56// Help message if the --help option is given by the user
57void usage(){
58 std::cout << "Running options:" << std::endl;
59 std::cout << "YOU HAVE TO ADAPT THE OPTIONS TO FFJETSMEARINGCORRECTION" << std::endl;
60 std::cout << " --help : To get the help you're reading" << std::endl;
61 std::cout << " --jetColl= : Specify the jet collection (e.g. UFO+CSSK jets)" << std::endl;
62 std::cout << " --MassDef= : Specify the kind of jet mass used (e.g. UFO)" << std::endl;
63 std::cout << " --MCType= : Specify the MC campaign (e.g. MC20, MC21, MC23, MC20AF3, MC23AF3)" << std::endl;
64 std::cout << " --sample= : Specify input xAOD" << std::endl;
65 std::cout << " Example: FFJetSmearingTool_MyExample --truth_jetColl=AntiKt10TruthSoftDropBeta100Zcut10Jets --reco_jetColl=AntiKt10UFOCSSKSoftDropBeta100Zcut10Jets --MassDef=UFO --MCType=MC23 --eventsMax=30000 --sample=/eos/atlas/atlascerngroupdisk/perf-jets/INSITU/TestFiles/dcamarer/DAOD_PHYS_MC_Zqqjets/mc20_13TeV/DAOD_PHYS.40038344._000011.pool.root.1 --output=.root --ConfigFile=/afs/cern.ch/user/d/dcamarer/private/PostDoc/JETM/JMRcombination_PhaseILargeR/test-jetcalibtools-configs/R22/source/JetUncertainties/share/DCM_240502_new/R10_FullJMR_Phase1.config --DebugTool=false" << std::endl;
66}
67
69// Main Function //
71
72int main(int argc, char* argv[]){
73
75
77 std::string sample = "";
78 std::string truth_jetColl = "";
79 std::string reco_jetColl = "";
80 std::string string_kindofmass = "";
81 std::string string_kindofmc = "";
82 std::string output = "";
83 std::string string_configfile_path = "";
84 std::string string_debugtool = "";
85 int eventsMax = 30000;
86
88 for (int i=1; i< argc; i++){
89
90 std::string opt(argv[i]); std::vector< std::string > v;
91 std::istringstream iss(opt);
92 std::string item;
93 char delim = '=';
94
95 while (std::getline(iss, item, delim)){
96 v.push_back(item);
97 }
98
99 if ( opt.find("--help") != std::string::npos ){
100 usage();
101 return 0;
102 }
103
104 if ( opt.find("--sample=") != std::string::npos ) sample = v[1];
105 if ( opt.find("--truth_jetColl=") != std::string::npos ) truth_jetColl = v[1];
106 if ( opt.find("--reco_jetColl=") != std::string::npos ) reco_jetColl = v[1];
107 if ( opt.find("--MassDef=") != std::string::npos ) string_kindofmass = v[1];
108 if ( opt.find("--MCType=") != std::string::npos ) string_kindofmc = v[1];
109 if ( opt.find("--ConfigFile=") != std::string::npos ) string_configfile_path = v[1];
110 if ( opt.find("--DebugTool=") != std::string::npos ) string_debugtool = v[1];
111 if ( opt.find("--output=") != std::string::npos ) output = v[1];
112 if ( opt.find("--eventsMax=") != std::string::npos ) eventsMax = std::atoi(v[1].data());
113
114 }
115
116 if (sample==""){
117 std::cout << "No input xAOD file specified, exiting" << std::endl;
118 return 1;
119 }
120 if (truth_jetColl==""){
121 std::cout << "No truth jet collection specified, exiting" << std::endl;
122 return 1;
123 }
124 if (reco_jetColl==""){
125 std::cout << "No truth jet collection specified, exiting" << std::endl;
126 return 1;
127 }
128 if (string_kindofmass==""){
129 std::cout << "No kind of jet mass specified, exiting" << std::endl;
130 return 1;
131 }
132 if (string_kindofmc==""){
133 std::cout << "No kind of MC specified, exiting" << std::endl;
134 return 1;
135 }
136 if (string_configfile_path==""){
137 std::cout << "No ConfigFile specified, exiting" << std::endl;
138 return 1;
139 }
140 if (string_debugtool==""){
141 std::cout << "No debugtool specified, exiting" << std::endl;
142 return 1;
143 }
144 if (output==""){
145 std::cout << "Output not specified, exiting" << std::endl;
146 return 1;
147 }
148 std::cout << sample << truth_jetColl << reco_jetColl << output << string_debugtool << std::endl;
149
150 TString kindofmass = string_kindofmass;
151 TString kindofmc = string_kindofmc;
152
153 bool want_to_debug = false;
154 if (string_debugtool == "true"){
155 want_to_debug = true;
156 }else if (string_debugtool == "false"){
157 want_to_debug = false;
158 }
159
161
162 std::unique_ptr< TFile > ifile( TFile::Open( sample.c_str(), "READ" ) );
163
164 // Create a Event object
165 auto event = xAOD::Event::createAndReadFrom(*ifile);
166 if (!event) {
167 std::cout << "Failed to open file " << std::endl;
168 return 1;
169 }
170
172
173 const std::string name_FFJetSmearingTool = "FFJetSmearing_Example";
174 CP::FFJetSmearingTool ffjetsmearingtool(name_FFJetSmearingTool.c_str());
175 CHECK(ffjetsmearingtool.setProperty("MassDef", kindofmass));
176 CHECK(ffjetsmearingtool.setProperty("MCType", kindofmc));
177 CHECK(ffjetsmearingtool.setProperty("ConfigFile", string_configfile_path));
178 CHECK(ffjetsmearingtool.setProperty("OutputLevel", MSG::ERROR)); // Set output level threshold (1=VERBOSE, 2=DEBUG, 3=INFO, 4=WARNING, 5=ERROR, 6=FATAL, 7=ALWAYS)
179 if (want_to_debug){
180 CHECK(ffjetsmearingtool.setProperty("OutputLevel", MSG::VERBOSE));
181 }
182 CHECK(ffjetsmearingtool.initialize());
183
185
186 const CP::SystematicSet& recommendedSysts = ffjetsmearingtool.recommendedSystematics();
187 std::cout << "\nRecommended systematics:\n" << std::endl;
188 for (auto sysItr = recommendedSysts.begin(); sysItr != recommendedSysts.end(); ++sysItr){
189 std::cout << sysItr->name().c_str() << std::endl;
190 }
191
192 std::vector<CP::SystematicSet> sysList;
193
195 // Note: in principle not needed for most of the derivation formats, as we save this information in DAODs
196 /*
197 JetTruthLabelingTool m_JetTruthLabelingTool("JetTruthLabelingTool");
198 CHECK(m_JetTruthLabelingTool.setProperty("TruthLabelName", "R10TruthLabel_R22v1"));
199 CHECK(m_JetTruthLabelingTool.setProperty("UseTRUTH3", false)); // Set this to false only if you have the FULL !TruthParticles container in your input file
200 CHECK(m_JetTruthLabelingTool.setProperty("TruthParticleContainerName", "TruthParticles")); // Set this if you have the FULL !TruthParticles container but have named it something else
201 CHECK(m_JetTruthLabelingTool.initialize());
202 */
203
205
206 const std::string name_JetCalibTools = "JetCalibrationTool";
207 JetCalibrationTool m_JetCalibrationTool_handle(name_JetCalibTools.c_str());
208
209 CHECK(m_JetCalibrationTool_handle.setProperty("JetCollection", "AntiKt10UFOCSSKSoftDropBeta100Zcut10"));
210 CHECK(m_JetCalibrationTool_handle.setProperty("ConfigFile", "JES_MC20PreRecommendation_R10_UFO_CSSK_SoftDrop_JMS_R21Insitu_26Nov2024.config"));
211 //CHECK(m_JetCalibrationTool_handle.setProperty("ConfigFile", "share/JES_MC20PreRecommendation_R10_UFO_CSSK_SoftDrop_JMS_R21Insitu_02Aug2024.config")); // When using DEVmode
212 CHECK(m_JetCalibrationTool_handle.setProperty("CalibSequence", "EtaJES_JMS"));
213 CHECK(m_JetCalibrationTool_handle.setProperty("IsData", false));
214 CHECK(m_JetCalibrationTool_handle.setProperty("CalibArea", "00-04-82"));
215 CHECK(m_JetCalibrationTool_handle.setProperty("DEVmode", false));
216 //CHECK(m_JetCalibrationTool_handle.setProperty("DEVmode", true)); // When using DEVmode
217 CHECK(m_JetCalibrationTool_handle.initialize());
218
220
222 /*
223 sysList.insert( sysList.begin(), CP::SystematicSet() );
224 const CP::SystematicVariation nullVar = CP::SystematicVariation("Nominal");
225 sysList.back().insert(nullVar);
226 */
227
228 for (auto sys : recommendedSysts){
229 sysList.push_back(CP::SystematicSet({sys}));
230 }
231
232 std::cout << "\n=============SYSTEMATICS CHECK NOW" << std::endl;
233 for (auto sys : sysList){
234
235 std::cout << "\nRunning over the systematic " << sys.name().c_str() << std::endl;
236 static constexpr float MeVtoGeV = 1.e-3;
237
238 // Tell the calibration tool which variation to apply
239 if (ffjetsmearingtool.applySystematicVariation(sys) != StatusCode::SUCCESS){
240 std::cout << "Error, Cannot configure calibration tool for systematics" << std::endl;
241 }
242
243 // Define histograms
244 TH1F* reco_jet_mass_hist = new TH1F("reco_jet_mass_hist","reco_jet_mass_hist", 300, 0, 600);
245 TH1F* matched_truth_jet_mass_hist = new TH1F("matched_truth_jet_mass_hist","matched_truth_jet_mass_hist", 300, 0, 600);
246 TH1F* smeared_reco_jet_mass_hist = new TH1F("smeared_reco_jet_mass_hist","smeared_reco_jet_mass_hist", 300, 0, 600);
247 //
248 TH1F* reco_jet_pt_hist = new TH1F("reco_jet_pt_hist","reco_jet_pt_hist", 800, 0, 4000);
249 TH1F* matched_truth_jet_pt_hist = new TH1F("matched_truth_jet_pt_hist","matched_truth_jet_pt_hist", 800, 0, 4000);
250 TH1F* smeared_reco_jet_pt_hist = new TH1F("smeared_reco_jet_pt_hist","smeared_reco_jet_pt_hist", 800, 0, 4000);
251 //
252 TH1F* reco_jet_rapidity_hist = new TH1F("reco_jet_rapidity_hist","reco_jet_rapidity_hist", 50, -2.5, 2.5);
253 TH1F* matched_truth_jet_rapidity_hist = new TH1F("matched_truth_jet_rapidity_hist","matched_truth_jet_rapidity_hist", 50, -2.5, 2.5);
254 TH1F* smeared_reco_jet_rapidity_hist = new TH1F("smeared_reco_jet_rapidity_hist","smeared_reco_jet_rapidity_hist", 50, -2.5, 2.5);
255
256 reco_jet_mass_hist->Sumw2();
257 matched_truth_jet_mass_hist->Sumw2();
258 smeared_reco_jet_mass_hist->Sumw2();
259
260 reco_jet_pt_hist->Sumw2();
261 matched_truth_jet_pt_hist->Sumw2();
262 smeared_reco_jet_pt_hist->Sumw2();
263
264 reco_jet_rapidity_hist->Sumw2();
265 matched_truth_jet_rapidity_hist->Sumw2();
266 smeared_reco_jet_rapidity_hist->Sumw2();
267
268 int upperlimit1 = 600;
269 int upperlimit2 = 1000;
270
271 int numBinsMass = 120;
272 int numBinsPt = 100;
273
274 std::unique_ptr<TH3F> hist_jet_mass_scale_change_3D;
275 hist_jet_mass_scale_change_3D = std::make_unique<TH3F>("hist_jet_mass_scale_change_3D","hist_jet_mass_scale_change_3D",numBinsPt,0,upperlimit2,numBinsMass,0,upperlimit1,numBinsMass,0,upperlimit1);
276
277 float lowerlimit3 = -0.5;
278 float upperlimit3 = 0.5;
279 int numBinsDiff = 100;
280
281 std::unique_ptr<TH3F> hist_jet_mass_resolution_change_3D;
282 hist_jet_mass_resolution_change_3D = std::make_unique<TH3F>("hist_jet_mass_resolution_change_3D","hist_jet_mass_resolution_change_3D",numBinsPt,0,upperlimit2,numBinsMass,0,upperlimit1,numBinsDiff,lowerlimit3,upperlimit3);
283
285
286 Long64_t nevents = event->getEntries();
287
288 if (eventsMax < nevents){
289 nevents = eventsMax;
290 }
291
292 for (Long64_t ievent = 0; ievent < nevents; ++ievent){
293
294 // Load the event:
295 if ( event->getEntry( ievent ) < 0 ){
296 std::cout << "Failed to load entry " << ievent << std::endl;
297 return 1;
298 }
299
300 // Show status
301 if (ievent % 1000==0){
302 std::cout << "Event " << ievent << " of " << nevents << std::endl;
303 }
304
305 // Print some event information for fun
306 if (want_to_debug){
307 const xAOD::EventInfo* ei = nullptr;
308 CHECK( event->retrieve(ei, "EventInfo") );
309
310 std::cout << "===>>> start processing event " << ei->eventNumber() << ", run " << ei->runNumber() << " - Events processed so far: " << ievent << std::endl;
311
312 // Get the truth jets from the event
313 const xAOD::JetContainer* jets_truth = nullptr;
314 CHECK( event->retrieve(jets_truth, truth_jetColl) );
315 std::cout << "Number of truth jets: " << jets_truth->size() << std::endl;
316
317 // Loop over the truth jets in the event
318 for (const xAOD::Jet* jet_truth : *jets_truth){
319 // Print basic info about this jet
320 std::cout << "Truth Jet: pt = " << jet_truth->pt()*MeVtoGeV << ", mass = " << jet_truth->m()*MeVtoGeV << ", eta = " << jet_truth->eta() << std::endl;
321 }
322
323 // Get the reco jets from the event
324 const xAOD::JetContainer* jets_reco = nullptr;
325 CHECK( event->retrieve(jets_reco, reco_jetColl) );
326 std::cout << "Number of reco jets: " << jets_reco->size() << std::endl;
327
328 //Loop over the reco jets in the event
329 for (const xAOD::Jet* jet_reco : *jets_reco){
330 // Print basic info about this jet
331 std::cout << "Reco Jet: pt = " << jet_reco->pt()*MeVtoGeV << ", mass = " << jet_reco->m()*MeVtoGeV << ", eta = " << jet_reco->eta() << std::endl;
332 }
333 }
334
335 xAOD::Jet jet_truth_matched;
336 jet_truth_matched.makePrivateStore();
337
338 // Retrieve jet container
339 const xAOD::JetContainer* jets = nullptr;
340 CHECK( event->retrieve( jets, reco_jetColl ) );
341
342 // Shallow copy
343 auto jets_shallowCopy = xAOD::shallowCopy( *jets );
344
346 // Note: in principle not needed for most of the derivaiton formats, as we save this information in DAODs
347 //m_JetTruthLabelingTool.modify(*(jets_shallowCopy.first));
348
349 if (want_to_debug){
350 std::cout << "Start the loop over the jets" << std::endl;
351 }
352
353 bool lead_jet = true;
354
355 // Apply the jet calibration
356 if ( m_JetCalibrationTool_handle.applyCalibration( *(jets_shallowCopy.first) ) != StatusCode::SUCCESS ){
357 std::cout << "JetCalibration tool reported a EL::StatusCode::FAILURE" << std::endl;
358 return 1;
359 }
360
362
363 for ( xAOD::Jet* jet_reco : *(jets_shallowCopy.first) ){
364
365 double jetpt = jet_reco->pt() * MeVtoGeV;
366 double jetm = jet_reco->m() * MeVtoGeV;
367 double jetrapidity = jet_reco->rapidity();
368 // Jet pT cut
369 if (jetpt < 200 or jetpt > 3000) continue;
370 // Jet mass cut
371 if (jetm > 600) continue;
372 // Jet rapidity cut
373 if (abs(jetrapidity) > 2) continue;
374
375 if ( ffjetsmearingtool.getMatchedTruthJet(*jet_reco, jet_truth_matched) != StatusCode::SUCCESS ){
376 continue;
377 }
378
379 double aux_original_jet_mass = jet_reco->m() * MeVtoGeV;
380
381 if (lead_jet == true && aux_original_jet_mass > 0){
382
383 reco_jet_mass_hist->Fill(jet_reco->m() * MeVtoGeV);
384 reco_jet_pt_hist->Fill(jet_reco->pt() * MeVtoGeV);
385 reco_jet_rapidity_hist->Fill(jet_reco->rapidity());
386 //
387 matched_truth_jet_mass_hist->Fill(jet_truth_matched.m() * MeVtoGeV);
388 matched_truth_jet_pt_hist->Fill(jet_truth_matched.pt() * MeVtoGeV);
389 matched_truth_jet_rapidity_hist->Fill(jet_truth_matched.rapidity());
390
391 // Smear the jets for each systematic variation
392 if ( ffjetsmearingtool.applyCorrection(*jet_reco) == CP::CorrectionCode::Error ) {
393 std::cout << "FFJetSmearingTool reported a EL::StatusCode::FAILURE" << std::endl;
394 return 1;
395 }
396
397 smeared_reco_jet_mass_hist->Fill(jet_reco->m() * MeVtoGeV);
398 smeared_reco_jet_pt_hist->Fill(jet_reco->pt() * MeVtoGeV);
399 smeared_reco_jet_rapidity_hist->Fill(jet_reco->rapidity());
400
401 hist_jet_mass_scale_change_3D->Fill(jet_reco->pt() * MeVtoGeV, aux_original_jet_mass, jet_reco->m() * MeVtoGeV);
402 hist_jet_mass_resolution_change_3D->Fill(jet_reco->pt() * MeVtoGeV, aux_original_jet_mass, TMath::Abs( (jet_reco->m()*MeVtoGeV) - (aux_original_jet_mass) )/aux_original_jet_mass);
403
404 lead_jet = false;
405
406 }
407
408 } // Loop over reco jets
409
410 } // Loop over number of events
411
412 TString output_path = "output/" + sys.name() + output ;
413 TFile *jetmass_histograms = new TFile(output_path,"recreate");
414
415 reco_jet_mass_hist->Write();
416 matched_truth_jet_mass_hist->Write();
417 smeared_reco_jet_mass_hist->Write();
418
419 reco_jet_pt_hist->Write();
420 matched_truth_jet_pt_hist->Write();
421 smeared_reco_jet_pt_hist->Write();
422
423 reco_jet_rapidity_hist->Write();
424 matched_truth_jet_rapidity_hist->Write();
425 smeared_reco_jet_rapidity_hist->Write();
426
428
429 TH2F* hist_jet_mass_scale_change_2D = new TH2F("hist_jet_mass_scale_change_2D","hist_jet_mass_scale_change_2D",numBinsPt,0,upperlimit2,numBinsMass,0,upperlimit1);
430 for (int i=1; i<=hist_jet_mass_scale_change_3D->GetNbinsX(); i++){
431 for (int j=1; j<= hist_jet_mass_scale_change_3D->GetNbinsY(); j++){
432 TH1F* hold = new TH1F("","",numBinsMass,0,upperlimit1);
433 for (int k=1; k<= hist_jet_mass_scale_change_3D->GetNbinsZ(); k++){
434 hold->SetBinContent(k,hist_jet_mass_scale_change_3D->GetBinContent(i,j,k)); // A 2D projection of the TH3
435 }
436 hist_jet_mass_scale_change_2D->SetBinContent(i,j,hold->GetMean()/hist_jet_mass_scale_change_3D->GetYaxis()->GetBinCenter(j));
437 delete hold;
438 }
439 }
440
441 TH2F* hist_jet_mass_resolution_change_2D = new TH2F("hist_jet_mass_resolution_change_2D","hist_jet_mass_resolution_change_2D",numBinsPt,0,upperlimit2,numBinsMass,0,upperlimit1);
442 for (int i=1; i<=hist_jet_mass_resolution_change_3D->GetNbinsX(); i++){
443 for (int j=1; j<= hist_jet_mass_resolution_change_3D->GetNbinsY(); j++){
444 TH1F* hold = new TH1F("","",numBinsDiff,lowerlimit3,upperlimit3);
445 for (int k=1; k<= hist_jet_mass_resolution_change_3D->GetNbinsZ(); k++){
446 hold->SetBinContent(k,hist_jet_mass_resolution_change_3D->GetBinContent(i,j,k)); // A 2D projection of the TH3
447 }
448 hist_jet_mass_resolution_change_2D->SetBinContent(i,j,hold->GetMean()/hist_jet_mass_resolution_change_3D->GetYaxis()->GetBinCenter(j));
449 delete hold;
450 }
451 }
452
454
455 double w = 650;
456 double h = 650;
457
458 TCanvas *c1 = new TCanvas("c1","c1",w,h);
459 c1->SetWindowSize(w+(w-c1->GetWw()),h+(h-c1->GetWh()));
460 c1->cd();
461
462 TPad *canvas_1 = new TPad("canvas_1", "canvas_1",0.0,0.0,1.0,1.0);
463 canvas_1->SetRightMargin(0.10);
464 canvas_1->SetFillStyle(4000);
465 canvas_1->Draw();
466 canvas_1->cd();
467
468 hist_jet_mass_scale_change_2D->GetXaxis()->SetTitleOffset(1.5);
469 hist_jet_mass_scale_change_2D->GetXaxis()->SetNdivisions(5);
470 hist_jet_mass_scale_change_2D->GetYaxis()->SetTitleOffset(1.5);
471 hist_jet_mass_scale_change_2D->GetZaxis()->SetTitleOffset(1.5);
472 hist_jet_mass_scale_change_2D->GetYaxis()->SetTitle("Initial Reco Mass [GeV]");
473 hist_jet_mass_scale_change_2D->GetXaxis()->SetTitle("Reco p_{T} [GeV]");
474 hist_jet_mass_scale_change_2D->GetZaxis()->SetTitle("Average Mass smearing (Initial_reco_mass/smeared_reco_mass)");
475 hist_jet_mass_scale_change_2D->GetZaxis()->SetRangeUser(0.95,1.05);
476 hist_jet_mass_scale_change_2D->GetZaxis()->SetLabelSize(0.035);
477 hist_jet_mass_scale_change_2D->Draw("colz");
478 gPad->RedrawAxis();
479
480 TString output_path_scale_debug = "output/debug_plots/scale_variations/" + sys.name() + "_scaleDebug.pdf" ;
481 c1->Print(output_path_scale_debug);
482
483 delete hist_jet_mass_scale_change_2D;
484 delete c1;
485
486 TCanvas *c2 = new TCanvas("c2","c2",w,h);
487 c2->SetWindowSize(w+(w-c2->GetWw()),h+(h-c2->GetWh()));
488 c2->cd();
489
490 TPad *canvas_2 = new TPad("canvas_2", "canvas_2",0.0,0.0,1.0,1.0);
491 canvas_2->SetRightMargin(0.10);
492 canvas_2->SetFillStyle(4000);
493 canvas_2->Draw();
494 canvas_2->cd();
495
496 hist_jet_mass_resolution_change_2D->GetXaxis()->SetTitleOffset(1.5);
497 hist_jet_mass_resolution_change_2D->GetXaxis()->SetNdivisions(5);
498 hist_jet_mass_resolution_change_2D->GetYaxis()->SetTitleOffset(1.5);
499 hist_jet_mass_resolution_change_2D->GetZaxis()->SetTitleOffset(1.5);
500 hist_jet_mass_resolution_change_2D->GetYaxis()->SetTitle("Initial Reco Mass [GeV]");
501 hist_jet_mass_resolution_change_2D->GetXaxis()->SetTitle("Reco p_{T} [GeV]");
502 hist_jet_mass_resolution_change_2D->GetZaxis()->SetTitle("Average Mass smearing (Initial_reco_mass/smeared_reco_mass)");
503 hist_jet_mass_resolution_change_2D->GetZaxis()->SetRangeUser(-0.1,0.1);
504 hist_jet_mass_resolution_change_2D->GetZaxis()->SetLabelSize(0.035);
505 hist_jet_mass_resolution_change_2D->Draw("colz");
506 gPad->RedrawAxis();
507
508 TString output_path_resolution_debug = "output/debug_plots/resolution_variations/" + sys.name() + "_resolutionDebug.pdf" ;
509 c2->Print(output_path_resolution_debug);
510
511 delete hist_jet_mass_resolution_change_2D;
512 delete c2;
513
514 jetmass_histograms->Close();
515 delete jetmass_histograms;
516
517 delete reco_jet_mass_hist;
518 delete matched_truth_jet_mass_hist;
519 delete smeared_reco_jet_mass_hist;
520
521 delete reco_jet_pt_hist;
522 delete matched_truth_jet_pt_hist;
523 delete smeared_reco_jet_pt_hist;
524
525 delete reco_jet_rapidity_hist;
526 delete matched_truth_jet_rapidity_hist;
527 delete smeared_reco_jet_rapidity_hist;
528
529 } // Loop over systematic uncertainties
530
531 return 0;
532
533}
#define CHECK(ARG)
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
Header file for AthHistogramAlgorithm.
@ Error
Some error happened during the object correction.
virtual CP::SystematicSet recommendedSystematics() const override
List of all systematics recommended for this tool.
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
virtual StatusCode applySystematicVariation(const CP::SystematicSet &systematics) override
Configure tool to apply systematic variation.
StatusCode getMatchedTruthJet(xAOD::Jet &jet_reco, xAOD::Jet &jet_truth_matched) const
virtual CP::CorrectionCode applyCorrection(xAOD::Jet &jet_reco) const override
Apply a systematic variation of get a new copy.
Class to wrap a set of SystematicVariations.
const_iterator end() const
description: const iterator to the end of the set
const_iterator begin() const
description: const iterator to the beginning of the set
size_type size() const noexcept
Returns the number of elements in the collection.
StatusCode initialize() override
Dummy implementation of the initialisation function.
StatusCode applyCalibration(xAOD::JetContainer &) const override
Apply calibration to a jet container.
uint32_t runNumber() const
The current event's run number.
uint64_t eventNumber() const
The current event's event number.
static std::unique_ptr< Event > createAndReadFrom(TFile &file)
static method to create an Event object and readFrom a file, given by a TFile.
Definition EventIO.cxx:43
virtual double pt() const
The transverse momentum ( ) of the particle.
Definition Jet_v1.cxx:44
virtual double rapidity() const
The true rapidity (y) of the particle.
Definition Jet_v1.cxx:67
virtual double m() const
The invariant mass of the particle.
Definition Jet_v1.cxx:59
int main()
Definition hello.cxx:18
Jet_v1 Jet
Definition of the current "jet version".
EventInfo_v1 EventInfo
Definition of the latest event info version.
ShallowCopyResult_t< T > shallowCopy(const T &cont, const EventContext &ctx)
Create a shallow copy of an existing container.
JetContainer_v1 JetContainer
Definition of the current "jet container version".