ATLAS Offline Software
Loading...
Searching...
No Matches
JetCalibTools_PlotJESFactors.cxx File Reference
#include "JetCalibTools/JetCalibrationTool.h"
#include "xAODJet/Jet.h"
#include "xAODJet/JetContainer.h"
#include "xAODJet/JetAuxContainer.h"
#include "xAODRootAccess/TEvent.h"
#include "xAODRootAccess/TStore.h"
#include "AsgTools/StandaloneToolHandle.h"
#include <vector>
#include "TString.h"
#include "TH2D.h"
#include "TCanvas.h"

Go to the source code of this file.

Classes

class  jet::JetFourMomAccessor
 JetFourMomAccessor is an extension of JetAttributeAccessor::AccessorWrapper<xAOD::JetFourMom_t> AccessorWrapper<xAOD::JetFourMom_t> purpose is to provide a direct and simple access to JetFourMom_t attributes (which are internally saved as 4 floats inside jets). More...

Namespaces

namespace  jet

Functions

int main (int argc, char *argv[])

Function Documentation

◆ main()

int main ( int argc,
char * argv[] )

Definition at line 41 of file JetCalibTools_PlotJESFactors.cxx.

42{
43 // Check argument usage
44 if (argc < 4 || argc > 6)
45 {
46 std::cout << "USAGE: " << argv[0] << " <JetCollection> <ConfigFile> <OutputFile> (dev mode switch) (JES_vs_E switch)" << std::endl;
47 return 1;
48 }
49
50 // Store the arguments
51 const TString jetAlgo = argv[1];
52 const TString config = argv[2];
53 const TString outFile = argv[3];
54 const bool isDevMode = argc > 4 && (TString(argv[4]) == "true" || TString(argv[4]) == "dev");
55 const bool vsE = argc > 5 && (TString(argv[5]) == "true");
56
57 // Derived information
58 const bool outFileIsExtensible = outFile.EndsWith(".pdf") || outFile.EndsWith(".ps") || outFile.EndsWith(".root");
59
60 // Assumed constants
61 const TString calibSeq = "EtaJES"; // only want to apply the JES here
62 const bool isData = false; // doesn't actually matter for JES, which is always active
63 float massForScan = 80.385e3; // W-boson
64 if ( !jetAlgo.Contains("AntiKt10") ) massForScan = 0;
65
66 // Accessor strings
67 TString startingScaleString = "JetConstitScaleMomentum";
68 if ( !jetAlgo.Contains("AntiKt10") ) startingScaleString = "JetPileupScaleMomentum";
69 const TString endingScaleString = "JetEtaJESScaleMomentum";
70 const TString detectorEtaString = "DetectorEta";
71
72 // Accessors
73 jet::JetFourMomAccessor startingScale(startingScaleString.Data());
74 jet::JetFourMomAccessor endingScale(endingScaleString.Data());
75 SG::AuxElement::Accessor<float> detectorEta(detectorEtaString.Data());
76
77
78 // Create the calib tool
80 calibTool.setTypeAndName("JetCalibrationTool/MyJetCalibTool");
81 {
82 if (calibTool.initialize().isFailure())
83 {
84 std::cout << "Failed to make ana tool" << std::endl;
85 return 2;
86 }
87 if (calibTool.setProperty("JetCollection",jetAlgo.Data()).isFailure())
88 {
89 std::cout << "Failed to set JetCollection: " << jetAlgo.Data() << std::endl;
90 return 3;
91 }
92 if (calibTool.setProperty("ConfigFile",config.Data()).isFailure())
93 {
94 std::cout << "Failed to set ConfigFile: " << config.Data() << std::endl;
95 return 3;
96 }
97 if (calibTool.setProperty("CalibSequence",calibSeq.Data()).isFailure())
98 {
99 std::cout << "Failed to set CalibSequence: " << calibSeq.Data() << std::endl;
100 return 3;
101 }
102 if (calibTool.setProperty("IsData",isData).isFailure())
103 {
104 std::cout << "Failed to set IsData: " << (isData ? std::string("true") : std::string("false")) << std::endl;
105 return 3;
106 }
107 if (isDevMode && calibTool.setProperty("DEVmode",isDevMode).isFailure())
108 {
109 std::cout << "Failed to set DEVmode" << std::endl;
110 return 4;
111 }
112 if (calibTool.retrieve().isFailure())
113 {
114 std::cout << "Failed to initialize the JetCalibTool" << std::endl;
115 return 5;
116 }
117 }
118
119
120 // Build a jet container and a jet for us to manipulate later
124 jets->setStore(new xAOD::JetAuxContainer());
125 jets->push_back(new xAOD::Jet());
126 xAOD::Jet* jet = jets->at(0);
127
128
129 // Make the histogram to fill
130 TH2D* hist_pt_eta;
131 if ( jetAlgo.Contains("AntiKt10") ) hist_pt_eta = new TH2D("JES_pt_eta",Form("JES for jets with mass=%.1f GeV",massForScan/1.e3),1200,100,2500,60,-3,3);
132 else { hist_pt_eta = new TH2D("JES_pt_eta",Form("JES for jets with mass=%.1f GeV",massForScan/1.e3),2500,20,5000,90,-4.5,4.5); }
133
134 // Fill the histogram
135 for (int xBin = 1; xBin <= hist_pt_eta->GetNbinsX(); ++xBin)
136 {
137 const double pt = hist_pt_eta->GetXaxis()->GetBinCenter(xBin)*1.e3; // E if vsE
138 for (int yBin = 1; yBin <= hist_pt_eta->GetNbinsY(); ++yBin)
139 {
140 const double eta = hist_pt_eta->GetYaxis()->GetBinCenter(yBin);
141
142 // Set the main 4-vector and scale 4-vector
143 if ( !vsE ){
144 jet->setJetP4(xAOD::JetFourMom_t(pt,eta,0,massForScan));
145 detectorEta(*jet) = eta;
146 startingScale.setAttribute(*jet,xAOD::JetFourMom_t(pt,eta,0,massForScan));
147 } else {
148 const double E = pt; // pt is actually E if vsE
149 const double pT = sqrt((E*E)-(massForScan*massForScan))/cosh(eta);
150 jet->setJetP4(xAOD::JetFourMom_t(pT,eta,0,massForScan));
151 detectorEta(*jet) = eta;
152 startingScale.setAttribute(*jet,xAOD::JetFourMom_t(pT,eta,0,massForScan));
153 }
154 const double startingE = startingScale(*jet).e();
155
156 // Jet kinematics set, now apply calibration
157 if(calibTool->modify(*jets).isFailure()){
158 std::cout << "Failed to apply jet calibration" << std::endl;
159 return 6;
160 }
161 // Calculate the scale factors
162 const double JES = jet->e()/startingE;
163
164 // JMS retrieved, fill the plot(s)
165 hist_pt_eta->SetBinContent(xBin,yBin,JES);
166 }
167 }
168
169
170
171 // Make the plots
172
173 // First the canvas
174 TCanvas* canvas = new TCanvas("canvas");
175 canvas->SetMargin(0.07,0.13,0.1,0.10);
176 canvas->SetFillStyle(4000);
177 canvas->SetFillColor(0);
178 canvas->SetFrameBorderMode(0);
179 canvas->cd();
180 canvas->SetLogx(true);
181
182 // Now labels/etc
183 hist_pt_eta->SetStats(false);
184 if ( !vsE ) hist_pt_eta->GetXaxis()->SetTitle("Jet #it{p}_{T} [GeV]");
185 else { hist_pt_eta->GetXaxis()->SetTitle("Jet E [GeV]"); }
186 hist_pt_eta->GetXaxis()->SetTitleOffset(1.35);
187 hist_pt_eta->GetXaxis()->SetMoreLogLabels();
188 hist_pt_eta->GetYaxis()->SetTitle("#eta");
189 hist_pt_eta->GetYaxis()->SetTitleOffset(0.9);
190 hist_pt_eta->GetZaxis()->SetTitle("JES_{Factor}");
191
192 // Now write them out
193 if (outFileIsExtensible)
194 {
195
196 if ( outFile.EndsWith(".pdf") || outFile.EndsWith(".ps") ){
197 canvas->Print(outFile+"[");
198 hist_pt_eta->Draw("colz");
199 canvas->Print(outFile);
200 canvas->Print(outFile+"]");
201 } else if ( outFile.EndsWith(".root") ){
202 TFile *fout = new TFile(outFile.Data(),"RECREATE");
203 hist_pt_eta->Write();
204 fout->Close();
205 delete fout;
206 }
207
208 }
209 else
210 {
211 unsigned int counter = 1;
212 hist_pt_eta->Draw("colz");
213 canvas->Print(Form("%u-fixMass-%s",counter++,outFile.Data()));
214 }
215
216 delete hist_pt_eta;
217
218 return 0;
219}
Scalar eta() const
pseudorapidity method
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
an "initializing" ToolHandle for stand-alone applications
StatusCode setProperty(const std::string &name, T2 &&value)
StatusCode retrieve()
initialize the tool, will succeed if the tool was already initialized
void setTypeAndName(const std::string &typeAndName)
StatusCode initialize()
initialize the tool, will fail if the tool was already initialized
JetFourMomAccessor is an extension of JetAttributeAccessor::AccessorWrapper<xAOD::JetFourMom_t> Acces...
STL class.
Tool for accessing xAOD files outside of Athena.
A relatively simple transient store for objects created in analysis.
Definition TStore.h:45
static TFile * fout
Definition listroot.cxx:40
outFile
Comment Out Those You do not wish to run.
TestStore store
Definition TestStore.cxx:23
STL namespace.
Jet_v1 Jet
Definition of the current "jet version".
JetAuxContainer_v1 JetAuxContainer
Definition of the current jet auxiliary container.
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