[cf88574] | 1 | #ifdef __CLING__
|
---|
| 2 | R__LOAD_LIBRARY(libDelphes)
|
---|
| 3 | #include "classes/DelphesClasses.h"
|
---|
| 4 | #include "external/ExRootAnalysis/ExRootTreeReader.h"
|
---|
| 5 | #include "external/ExRootAnalysis/ExRootResult.h"
|
---|
| 6 | #else
|
---|
| 7 | class ExRootTreeReader;
|
---|
| 8 | class ExRootResult;
|
---|
| 9 | #endif
|
---|
| 10 |
|
---|
| 11 | #include "TCanvas.h"
|
---|
| 12 | #include "TSystem.h"
|
---|
| 13 | #include "TStyle.h"
|
---|
| 14 | #include "TLegend.h"
|
---|
| 15 | #include <TH1.h>
|
---|
| 16 | #include "TString.h"
|
---|
| 17 | #include "vector"
|
---|
| 18 | #include <TMath.h>
|
---|
| 19 | #include <iostream>
|
---|
| 20 | #include "TGraph.h"
|
---|
| 21 | #include "TGraphErrors.h"
|
---|
| 22 | #include "TMultiGraph.h"
|
---|
| 23 | #include <typeinfo>
|
---|
| 24 | #include "TLorentzVector.h"
|
---|
| 25 |
|
---|
| 26 | //------------------------------------------------------------------------------
|
---|
| 27 |
|
---|
| 28 | double ptrangemin = 10;
|
---|
| 29 | double ptrangemax = 10000;
|
---|
| 30 | static const int Nbins = 20;
|
---|
| 31 |
|
---|
| 32 | int objStyle = 1;
|
---|
| 33 | int trackStyle = 7;
|
---|
| 34 | int towerStyle = 5;
|
---|
| 35 |
|
---|
| 36 | Color_t objColor = kBlack;
|
---|
| 37 | Color_t trackColor = kBlack;
|
---|
| 38 | Color_t towerColor = kBlack;
|
---|
| 39 |
|
---|
| 40 | double effLegXmin = 0.22;
|
---|
| 41 | double effLegXmax = 0.5;
|
---|
| 42 | double effLegYmin = 0.22;
|
---|
| 43 | double effLegYmax = 0.4;
|
---|
| 44 |
|
---|
| 45 | struct resolPlot
|
---|
| 46 | {
|
---|
| 47 | TH1 *cenResolHist;
|
---|
| 48 | TH1 *fwdResolHist;
|
---|
| 49 | double ptmin;
|
---|
| 50 | double ptmax;
|
---|
[6563d4d] | 51 | double xmin;
|
---|
| 52 | double xmax;
|
---|
[cf88574] | 53 | TString obj;
|
---|
| 54 |
|
---|
| 55 | resolPlot();
|
---|
| 56 | resolPlot(double ptdown, double ptup, TString object);
|
---|
[6563d4d] | 57 | void set(double ptdown, double ptup, TString object, double xmin = 0, double xmax = 2);
|
---|
[cf88574] | 58 | void print(){std::cout << ptmin << std::endl;}
|
---|
| 59 | };
|
---|
| 60 |
|
---|
| 61 |
|
---|
| 62 | resolPlot::resolPlot()
|
---|
| 63 | {
|
---|
| 64 | }
|
---|
| 65 |
|
---|
| 66 | resolPlot::resolPlot(double ptdown, double ptup, TString object)
|
---|
| 67 | {
|
---|
| 68 | this->set(ptdown,ptup,object);
|
---|
| 69 | }
|
---|
| 70 |
|
---|
[6563d4d] | 71 | void resolPlot::set(double ptdown, double ptup, TString object, double xmin, double xmax){
|
---|
[cf88574] | 72 | ptmin = ptdown;
|
---|
| 73 | ptmax = ptup;
|
---|
| 74 | obj = object;
|
---|
| 75 |
|
---|
[6563d4d] | 76 | cenResolHist = new TH1D(obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_cen", obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_cen", 500, xmin, xmax);
|
---|
| 77 | fwdResolHist = new TH1D(obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_fwd", obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_fwd", 500, xmin, xmax);
|
---|
[cf88574] | 78 |
|
---|
| 79 | }
|
---|
| 80 |
|
---|
[6563d4d] | 81 | void HistogramsCollection(std::vector<resolPlot> *histos, double ptmin, double ptmax, TString obj, double xmin = 0, double xmax = 2)
|
---|
[cf88574] | 82 | {
|
---|
| 83 | double width;
|
---|
| 84 | double ptdown;
|
---|
| 85 | double ptup;
|
---|
| 86 | resolPlot ptemp;
|
---|
| 87 |
|
---|
| 88 | for (int i = 0; i < Nbins; i++)
|
---|
| 89 | {
|
---|
| 90 | width = (ptmax - ptmin) / Nbins;
|
---|
| 91 | ptdown = TMath::Power(10,ptmin + i * width );
|
---|
| 92 | ptup = TMath::Power(10,ptmin + (i+1) * width );
|
---|
[6563d4d] | 93 | ptemp.set(ptdown, ptup, obj, xmin, xmax);
|
---|
[cf88574] | 94 | histos->push_back(ptemp);
|
---|
| 95 | }
|
---|
| 96 | }
|
---|
| 97 |
|
---|
| 98 | //------------------------------------------------------------------------------
|
---|
| 99 |
|
---|
| 100 | class ExRootResult;
|
---|
| 101 | class ExRootTreeReader;
|
---|
| 102 |
|
---|
| 103 | //------------------------------------------------------------------------------
|
---|
| 104 |
|
---|
| 105 | void BinLogX(TH1*h)
|
---|
| 106 | {
|
---|
| 107 |
|
---|
| 108 | TAxis *axis = h->GetXaxis();
|
---|
| 109 | int bins = axis->GetNbins();
|
---|
| 110 |
|
---|
| 111 | Axis_t from = axis->GetXmin();
|
---|
| 112 | Axis_t to = axis->GetXmax();
|
---|
| 113 | Axis_t width = (to - from) / bins;
|
---|
| 114 | Axis_t *new_bins = new Axis_t[bins + 1];
|
---|
| 115 |
|
---|
| 116 | for (int i = 0; i <= bins; i++) {
|
---|
| 117 | new_bins[i] = TMath::Power(10, from + i * width);
|
---|
| 118 |
|
---|
| 119 | }
|
---|
| 120 | axis->Set(bins, new_bins);
|
---|
| 121 | delete new_bins;
|
---|
| 122 | }
|
---|
| 123 |
|
---|
| 124 |
|
---|
| 125 | //------------------------------------------------------------------------------
|
---|
| 126 |
|
---|
| 127 | template<typename T>
|
---|
| 128 | std::pair<TH1D*, TH1D*> GetEff(TClonesArray *branchReco, TClonesArray *branchParticle, TString name, int pdgID, ExRootTreeReader *treeReader)
|
---|
| 129 | {
|
---|
| 130 |
|
---|
| 131 | cout << "** Computing Efficiency of reconstructing "<< branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
|
---|
| 132 |
|
---|
| 133 | Long64_t allEntries = treeReader->GetEntries();
|
---|
| 134 |
|
---|
| 135 | GenParticle *particle;
|
---|
| 136 | T *recoObj;
|
---|
| 137 |
|
---|
| 138 | TLorentzVector recoMomentum, genMomentum, bestRecoMomentum;
|
---|
| 139 |
|
---|
| 140 | Float_t deltaR;
|
---|
| 141 | Float_t pt, eta;
|
---|
| 142 | Long64_t entry;
|
---|
| 143 |
|
---|
| 144 | Int_t i, j;
|
---|
| 145 |
|
---|
[20ca0cd] | 146 | TH1D *histGenPtcen = new TH1D(name+" gen spectra Pt",name+" gen spectra cen", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
|
---|
| 147 | TH1D *histRecoPtcen = new TH1D(name+" reco spectra Pt",name+" reco spectra cen", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
|
---|
| 148 | TH1D *histGenPtfwd = new TH1D(name+" gen spectra Eta",name+" gen spectra fwd", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
|
---|
| 149 | TH1D *histRecoPtfwd = new TH1D(name+" reco spectra Eta",name+" reco spectra fwd", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
|
---|
[cf88574] | 150 |
|
---|
| 151 |
|
---|
[20ca0cd] | 152 | BinLogX(histGenPtcen);
|
---|
| 153 | BinLogX(histRecoPtcen);
|
---|
| 154 | BinLogX(histGenPtfwd);
|
---|
| 155 | BinLogX(histRecoPtfwd);
|
---|
[cf88574] | 156 |
|
---|
| 157 | // Loop over all events
|
---|
| 158 | for(entry = 0; entry < allEntries; ++entry)
|
---|
| 159 | {
|
---|
| 160 | // Load selected branches with data from specified event
|
---|
| 161 | treeReader->ReadEntry(entry);
|
---|
| 162 |
|
---|
| 163 | // Loop over all generated particle in event
|
---|
| 164 | for(i = 0; i < branchParticle->GetEntriesFast(); ++i)
|
---|
| 165 | {
|
---|
| 166 |
|
---|
| 167 | particle = (GenParticle*) branchParticle->At(i);
|
---|
| 168 | genMomentum = particle->P4();
|
---|
| 169 |
|
---|
| 170 | deltaR = 999;
|
---|
| 171 |
|
---|
| 172 | if (particle->PID == pdgID && genMomentum.Pt() > ptrangemin && genMomentum.Pt() < ptrangemax )
|
---|
| 173 | {
|
---|
| 174 |
|
---|
| 175 | // Loop over all reco object in event
|
---|
| 176 | for(j = 0; j < branchReco->GetEntriesFast(); ++j)
|
---|
| 177 | {
|
---|
| 178 | recoObj = (T*)branchReco->At(j);
|
---|
| 179 | recoMomentum = recoObj->P4();
|
---|
| 180 | // this is simply to avoid warnings from initial state particle
|
---|
| 181 | // having infite rapidity ...
|
---|
| 182 | //if(Momentum.Px() == 0 && genMomentum.Py() == 0) continue;
|
---|
| 183 |
|
---|
| 184 | // take the closest parton candidate
|
---|
| 185 | if(TMath::Abs(pdgID) == 5)
|
---|
| 186 | {
|
---|
| 187 | Jet *jet = (Jet *)recoObj;
|
---|
| 188 | if(jet->BTag != 1) continue;
|
---|
| 189 | }
|
---|
| 190 | if(TMath::Abs(pdgID) == 15)
|
---|
| 191 | {
|
---|
| 192 | Jet *jet = (Jet *)recoObj;
|
---|
| 193 | if(jet->TauTag != 1) continue;
|
---|
| 194 | }
|
---|
| 195 | if(genMomentum.DeltaR(recoMomentum) < deltaR)
|
---|
| 196 | {
|
---|
| 197 | deltaR = genMomentum.DeltaR(recoMomentum);
|
---|
| 198 | bestRecoMomentum = recoMomentum;
|
---|
| 199 | }
|
---|
| 200 | }
|
---|
| 201 |
|
---|
| 202 | pt = genMomentum.Pt();
|
---|
| 203 | eta = genMomentum.Eta();
|
---|
| 204 |
|
---|
[20ca0cd] | 205 | if (TMath::Abs(eta) < 1.5)
|
---|
| 206 | {
|
---|
| 207 | histGenPtcen->Fill(pt);
|
---|
| 208 | if(deltaR < 0.3) { histRecoPtcen->Fill(pt); }
|
---|
| 209 | }
|
---|
| 210 | else if (TMath::Abs(eta) < 2.5)
|
---|
[cf88574] | 211 | {
|
---|
[20ca0cd] | 212 | histGenPtfwd->Fill(pt);
|
---|
| 213 | if(deltaR < 0.3) { histRecoPtfwd->Fill(pt); }
|
---|
| 214 |
|
---|
[cf88574] | 215 | }
|
---|
| 216 | }
|
---|
| 217 | }
|
---|
| 218 | }
|
---|
| 219 |
|
---|
| 220 |
|
---|
| 221 | std::pair<TH1D*,TH1D*> histos;
|
---|
| 222 |
|
---|
[20ca0cd] | 223 | histRecoPtcen->Divide(histGenPtcen);
|
---|
| 224 | histRecoPtfwd->Divide(histGenPtfwd);
|
---|
[cf88574] | 225 |
|
---|
[20ca0cd] | 226 | histos.first = histRecoPtcen;
|
---|
| 227 | histos.second = histRecoPtfwd;
|
---|
[cf88574] | 228 |
|
---|
| 229 | return histos;
|
---|
| 230 | }
|
---|
| 231 |
|
---|
| 232 | template<typename T>
|
---|
| 233 | void GetEres(std::vector<resolPlot> *histos, TClonesArray *branchReco, TClonesArray *branchParticle, int pdgID, ExRootTreeReader *treeReader)
|
---|
| 234 | {
|
---|
| 235 | Long64_t allEntries = treeReader->GetEntries();
|
---|
| 236 |
|
---|
| 237 | cout << "** Computing resolution of " << branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
|
---|
| 238 |
|
---|
| 239 | GenParticle *particle;
|
---|
| 240 | T* recoObj;
|
---|
| 241 |
|
---|
| 242 | TLorentzVector recoMomentum, genMomentum, bestGenMomentum;
|
---|
| 243 |
|
---|
| 244 | Float_t deltaR;
|
---|
| 245 | Float_t pt, eta;
|
---|
| 246 | Long64_t entry;
|
---|
| 247 |
|
---|
| 248 | Int_t i, j, bin;
|
---|
| 249 |
|
---|
| 250 | // Loop over all events
|
---|
| 251 | for(entry = 0; entry < allEntries; ++entry)
|
---|
| 252 | {
|
---|
| 253 | // Load selected branches with data from specified event
|
---|
| 254 | treeReader->ReadEntry(entry);
|
---|
| 255 |
|
---|
| 256 | // Loop over all reconstructed jets in event
|
---|
| 257 | for(i = 0; i < branchReco->GetEntriesFast(); ++i)
|
---|
| 258 | {
|
---|
| 259 | recoObj = (T*) branchReco->At(i);
|
---|
| 260 | recoMomentum = recoObj->P4();
|
---|
| 261 |
|
---|
| 262 | deltaR = 999;
|
---|
| 263 |
|
---|
| 264 | // Loop over all hard partons in event
|
---|
| 265 | for(j = 0; j < branchParticle->GetEntriesFast(); ++j)
|
---|
| 266 | {
|
---|
| 267 | particle = (GenParticle*) branchParticle->At(j);
|
---|
| 268 | if (particle->PID == pdgID && particle->Status == 1)
|
---|
| 269 | {
|
---|
| 270 | genMomentum = particle->P4();
|
---|
| 271 |
|
---|
| 272 | // this is simply to avoid warnings from initial state particle
|
---|
| 273 | // having infite rapidity ...
|
---|
| 274 | if(genMomentum.Px() == 0 && genMomentum.Py() == 0) continue;
|
---|
| 275 |
|
---|
| 276 | // take the closest parton candidate
|
---|
| 277 | if(genMomentum.DeltaR(recoMomentum) < deltaR)
|
---|
| 278 | {
|
---|
| 279 | deltaR = genMomentum.DeltaR(recoMomentum);
|
---|
| 280 | bestGenMomentum = genMomentum;
|
---|
| 281 | }
|
---|
| 282 | }
|
---|
| 283 | }
|
---|
| 284 |
|
---|
| 285 | if(deltaR < 0.3)
|
---|
| 286 | {
|
---|
| 287 | pt = bestGenMomentum.Pt();
|
---|
| 288 | eta = TMath::Abs(bestGenMomentum.Eta());
|
---|
| 289 |
|
---|
| 290 | for (bin = 0; bin < Nbins; bin++)
|
---|
| 291 | {
|
---|
| 292 | if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta > 0.0 && eta < 2.5)
|
---|
| 293 | {
|
---|
| 294 | if (eta < 1.5) {histos->at(bin).cenResolHist->Fill(recoMomentum.Pt()/bestGenMomentum.Pt());}
|
---|
| 295 | else if (eta < 2.5) {histos->at(bin).fwdResolHist->Fill(recoMomentum.Pt()/bestGenMomentum.Pt());}
|
---|
| 296 | }
|
---|
| 297 | }
|
---|
| 298 | }
|
---|
| 299 | }
|
---|
| 300 | }
|
---|
| 301 | }
|
---|
| 302 | void GetJetsEres(std::vector<resolPlot> *histos, TClonesArray *branchJet, TClonesArray *branchGenJet, ExRootTreeReader *treeReader)
|
---|
| 303 | {
|
---|
| 304 |
|
---|
| 305 | Long64_t allEntries = treeReader->GetEntries();
|
---|
| 306 |
|
---|
| 307 | cout << "** Computing resolution of " << branchJet->GetName() << " induced by " << branchGenJet->GetName() << endl;
|
---|
| 308 |
|
---|
| 309 | Jet *jet, *genjet;
|
---|
| 310 |
|
---|
| 311 | TLorentzVector jetMomentum, genJetMomentum, bestGenJetMomentum;
|
---|
| 312 |
|
---|
| 313 | Float_t deltaR;
|
---|
| 314 | Float_t pt, eta;
|
---|
| 315 | Long64_t entry;
|
---|
| 316 |
|
---|
| 317 | Int_t i, j, bin;
|
---|
| 318 |
|
---|
| 319 | // Loop over all events
|
---|
| 320 | for(entry = 0; entry < allEntries; ++entry)
|
---|
| 321 | {
|
---|
| 322 | // Load selected branches with data from specified event
|
---|
| 323 | treeReader->ReadEntry(entry);
|
---|
| 324 |
|
---|
| 325 | if(entry%10000 == 0) cout << "Event number: "<< entry <<endl;
|
---|
| 326 |
|
---|
| 327 | // Loop over all reconstructed jets in event
|
---|
| 328 | for(i = 0; i < TMath::Min(2,branchJet->GetEntriesFast()); ++i) //branchJet->GetEntriesFast(); ++i)
|
---|
| 329 | {
|
---|
| 330 |
|
---|
| 331 | jet = (Jet*) branchJet->At(i);
|
---|
| 332 | jetMomentum = jet->P4();
|
---|
| 333 |
|
---|
| 334 | deltaR = 999;
|
---|
| 335 |
|
---|
| 336 | // Loop over all hard partons in event
|
---|
| 337 | for(j = 0; j < TMath::Min(2,branchGenJet->GetEntriesFast()); ++j)
|
---|
| 338 | {
|
---|
| 339 | genjet = (Jet*) branchGenJet->At(j);
|
---|
| 340 |
|
---|
| 341 | genJetMomentum = genjet->P4();
|
---|
| 342 |
|
---|
| 343 | // this is simply to avoid warnings from initial state particle
|
---|
| 344 | // having infite rapidity ...
|
---|
| 345 | if(genJetMomentum.Px() == 0 && genJetMomentum.Py() == 0) continue;
|
---|
| 346 |
|
---|
| 347 | // take the closest parton candidate
|
---|
| 348 | if(genJetMomentum.DeltaR(jetMomentum) < deltaR)
|
---|
| 349 | {
|
---|
| 350 | deltaR = genJetMomentum.DeltaR(jetMomentum);
|
---|
| 351 | bestGenJetMomentum = genJetMomentum;
|
---|
| 352 | }
|
---|
| 353 | }
|
---|
| 354 |
|
---|
| 355 | if(deltaR < 0.25)
|
---|
| 356 | {
|
---|
| 357 | pt = genJetMomentum.Pt();
|
---|
| 358 | eta = TMath::Abs(genJetMomentum.Eta());
|
---|
| 359 |
|
---|
| 360 | for (bin = 0; bin < Nbins; bin++)
|
---|
| 361 | {
|
---|
[20ca0cd] | 362 | if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta < 1.5)
|
---|
[cf88574] | 363 | {
|
---|
| 364 | histos->at(bin).cenResolHist->Fill(jetMomentum.Pt()/bestGenJetMomentum.Pt());
|
---|
| 365 | }
|
---|
[20ca0cd] | 366 | else if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta < 2.5)
|
---|
| 367 | {
|
---|
| 368 | histos->at(bin).fwdResolHist->Fill(jetMomentum.Pt()/bestGenJetMomentum.Pt());
|
---|
| 369 | }
|
---|
| 370 |
|
---|
[cf88574] | 371 | }
|
---|
| 372 | }
|
---|
| 373 | }
|
---|
| 374 | }
|
---|
| 375 | }
|
---|
[6563d4d] | 376 |
|
---|
| 377 | void GetMetres(std::vector<resolPlot> *histos, TClonesArray *branchScalarHT, TClonesArray *branchMet, TClonesArray *branchJet, ExRootTreeReader *treeReader)
|
---|
| 378 | {
|
---|
| 379 |
|
---|
| 380 | Long64_t allEntries = treeReader->GetEntries();
|
---|
| 381 |
|
---|
| 382 | cout << "** Computing resolution of " << branchMet->GetName() << " vs " << branchScalarHT->GetName() << endl;
|
---|
| 383 |
|
---|
| 384 | MissingET *met;
|
---|
| 385 | ScalarHT *scalarHT;
|
---|
| 386 |
|
---|
| 387 | Long64_t entry;
|
---|
| 388 |
|
---|
| 389 | Int_t bin;
|
---|
| 390 | Double_t ht;
|
---|
| 391 |
|
---|
| 392 | // Loop over all events
|
---|
| 393 | for(entry = 0; entry < allEntries; ++entry)
|
---|
| 394 | {
|
---|
| 395 | // Load selected branches with data from specified event
|
---|
| 396 | treeReader->ReadEntry(entry);
|
---|
| 397 |
|
---|
| 398 | if(entry%10000 == 0) cout << "Event number: "<< entry <<endl;
|
---|
| 399 |
|
---|
| 400 | if (branchJet->GetEntriesFast() > 1)
|
---|
| 401 | {
|
---|
| 402 | met = (MissingET*) branchMet->At(0);
|
---|
| 403 | scalarHT = (ScalarHT*) branchScalarHT->At(0);
|
---|
| 404 | ht = scalarHT->HT;
|
---|
| 405 |
|
---|
| 406 | for (bin = 0; bin < Nbins; bin++)
|
---|
| 407 | {
|
---|
| 408 | if(ht > histos->at(bin).ptmin && ht < histos->at(bin).ptmax )
|
---|
| 409 | {
|
---|
| 410 | histos->at(bin).cenResolHist->Fill(met->P4().Px());
|
---|
| 411 | histos->at(bin).fwdResolHist->Fill(met->P4().Py());
|
---|
| 412 | }
|
---|
| 413 | }
|
---|
| 414 | }
|
---|
| 415 | }
|
---|
| 416 | }
|
---|
| 417 |
|
---|
| 418 |
|
---|
[cf88574] | 419 | std::pair<Double_t, Double_t> GausFit(TH1* hist)
|
---|
| 420 | {
|
---|
| 421 | TF1 *f1 = new TF1("f1", "gaus", hist->GetMean()-2*hist->GetRMS(), hist->GetMean()+2*hist->GetRMS());
|
---|
| 422 | hist->Fit("f1","RQ");
|
---|
| 423 | Double_t sig = f1->GetParameter(2);
|
---|
| 424 | Double_t sigErr = f1->GetParError(2);
|
---|
| 425 | delete f1;
|
---|
| 426 | return make_pair (sig, sigErr);
|
---|
| 427 | //return make_pair (hist->GetRMS(), hist->GetRMSError());
|
---|
| 428 | }
|
---|
| 429 |
|
---|
| 430 |
|
---|
[6563d4d] | 431 | TGraphErrors EresGraph(std::vector<resolPlot> *histos, bool central, bool rms = false)
|
---|
[cf88574] | 432 | {
|
---|
| 433 | Int_t bin;
|
---|
| 434 | Int_t count = 0;
|
---|
| 435 | TGraphErrors gr = TGraphErrors(Nbins/2);
|
---|
| 436 | Double_t sig = 0;
|
---|
| 437 | Double_t sigErr = 0;
|
---|
| 438 | for (bin = 0; bin < Nbins; bin++)
|
---|
| 439 | {
|
---|
| 440 | if (central == true && histos->at(bin).cenResolHist->GetEntries() > 100)
|
---|
| 441 | {
|
---|
| 442 | std::pair<Double_t, Double_t> sigvalues = GausFit(histos->at(bin).cenResolHist);
|
---|
[6563d4d] | 443 | if (rms == true)
|
---|
| 444 | {
|
---|
| 445 | gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.second);
|
---|
| 446 | gr.SetPointError(count,0, sigvalues.second); // to correct
|
---|
| 447 | }
|
---|
| 448 | else
|
---|
| 449 | {
|
---|
| 450 | gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.first);
|
---|
| 451 | gr.SetPointError(count,0, sigvalues.second);
|
---|
| 452 | }
|
---|
[cf88574] | 453 | count++;
|
---|
| 454 | }
|
---|
[20ca0cd] | 455 |
|
---|
| 456 | else if (central == false && histos->at(bin).fwdResolHist->GetEntries() > 10)
|
---|
[cf88574] | 457 | {
|
---|
[20ca0cd] | 458 | std::pair<Double_t, Double_t> sigvalues = GausFit(histos->at(bin).fwdResolHist);
|
---|
| 459 | if (rms == true)
|
---|
| 460 | {
|
---|
| 461 | gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.second);
|
---|
| 462 | gr.SetPointError(count,0, sigvalues.second); // to correct
|
---|
| 463 | }
|
---|
| 464 | else
|
---|
| 465 | {
|
---|
| 466 | gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.first);
|
---|
| 467 | gr.SetPointError(count,0, sigvalues.second);
|
---|
| 468 | }
|
---|
| 469 | count++;
|
---|
[cf88574] | 470 | }
|
---|
[20ca0cd] | 471 |
|
---|
[cf88574] | 472 | }
|
---|
| 473 | return gr;
|
---|
| 474 | }
|
---|
| 475 |
|
---|
| 476 |
|
---|
| 477 | //------------------------------------------------------------------------------
|
---|
| 478 |
|
---|
| 479 |
|
---|
| 480 | // type 1 : object, 2 : track, 3 : tower
|
---|
| 481 |
|
---|
| 482 | void addGraph(TMultiGraph *mg, TGraphErrors *gr, TLegend *leg, int type)
|
---|
| 483 | {
|
---|
| 484 |
|
---|
| 485 | gr->SetLineWidth(3);
|
---|
| 486 |
|
---|
| 487 | switch ( type )
|
---|
| 488 | {
|
---|
| 489 | case 1:
|
---|
| 490 | gr->SetLineColor(objColor);
|
---|
| 491 | gr->SetLineStyle(objStyle);
|
---|
| 492 | std::cout << "Adding " << gr->GetName() << std::endl;
|
---|
| 493 | mg->Add(gr);
|
---|
| 494 | leg->AddEntry(gr,"Reco","l");
|
---|
| 495 | break;
|
---|
| 496 |
|
---|
| 497 | case 2:
|
---|
| 498 | gr->SetLineColor(trackColor);
|
---|
| 499 | gr->SetLineStyle(trackStyle);
|
---|
| 500 | mg->Add(gr);
|
---|
| 501 | leg->AddEntry(gr,"Track","l");
|
---|
| 502 | break;
|
---|
| 503 |
|
---|
| 504 | case 3:
|
---|
| 505 | gr->SetLineColor(towerColor);
|
---|
| 506 | gr->SetLineStyle(towerStyle);
|
---|
| 507 | mg->Add(gr);
|
---|
| 508 | leg->AddEntry(gr,"Tower","l");
|
---|
| 509 | break;
|
---|
| 510 |
|
---|
| 511 | default:
|
---|
| 512 | std::cout << "wrong type, possibles choices are Object, Track and Tower" << std::endl;
|
---|
| 513 | break;
|
---|
| 514 | }
|
---|
| 515 | }
|
---|
| 516 |
|
---|
| 517 | void addHist(TH1D *h, TLegend *leg, int type)
|
---|
| 518 | {
|
---|
| 519 | h->SetLineWidth(3);
|
---|
| 520 |
|
---|
| 521 | switch ( type )
|
---|
| 522 | {
|
---|
| 523 | case 1:
|
---|
| 524 | h->SetLineColor(objColor);
|
---|
| 525 | h->SetLineStyle(objStyle);
|
---|
| 526 | leg->AddEntry(h,"Reco","l");
|
---|
| 527 | break;
|
---|
| 528 |
|
---|
| 529 | case 2:
|
---|
| 530 | h->SetLineColor(trackColor);
|
---|
| 531 | h->SetLineStyle(trackStyle);
|
---|
| 532 | leg->AddEntry(h,"Track","l");
|
---|
| 533 | break;
|
---|
| 534 |
|
---|
| 535 | case 3:
|
---|
| 536 | h->SetLineColor(towerColor);
|
---|
| 537 | h->SetLineStyle(towerStyle);
|
---|
| 538 | leg->AddEntry(h,"Tower","l");
|
---|
| 539 | break;
|
---|
| 540 |
|
---|
| 541 | default:
|
---|
| 542 | std::cout << "wrong type, possibles choices are Object, Track and Tower" << std::endl;
|
---|
| 543 | break;
|
---|
| 544 | }
|
---|
| 545 | }
|
---|
| 546 |
|
---|
| 547 | void DrawAxis(TMultiGraph *mg, TLegend *leg, double max)
|
---|
| 548 | {
|
---|
| 549 | mg->SetMinimum(0.);
|
---|
| 550 | mg->SetMaximum(max);
|
---|
| 551 | mg->GetXaxis()->SetLimits(ptrangemin,ptrangemax);
|
---|
| 552 | mg->GetYaxis()->SetTitle("#sigma (E) / E_{gen}");
|
---|
| 553 | mg->GetXaxis()->SetTitle("E_{gen}");
|
---|
| 554 | mg->GetYaxis()->SetTitleSize(0.07);
|
---|
| 555 | mg->GetXaxis()->SetTitleSize(0.07);
|
---|
| 556 | mg->GetYaxis()->SetLabelSize(0.07);
|
---|
| 557 | mg->GetXaxis()->SetLabelSize(0.07);
|
---|
| 558 |
|
---|
| 559 | leg->SetLineStyle(0);
|
---|
| 560 | leg->SetFillStyle(0);
|
---|
| 561 | leg->SetLineWidth(0);
|
---|
| 562 | leg->SetLineColor(0);
|
---|
| 563 |
|
---|
| 564 | gStyle->SetOptTitle(0);
|
---|
| 565 | gPad->SetLogx();
|
---|
| 566 | gPad->SetBottomMargin(0.2);
|
---|
| 567 | gPad->SetLeftMargin(0.2);
|
---|
| 568 | gPad->Modified();
|
---|
| 569 | gPad->Update();
|
---|
| 570 |
|
---|
| 571 | }
|
---|
| 572 |
|
---|
| 573 | void DrawAxis(TH1D *h, TLegend *leg, int type)
|
---|
| 574 | {
|
---|
| 575 |
|
---|
| 576 | h->GetYaxis()->SetRangeUser(0,1.2);
|
---|
| 577 | if (type == 0) h->GetXaxis()->SetTitle("E_{gen}");
|
---|
| 578 | else h->GetXaxis()->SetTitle("#eta");
|
---|
| 579 | h->GetYaxis()->SetTitle("#epsilon");
|
---|
| 580 | h->GetYaxis()->SetTitleSize(0.07);
|
---|
| 581 | h->GetXaxis()->SetTitleSize(0.07);
|
---|
| 582 | h->GetYaxis()->SetLabelSize(0.07);
|
---|
| 583 | h->GetXaxis()->SetLabelSize(0.07);
|
---|
| 584 | leg->SetLineStyle(0);
|
---|
| 585 | leg->SetFillStyle(0);
|
---|
| 586 | leg->SetLineWidth(0);
|
---|
| 587 | leg->SetLineColor(0);
|
---|
| 588 |
|
---|
| 589 | gStyle->SetOptTitle(0);
|
---|
| 590 | gStyle->SetOptStat(0);
|
---|
| 591 | gPad->SetBottomMargin(0.2);
|
---|
| 592 | gPad->SetLeftMargin(0.2);
|
---|
| 593 |
|
---|
| 594 | gPad->Modified();
|
---|
| 595 | gPad->Update();
|
---|
| 596 |
|
---|
| 597 | }
|
---|
| 598 |
|
---|
| 599 |
|
---|
| 600 | void Validation(const char *inputFile, const char *outputFile)
|
---|
| 601 | {
|
---|
| 602 | //gSystem->Load("libDelphes");
|
---|
| 603 |
|
---|
| 604 | std::cout << "input file : " << inputFile << " " << " , output file : " << outputFile << std::endl;
|
---|
| 605 |
|
---|
| 606 | TChain *chain = new TChain("Delphes");
|
---|
| 607 | chain->Add(inputFile);
|
---|
| 608 |
|
---|
| 609 | ExRootTreeReader *treeReader = new ExRootTreeReader(chain);
|
---|
| 610 | TClonesArray *branchParticle = treeReader->UseBranch("Particle");
|
---|
| 611 | TClonesArray *branchElectron = treeReader->UseBranch("Electron");
|
---|
| 612 | TClonesArray *branchMuon = treeReader->UseBranch("Muon");
|
---|
| 613 | TClonesArray *branchPhoton = treeReader->UseBranch("Photon");
|
---|
| 614 | TClonesArray *branchTrack = treeReader->UseBranch("Track");
|
---|
| 615 | TClonesArray *branchTower = treeReader->UseBranch("Tower");
|
---|
| 616 | TClonesArray *branchGenJet = treeReader->UseBranch("GenJet");
|
---|
| 617 | TClonesArray *branchPFJet = treeReader->UseBranch("Jet");
|
---|
| 618 | TClonesArray *branchCaloJet = treeReader->UseBranch("CaloJet");
|
---|
[6563d4d] | 619 | TClonesArray *branchScalarHT = treeReader->UseBranch("ScalarHT");
|
---|
| 620 | TClonesArray *branchMet = treeReader->UseBranch("MissingET");
|
---|
[cf88574] | 621 |
|
---|
| 622 | #ifdef ELECTRON
|
---|
| 623 |
|
---|
| 624 | ///////////////
|
---|
| 625 | // Electrons //
|
---|
| 626 | ///////////////
|
---|
| 627 |
|
---|
| 628 | // Reconstruction efficiency
|
---|
| 629 | TString elecs = "Electron";
|
---|
| 630 | int elID = 11;
|
---|
| 631 | std::pair<TH1D*,TH1D*> histos_el = GetEff<Electron>(branchElectron, branchParticle, "Electron", elID, treeReader);
|
---|
| 632 |
|
---|
| 633 | // tracking reconstruction efficiency
|
---|
| 634 | std::pair <TH1D*,TH1D*> histos_eltrack = GetEff<Track>(branchTrack, branchParticle, "electronTrack", elID, treeReader);
|
---|
| 635 |
|
---|
| 636 | // Tower reconstruction efficiency
|
---|
| 637 | std::pair <TH1D*,TH1D*> histos_eltower = GetEff<Tower>(branchTower, branchParticle, "electronTower", elID, treeReader);
|
---|
| 638 |
|
---|
| 639 | // Electron Energy Resolution
|
---|
| 640 | std::vector<resolPlot> plots_el;
|
---|
| 641 | HistogramsCollection(&plots_el, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electrons");
|
---|
| 642 | GetEres<Electron>( &plots_el, branchElectron, branchParticle, elID, treeReader);
|
---|
| 643 | TGraphErrors gr_el = EresGraph(&plots_el, true);
|
---|
[20ca0cd] | 644 | TGraphErrors gr_elFwd = EresGraph(&plots_el, false);
|
---|
[cf88574] | 645 | gr_el.SetName("Electron");
|
---|
[20ca0cd] | 646 | gr_elFwd.SetName("ElectronFwd");
|
---|
[cf88574] | 647 |
|
---|
| 648 | // Electron Track Energy Resolution
|
---|
| 649 | std::vector<resolPlot> plots_eltrack;
|
---|
| 650 | HistogramsCollection(&plots_eltrack, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electronsTracks");
|
---|
| 651 | GetEres<Track>( &plots_eltrack, branchTrack, branchParticle, elID, treeReader);
|
---|
| 652 | TGraphErrors gr_eltrack = EresGraph(&plots_eltrack, true);
|
---|
[20ca0cd] | 653 | TGraphErrors gr_eltrackFwd = EresGraph(&plots_eltrack, false);
|
---|
[cf88574] | 654 | gr_eltrack.SetName("ElectronTracks");
|
---|
[20ca0cd] | 655 | gr_eltrackFwd.SetName("ElectronTracksFwd");
|
---|
[cf88574] | 656 |
|
---|
| 657 | // Electron Tower Energy Resolution
|
---|
| 658 | std::vector<resolPlot> plots_eltower;
|
---|
| 659 | HistogramsCollection(&plots_eltower, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "electronsTower");
|
---|
| 660 | GetEres<Tower>( &plots_eltower, branchTower, branchParticle, elID, treeReader);
|
---|
| 661 | TGraphErrors gr_eltower = EresGraph(&plots_eltower, true);
|
---|
[20ca0cd] | 662 | TGraphErrors gr_eltowerFwd = EresGraph(&plots_eltower, false);
|
---|
[cf88574] | 663 | gr_eltower.SetName("ElectronTower");
|
---|
[20ca0cd] | 664 | gr_eltrackFwd.SetName("ElectronTracksFwd");
|
---|
[cf88574] | 665 |
|
---|
| 666 | // Canvases
|
---|
| 667 | TString elEff = "electronEff";
|
---|
| 668 | TCanvas *C_el1 = new TCanvas(elEff,elEff, 1000, 500);
|
---|
| 669 | C_el1->Divide(2);
|
---|
| 670 | C_el1->cd(1);
|
---|
| 671 | TLegend *leg_el1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
|
---|
| 672 |
|
---|
| 673 | gPad->SetLogx();
|
---|
| 674 | histos_eltrack.first->Draw();
|
---|
| 675 | addHist(histos_eltrack.first, leg_el1, 2);
|
---|
| 676 | histos_eltower.first->Draw("same");
|
---|
| 677 | addHist(histos_eltower.first, leg_el1, 3);
|
---|
| 678 | histos_el.first->Draw("same");
|
---|
| 679 | addHist(histos_el.first, leg_el1, 1);
|
---|
| 680 |
|
---|
| 681 | DrawAxis(histos_eltrack.first, leg_el1, 0);
|
---|
| 682 | leg_el1->Draw();
|
---|
| 683 |
|
---|
| 684 | C_el1->cd(2);
|
---|
| 685 | TLegend *leg_el2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
|
---|
| 686 |
|
---|
[20ca0cd] | 687 | gPad->SetLogx();
|
---|
[cf88574] | 688 | histos_eltrack.second->Draw();
|
---|
| 689 | addHist(histos_eltrack.second, leg_el2, 2);
|
---|
| 690 | histos_eltower.second->Draw("same");
|
---|
| 691 | addHist(histos_eltower.second, leg_el2, 3);
|
---|
| 692 | histos_el.second->Draw("same");
|
---|
| 693 | addHist(histos_el.second, leg_el2, 1);
|
---|
| 694 |
|
---|
[20ca0cd] | 695 | DrawAxis(histos_eltrack.second, leg_el2, 0);
|
---|
[cf88574] | 696 | delete(leg_el2);
|
---|
| 697 | C_el1->cd(0);
|
---|
| 698 |
|
---|
| 699 | TString elRes = "electronERes";
|
---|
[20ca0cd] | 700 | TString elResFwd = "electronEResForward";
|
---|
[cf88574] | 701 | TCanvas *C_el2 = new TCanvas(elRes,elRes, 1000, 500);
|
---|
[20ca0cd] | 702 | C_el2->Divide(2);
|
---|
| 703 | C_el2->cd(1);
|
---|
[cf88574] | 704 | TMultiGraph *mg_el = new TMultiGraph(elRes,elRes);
|
---|
| 705 | TLegend *leg_el = new TLegend(0.52,0.7,0.9,0.9);
|
---|
| 706 |
|
---|
| 707 | addGraph(mg_el, &gr_eltower, leg_el, 3);
|
---|
| 708 | addGraph(mg_el, &gr_eltrack, leg_el, 2);
|
---|
| 709 | addGraph(mg_el, &gr_el, leg_el, 1);
|
---|
| 710 |
|
---|
| 711 | mg_el->Draw("ACX");
|
---|
| 712 | leg_el->Draw();
|
---|
| 713 |
|
---|
| 714 | DrawAxis(mg_el, leg_el, 0.1);
|
---|
[20ca0cd] | 715 |
|
---|
| 716 | C_el2->cd(2);
|
---|
| 717 | TMultiGraph *mg_elFwd = new TMultiGraph(elResFwd,elResFwd);
|
---|
| 718 | TLegend *leg_elFwd = new TLegend(0.52,0.7,0.9,0.9);
|
---|
| 719 |
|
---|
| 720 | addGraph(mg_elFwd, &gr_eltowerFwd, leg_elFwd, 3);
|
---|
| 721 | addGraph(mg_elFwd, &gr_eltrackFwd, leg_elFwd, 2);
|
---|
| 722 | addGraph(mg_elFwd, &gr_elFwd, leg_elFwd, 1);
|
---|
| 723 |
|
---|
| 724 | mg_elFwd->Draw("ACX");
|
---|
| 725 |
|
---|
| 726 | DrawAxis(mg_elFwd, leg_elFwd, 0.1);
|
---|
| 727 |
|
---|
[cf88574] | 728 |
|
---|
| 729 | C_el1->SaveAs(elEff+".eps");
|
---|
| 730 | C_el2->SaveAs(elRes+".eps");
|
---|
| 731 |
|
---|
| 732 | #endif
|
---|
| 733 |
|
---|
| 734 | gDirectory->cd(0);
|
---|
| 735 |
|
---|
| 736 | #ifdef MUON
|
---|
| 737 |
|
---|
| 738 | ///////////
|
---|
| 739 | // Muons //
|
---|
| 740 | ///////////
|
---|
| 741 |
|
---|
| 742 | // Reconstruction efficiency
|
---|
| 743 | int muID = 13;
|
---|
| 744 | std::pair<TH1D*,TH1D*> histos_mu = GetEff<Muon>(branchMuon, branchParticle,"Muon", muID, treeReader);
|
---|
| 745 |
|
---|
| 746 | // muon tracking reconstruction efficiency
|
---|
| 747 | std::pair <TH1D*,TH1D*> histos_mutrack = GetEff<Track>(branchTrack, branchParticle, "muonTrack", muID, treeReader);
|
---|
| 748 |
|
---|
| 749 | // Muon Energy Resolution
|
---|
| 750 | std::vector<resolPlot> plots_mu;
|
---|
| 751 | HistogramsCollection(&plots_mu, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "muons");
|
---|
| 752 | GetEres<Muon>( &plots_mu, branchMuon, branchParticle, muID, treeReader);
|
---|
| 753 | TGraphErrors gr_mu = EresGraph(&plots_mu, true);
|
---|
[20ca0cd] | 754 | TGraphErrors gr_muFwd = EresGraph(&plots_mu, false);
|
---|
[cf88574] | 755 | gr_mu.SetName("Muon");
|
---|
[20ca0cd] | 756 | gr_muFwd.SetName("MuonFwd");
|
---|
[cf88574] | 757 |
|
---|
| 758 | // Muon Track Energy Resolution
|
---|
| 759 | std::vector<resolPlot> plots_mutrack;
|
---|
| 760 | HistogramsCollection(&plots_mutrack, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "muonsTracks");
|
---|
| 761 | GetEres<Track>( &plots_mutrack, branchTrack, branchParticle, muID, treeReader);
|
---|
| 762 | TGraphErrors gr_mutrack = EresGraph(&plots_mutrack, true);
|
---|
[20ca0cd] | 763 | TGraphErrors gr_mutrackFwd = EresGraph(&plots_mutrack, false);
|
---|
| 764 | gr_mutrackFwd.SetName("MuonTracksFwd");
|
---|
[cf88574] | 765 |
|
---|
| 766 | // Canvas
|
---|
| 767 |
|
---|
| 768 | TString muEff = "muonEff";
|
---|
| 769 | TCanvas *C_mu1 = new TCanvas(muEff,muEff, 1000, 500);
|
---|
| 770 | C_mu1->Divide(2);
|
---|
| 771 | C_mu1->cd(1);
|
---|
| 772 | TLegend *leg_mu1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
|
---|
| 773 |
|
---|
| 774 | gPad->SetLogx();
|
---|
| 775 | histos_mutrack.first->Draw();
|
---|
| 776 | addHist(histos_mutrack.first, leg_mu1, 2);
|
---|
| 777 | histos_mu.first->Draw("same");
|
---|
| 778 | addHist(histos_mu.first, leg_mu1, 1);
|
---|
| 779 |
|
---|
| 780 | DrawAxis(histos_mutrack.first, leg_mu1, 0);
|
---|
| 781 | leg_mu1->Draw();
|
---|
| 782 |
|
---|
| 783 | C_mu1->cd(2);
|
---|
| 784 | TLegend *leg_mu2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
|
---|
| 785 |
|
---|
[20ca0cd] | 786 | gPad->SetLogx();
|
---|
[cf88574] | 787 | histos_mutrack.second->Draw();
|
---|
| 788 | addHist(histos_mutrack.second, leg_mu2, 2);
|
---|
| 789 | histos_mu.second->Draw("same");
|
---|
| 790 | addHist(histos_mu.second, leg_mu2, 1);
|
---|
| 791 |
|
---|
| 792 | DrawAxis(histos_mutrack.second, leg_mu2, 1);
|
---|
| 793 |
|
---|
| 794 | TString muRes = "muonERes";
|
---|
[20ca0cd] | 795 | TString muResFwd = "muonEResFwd";
|
---|
| 796 |
|
---|
[cf88574] | 797 | TCanvas *C_mu = new TCanvas(muRes,muRes, 1000, 500);
|
---|
[20ca0cd] | 798 | C_mu->Divide(2);
|
---|
| 799 | C_mu->cd(1);
|
---|
[cf88574] | 800 | TMultiGraph *mg_mu = new TMultiGraph(muRes,muRes);
|
---|
| 801 | TLegend *leg_mu = new TLegend(0.52,0.7,0.9,0.9);
|
---|
| 802 |
|
---|
| 803 | addGraph(mg_mu, &gr_mutrack, leg_mu, 2);
|
---|
| 804 | addGraph(mg_mu, &gr_mu, leg_mu, 1);
|
---|
| 805 |
|
---|
| 806 | mg_mu->Draw("ACX");
|
---|
| 807 | leg_mu->Draw();
|
---|
| 808 |
|
---|
| 809 | DrawAxis(mg_mu, leg_mu, 0.3);
|
---|
| 810 |
|
---|
[20ca0cd] | 811 | C_mu->cd(2);
|
---|
| 812 | TMultiGraph *mg_muFwd = new TMultiGraph(muResFwd,muResFwd);
|
---|
| 813 | TLegend *leg_muFwd = new TLegend(0.52,0.7,0.9,0.9);
|
---|
| 814 |
|
---|
| 815 | addGraph(mg_muFwd, &gr_mutrackFwd, leg_muFwd, 2);
|
---|
| 816 | addGraph(mg_muFwd, &gr_muFwd, leg_muFwd, 1);
|
---|
| 817 |
|
---|
| 818 | mg_muFwd->Draw("ACX");
|
---|
| 819 |
|
---|
| 820 | DrawAxis(mg_muFwd, leg_muFwd, 0.3);
|
---|
| 821 |
|
---|
[cf88574] | 822 | C_mu1->SaveAs(muEff+".eps");
|
---|
| 823 | C_mu->SaveAs(muRes+".eps");
|
---|
| 824 |
|
---|
| 825 | #endif
|
---|
| 826 |
|
---|
| 827 | gDirectory->cd(0);
|
---|
| 828 |
|
---|
| 829 | #ifdef PHOTON
|
---|
| 830 |
|
---|
| 831 | /////////////
|
---|
| 832 | // Photons //
|
---|
| 833 | /////////////
|
---|
| 834 |
|
---|
| 835 | // Reconstruction efficiency
|
---|
| 836 | int phID = 22;
|
---|
| 837 | std::pair<TH1D*,TH1D*> histos_ph = GetEff<Electron>(branchPhoton, branchParticle, "Photon", phID, treeReader);
|
---|
| 838 | std::pair<TH1D*,TH1D*> histos_phtower = GetEff<Electron>(branchTower, branchParticle, "Photon", phID, treeReader);
|
---|
| 839 |
|
---|
| 840 | // Photon Energy Resolution
|
---|
| 841 | std::vector<resolPlot> plots_ph;
|
---|
| 842 | HistogramsCollection(&plots_ph, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "photons");
|
---|
| 843 | GetEres<Photon>( &plots_ph, branchPhoton, branchParticle, phID, treeReader);
|
---|
| 844 | TGraphErrors gr_ph = EresGraph(&plots_ph, true);
|
---|
[20ca0cd] | 845 | TGraphErrors gr_phFwd = EresGraph(&plots_ph, false);
|
---|
[cf88574] | 846 | gr_ph.SetName("Photon");
|
---|
[20ca0cd] | 847 | gr_phFwd.SetName("PhotonFwd");
|
---|
| 848 |
|
---|
[cf88574] | 849 |
|
---|
| 850 | // Photon Tower Energy Resolution
|
---|
| 851 | std::vector<resolPlot> plots_phtower;
|
---|
| 852 | HistogramsCollection(&plots_phtower, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "photonsTower");
|
---|
| 853 | GetEres<Tower>( &plots_phtower, branchTower, branchParticle, phID, treeReader);
|
---|
| 854 | TGraphErrors gr_phtower = EresGraph(&plots_phtower, true);
|
---|
[20ca0cd] | 855 | TGraphErrors gr_phtowerFwd = EresGraph(&plots_phtower, false);
|
---|
[cf88574] | 856 | gr_phtower.SetName("PhotonTower");
|
---|
[20ca0cd] | 857 | gr_phtowerFwd.SetName("PhotonTowerFwd");
|
---|
[cf88574] | 858 |
|
---|
| 859 | // Canvas
|
---|
| 860 |
|
---|
| 861 | TString phEff = "photonEff";
|
---|
| 862 | TCanvas *C_ph1 = new TCanvas(phEff,phEff, 1000, 500);
|
---|
| 863 | C_ph1->Divide(2);
|
---|
| 864 | C_ph1->cd(1);
|
---|
| 865 | TLegend *leg_ph1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
|
---|
| 866 |
|
---|
| 867 | gPad->SetLogx();
|
---|
| 868 | histos_phtower.first->Draw();
|
---|
| 869 | addHist(histos_phtower.first, leg_ph1, 3);
|
---|
| 870 | histos_ph.first->Draw("same");
|
---|
| 871 | addHist(histos_ph.first, leg_ph1, 1);
|
---|
| 872 |
|
---|
| 873 | DrawAxis(histos_phtower.first, leg_ph1, 0);
|
---|
| 874 | leg_ph1->Draw();
|
---|
| 875 |
|
---|
| 876 | C_ph1->cd(2);
|
---|
| 877 | TLegend *leg_ph2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
|
---|
| 878 |
|
---|
[20ca0cd] | 879 | gPad->SetLogx();
|
---|
[cf88574] | 880 | histos_phtower.second->Draw("same");
|
---|
| 881 | addHist(histos_phtower.second, leg_ph2, 3);
|
---|
| 882 | histos_ph.second->Draw("same");
|
---|
| 883 | addHist(histos_ph.second, leg_ph2, 1);
|
---|
| 884 |
|
---|
| 885 | DrawAxis(histos_phtower.second, leg_ph2, 1);
|
---|
| 886 |
|
---|
| 887 | C_ph1->SaveAs(phEff+".eps");
|
---|
| 888 |
|
---|
| 889 | TString phRes = "phERes";
|
---|
[20ca0cd] | 890 | TString phResFwd = "phEResFwd";
|
---|
| 891 |
|
---|
[cf88574] | 892 | TCanvas *C_ph = new TCanvas(phRes,phRes, 1000, 500);
|
---|
[20ca0cd] | 893 | C_ph->Divide(2);
|
---|
| 894 | C_ph->cd(1);
|
---|
[cf88574] | 895 | TMultiGraph *mg_ph = new TMultiGraph(phRes,phRes);
|
---|
| 896 | TLegend *leg_ph = new TLegend(0.52,0.7,0.9,0.9);
|
---|
| 897 |
|
---|
| 898 | addGraph(mg_ph, &gr_phtower, leg_ph, 3);
|
---|
| 899 | addGraph(mg_ph, &gr_ph, leg_ph, 1);
|
---|
| 900 |
|
---|
| 901 | mg_ph->Draw("ACX");
|
---|
| 902 | leg_ph->Draw();
|
---|
| 903 |
|
---|
| 904 | DrawAxis(mg_ph, leg_ph, 0.3);
|
---|
| 905 |
|
---|
[20ca0cd] | 906 | C_ph->cd(2);
|
---|
| 907 | TMultiGraph *mg_phFwd = new TMultiGraph(phResFwd,phResFwd);
|
---|
| 908 | TLegend *leg_phFwd = new TLegend(0.52,0.7,0.9,0.9);
|
---|
| 909 |
|
---|
| 910 | addGraph(mg_phFwd, &gr_phtowerFwd, leg_phFwd, 3);
|
---|
| 911 | addGraph(mg_phFwd, &gr_phFwd, leg_phFwd, 1);
|
---|
| 912 |
|
---|
| 913 | mg_phFwd->Draw("ACX");
|
---|
| 914 |
|
---|
| 915 | DrawAxis(mg_phFwd, leg_phFwd, 0.3);
|
---|
| 916 |
|
---|
[cf88574] | 917 | C_ph->SaveAs(phRes+".eps");
|
---|
| 918 |
|
---|
| 919 | #endif
|
---|
| 920 |
|
---|
| 921 | gDirectory->cd(0);
|
---|
| 922 |
|
---|
| 923 | #ifdef JET
|
---|
| 924 |
|
---|
| 925 | //////////
|
---|
| 926 | // Jets //
|
---|
| 927 | //////////
|
---|
| 928 |
|
---|
| 929 | // BJets Reconstruction efficiency
|
---|
| 930 | int bID = 5;
|
---|
| 931 | std::pair<TH1D*,TH1D*> histos_btag = GetEff<Jet>(branchPFJet, branchParticle,"BTag", bID, treeReader);
|
---|
| 932 |
|
---|
| 933 | // TauJets Reconstruction efficiency
|
---|
| 934 | int tauID = 15;
|
---|
| 935 | std::pair<TH1D*,TH1D*> histos_tautag = GetEff<Jet>(branchPFJet, branchParticle,"TauTag", tauID, treeReader);
|
---|
| 936 |
|
---|
| 937 | // PFJets Energy Resolution
|
---|
| 938 | std::vector<resolPlot> plots_pfjets;
|
---|
| 939 | HistogramsCollection(&plots_pfjets, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "PFJet");
|
---|
| 940 | GetJetsEres( &plots_pfjets, branchPFJet, branchGenJet, treeReader);
|
---|
| 941 | TGraphErrors gr_pfjets = EresGraph(&plots_pfjets, true);
|
---|
[20ca0cd] | 942 | TGraphErrors gr_pfjetsFwd = EresGraph(&plots_pfjets, false);
|
---|
[cf88574] | 943 | gr_pfjets.SetName("pfJet");
|
---|
[20ca0cd] | 944 | gr_pfjetsFwd.SetName("pfJetFwd");
|
---|
[cf88574] | 945 |
|
---|
[6563d4d] | 946 | // CaloJets Energy Resolution
|
---|
[cf88574] | 947 | std::vector<resolPlot> plots_calojets;
|
---|
| 948 | HistogramsCollection(&plots_calojets, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "CaloJet");
|
---|
| 949 | GetJetsEres( &plots_calojets, branchCaloJet, branchGenJet, treeReader);
|
---|
| 950 | TGraphErrors gr_calojets = EresGraph(&plots_calojets, true);
|
---|
[20ca0cd] | 951 | TGraphErrors gr_calojetsFwd = EresGraph(&plots_calojets, false);
|
---|
[cf88574] | 952 | gr_calojets.SetName("caloJet");
|
---|
[20ca0cd] | 953 | gr_calojetsFwd.SetName("caloJetFwd");
|
---|
[cf88574] | 954 |
|
---|
[6563d4d] | 955 | // MET Resolution vs HT
|
---|
| 956 | std::vector<resolPlot> plots_met;
|
---|
| 957 | HistogramsCollection(&plots_met, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "MET", -200, 200);
|
---|
| 958 | GetMetres( &plots_met, branchScalarHT, branchMet, branchPFJet, treeReader);
|
---|
| 959 | TGraphErrors gr_met = EresGraph(&plots_met, true);
|
---|
| 960 | gr_calojets.SetName("MET");
|
---|
| 961 |
|
---|
[cf88574] | 962 | // Canvas
|
---|
| 963 | TString btagEff = "btagEff";
|
---|
| 964 | TCanvas *C_btag1 = new TCanvas(btagEff,btagEff, 1000, 500);
|
---|
| 965 | C_btag1->Divide(2);
|
---|
| 966 | C_btag1->cd(1);
|
---|
| 967 | TLegend *leg_btag1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
|
---|
| 968 |
|
---|
| 969 | gPad->SetLogx();
|
---|
| 970 | histos_btag.first->Draw();
|
---|
| 971 | addHist(histos_btag.first, leg_btag1, 1);
|
---|
| 972 |
|
---|
| 973 | DrawAxis(histos_btag.first, leg_btag1, 0);
|
---|
| 974 | leg_btag1->Draw();
|
---|
| 975 |
|
---|
| 976 | C_btag1->cd(2);
|
---|
| 977 | TLegend *leg_btag2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
|
---|
| 978 |
|
---|
| 979 | histos_btag.second->Draw();
|
---|
| 980 | addHist(histos_btag.second, leg_btag2, 1);
|
---|
| 981 |
|
---|
| 982 | DrawAxis(histos_btag.second, leg_btag2, 1);
|
---|
| 983 | delete(leg_btag2);
|
---|
| 984 | C_btag1->cd(0);
|
---|
| 985 |
|
---|
| 986 | TString tautagEff = "tautagEff";
|
---|
| 987 | TCanvas *C_tautag1 = new TCanvas(tautagEff,tautagEff, 1000, 500);
|
---|
| 988 | C_tautag1->Divide(2);
|
---|
| 989 | C_tautag1->cd(1);
|
---|
| 990 | TLegend *leg_tautag1 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
|
---|
| 991 |
|
---|
| 992 | gPad->SetLogx();
|
---|
| 993 | histos_tautag.first->Draw();
|
---|
| 994 | addHist(histos_tautag.first, leg_tautag1, 1);
|
---|
| 995 |
|
---|
| 996 | DrawAxis(histos_tautag.first, leg_tautag1, 0);
|
---|
| 997 | leg_tautag1->Draw();
|
---|
| 998 |
|
---|
| 999 | C_tautag1->cd(2);
|
---|
| 1000 | TLegend *leg_tautag2 = new TLegend(effLegXmin,effLegYmin,effLegXmax,effLegYmax);
|
---|
| 1001 |
|
---|
| 1002 | histos_tautag.second->Draw();
|
---|
| 1003 | addHist(histos_tautag.second, leg_tautag2, 1);
|
---|
| 1004 |
|
---|
| 1005 | DrawAxis(histos_tautag.second, leg_tautag2, 1);
|
---|
| 1006 | delete(leg_tautag2);
|
---|
| 1007 | C_tautag1->cd(0);
|
---|
| 1008 |
|
---|
| 1009 | TString jetRes = "jetERes";
|
---|
[20ca0cd] | 1010 | TString jetResFwd = "jetEResFwd";
|
---|
[cf88574] | 1011 | TCanvas *C_jet = new TCanvas(jetRes,jetRes, 1000, 500);
|
---|
[20ca0cd] | 1012 | C_jet->Divide(2);
|
---|
| 1013 |
|
---|
| 1014 | C_jet->cd(1);
|
---|
[cf88574] | 1015 | TMultiGraph *mg_jet = new TMultiGraph(jetRes,jetRes);
|
---|
| 1016 | TLegend *leg_jet = new TLegend(0.52,0.7,0.9,0.9);
|
---|
| 1017 |
|
---|
| 1018 | addGraph(mg_jet, &gr_calojets, leg_jet, 3);
|
---|
| 1019 | addGraph(mg_jet, &gr_pfjets, leg_jet, 1);
|
---|
| 1020 |
|
---|
| 1021 | mg_jet->Draw("ACX");
|
---|
| 1022 | leg_jet->Draw();
|
---|
| 1023 |
|
---|
| 1024 | DrawAxis(mg_jet, leg_jet, 0.25);
|
---|
| 1025 |
|
---|
[20ca0cd] | 1026 | C_jet->cd(2);
|
---|
| 1027 | TMultiGraph *mg_jetFwd = new TMultiGraph(jetResFwd,jetResFwd);
|
---|
| 1028 | TLegend *leg_jetFwd = new TLegend(0.52,0.7,0.9,0.9);
|
---|
| 1029 |
|
---|
| 1030 | addGraph(mg_jetFwd, &gr_calojetsFwd, leg_jetFwd, 3);
|
---|
| 1031 | addGraph(mg_jetFwd, &gr_pfjetsFwd, leg_jetFwd, 1);
|
---|
| 1032 |
|
---|
| 1033 | mg_jetFwd->Draw("ACX");
|
---|
| 1034 |
|
---|
| 1035 | DrawAxis(mg_jetFwd, leg_jetFwd, 0.25);
|
---|
| 1036 |
|
---|
[cf88574] | 1037 | C_btag1->SaveAs(btagEff+".eps");
|
---|
| 1038 | C_jet->SaveAs(jetRes+".eps");
|
---|
| 1039 |
|
---|
[6563d4d] | 1040 | TString metRes = "MetRes";
|
---|
| 1041 | TCanvas *C_met = new TCanvas(metRes,metRes, 1000, 500);
|
---|
| 1042 |
|
---|
| 1043 | TMultiGraph *mg_met = new TMultiGraph(metRes,metRes);
|
---|
| 1044 | TLegend *leg_met = new TLegend(0.52,0.7,0.9,0.9);
|
---|
| 1045 |
|
---|
| 1046 | addGraph(mg_met, &gr_met, leg_met, 3);
|
---|
| 1047 |
|
---|
| 1048 | mg_met->Draw("ACX");
|
---|
| 1049 | leg_met->Draw();
|
---|
| 1050 |
|
---|
| 1051 | DrawAxis(mg_met, leg_met, 100);
|
---|
| 1052 |
|
---|
| 1053 | C_met->SaveAs(metRes+".eps");
|
---|
| 1054 |
|
---|
| 1055 |
|
---|
[cf88574] | 1056 | #endif
|
---|
| 1057 |
|
---|
| 1058 | /*
|
---|
| 1059 | // CaloJets Energy Resolution
|
---|
| 1060 | std::vector<resolPlot> plots_calojets;
|
---|
| 1061 | HistogramsCollection(&plots_calojets, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax), "caloJet");
|
---|
| 1062 | GetJetsEres( &plots_calojets, branchCaloJet, branchGenJet, treeReader);
|
---|
| 1063 | TGraphErrors gr_calojets = EresGraph(&plots_calojets, true);
|
---|
| 1064 | gr_calojets.SetName("caloJet");
|
---|
| 1065 | */
|
---|
| 1066 |
|
---|
| 1067 | TFile *fout = new TFile(outputFile,"recreate");
|
---|
| 1068 |
|
---|
| 1069 | #ifdef ELECTRON
|
---|
| 1070 | for (int bin = 0; bin < Nbins; bin++)
|
---|
| 1071 | {
|
---|
| 1072 | plots_el.at(bin).cenResolHist->Write();
|
---|
| 1073 | plots_eltrack.at(bin).cenResolHist->Write();
|
---|
| 1074 | plots_eltower.at(bin).cenResolHist->Write();
|
---|
[20ca0cd] | 1075 | plots_el.at(bin).fwdResolHist->Write();
|
---|
| 1076 | plots_eltrack.at(bin).fwdResolHist->Write();
|
---|
| 1077 | plots_eltower.at(bin).fwdResolHist->Write();
|
---|
[cf88574] | 1078 | }
|
---|
| 1079 | // gr.Write();
|
---|
| 1080 | histos_el.first->Write();
|
---|
| 1081 | //histos_el.second->Write();
|
---|
| 1082 | histos_eltrack.first->Write();
|
---|
| 1083 | //histos_eltrack.second->Write();
|
---|
| 1084 | histos_eltower.first->Write();
|
---|
| 1085 | C_el1->Write();
|
---|
| 1086 | C_el2->Write();
|
---|
| 1087 | #endif
|
---|
| 1088 |
|
---|
| 1089 | #ifdef MUON
|
---|
| 1090 | histos_mu.first->Write();
|
---|
| 1091 | histos_mu.second->Write();
|
---|
| 1092 | histos_mutrack.first->Write();
|
---|
| 1093 | histos_mutrack.second->Write();
|
---|
| 1094 | C_mu1->Write();
|
---|
| 1095 | C_mu->Write();
|
---|
| 1096 | #endif
|
---|
| 1097 |
|
---|
| 1098 | #ifdef PHOTON
|
---|
| 1099 | histos_ph.first->Write();
|
---|
| 1100 | histos_ph.second->Write();
|
---|
| 1101 | C_ph1->Write();
|
---|
| 1102 | C_ph->Write();
|
---|
| 1103 | #endif
|
---|
| 1104 |
|
---|
| 1105 | #ifdef JET
|
---|
| 1106 | for (int bin = 0; bin < Nbins; bin++)
|
---|
| 1107 | {
|
---|
| 1108 | plots_pfjets.at(bin).cenResolHist->Write();
|
---|
[20ca0cd] | 1109 | plots_pfjets.at(bin).fwdResolHist->Write();
|
---|
[cf88574] | 1110 | plots_calojets.at(bin).cenResolHist->Write();
|
---|
[20ca0cd] | 1111 | plots_calojets.at(bin).fwdResolHist->Write();
|
---|
[6563d4d] | 1112 | plots_met.at(bin).cenResolHist->Write();
|
---|
[cf88574] | 1113 | }
|
---|
| 1114 | histos_btag.first->Write();
|
---|
| 1115 | histos_btag.second->Write();
|
---|
| 1116 | histos_tautag.first->Write();
|
---|
| 1117 | histos_tautag.second->Write();
|
---|
| 1118 | C_btag1->Write();
|
---|
| 1119 | C_tautag1->Write();
|
---|
| 1120 | C_jet->Write();
|
---|
[6563d4d] | 1121 | C_met->Write();
|
---|
[cf88574] | 1122 | #endif
|
---|
| 1123 |
|
---|
| 1124 | fout->Write();
|
---|
| 1125 |
|
---|
| 1126 | cout << "** Exiting..." << endl;
|
---|
| 1127 |
|
---|
| 1128 | //delete plots;
|
---|
| 1129 | //delete result;
|
---|
| 1130 | delete treeReader;
|
---|
| 1131 | delete chain;
|
---|
| 1132 | }
|
---|
| 1133 |
|
---|
| 1134 | //------------------------------------------------------------------------------
|
---|