DESY Hbb Analysis Framework
Functions
AnalysisJetsExample.cc File Reference
#include "boost/program_options.hpp"
#include "boost/algorithm/string.hpp"
#include <string>
#include <iostream>
#include <vector>
#include "TFile.h"
#include "TFileCollection.h"
#include "TChain.h"
#include "TH1.h"
#include "Analysis/Tools/interface/Analysis.h"
#include "Analysis/Tools/bin/macro_config.h"

Go to the source code of this file.

Functions

float btagMin (const std::string &btagwp)
 
float btagMin (const string &wp)
 
int main (int argc, char *argv[])
 

Function Documentation

float btagMin ( const std::string &  btagwp)
float btagMin ( const string &  wp)

Definition at line 120 of file AnalysisJetsExample.cc.

References btagwploose_, btagwpmedium_, and btagwptight_.

121 {
122  float min = -1000;
123  if ( wp == "loose" ) min = btagwploose_;
124  if ( wp == "medium" ) min = btagwpmedium_;
125  if ( wp == "tight" ) min = btagwptight_;
126 
127  return min;
128 
129 }
float btagwploose_
Definition: macro_config.h:108
float btagwpmedium_
Definition: macro_config.h:109
float btagwptight_
Definition: macro_config.h:110
int main ( int  argc,
char *  argv[] 
)

Definition at line 22 of file AnalysisJetsExample.cc.

References analysis::tools::Analysis::addTree(), analysis::tools::Jet::btag(), btagalgo_, analysis::tools::Analysis::btagCalibration(), analysis::tools::Jet::btagSF(), btagsf_, analysis::tools::Jet::btagSFdown(), analysis::tools::Jet::btagSFup(), btagwp_, analysis::tools::Analysis::collection(), analysis::tools::Candidate::eta(), analysis::tools::Analysis::event(), analysis::tools::Jet::extendedFlavour(), analysis::tools::Jet::flavour(), genjetsCol_, genParticleCol_, inputlist_, analysis::tools::Jet::jerCorrection(), analysis::tools::Jet::jerInfo(), analysis::tools::Jet::jerMatch(), jerpt_, analysis::tools::Jet::jerSF(), jersf_, analysis::tools::Analysis::jetResolutionInfo(), jetsCol_, analysis::tools::Analysis::lumiSection(), macro_config(), nevtmax_, analysis::tools::Candidate::phi(), analysis::tools::Jet::pileupJetIdFullDiscriminant(), analysis::tools::Jet::pileupJetIdFullId(), analysis::tools::Candidate::pt(), analysis::tools::Jet::qgLikelihood(), analysis::tools::Analysis::run(), and analysis::tools::Analysis::size().

23 {
24  if ( macro_config(argc, argv) != 0 ) return -1;
25 
26  TH1::SetDefaultSumw2(); // proper treatment of errors when scaling histograms
27 
28  // Input files list
30 
31  // btag
32 // float btagmin = btagMin(btagwp_);
33  // b-tag scale factors
34  // inputs from the config file test/analysis_bjetsv2.cfg
35  // btagalgo_ = "deepcsv";
36  // btagsf_ = "../data/DeepCSV_94XSF_V3_B_F.csv";
37  // btagwp_ = "tight";
38  auto bsf_reader = analysis.btagCalibration(btagalgo_, btagsf_, btagwp_);
39 
40  // jer
41  // Jet energy resolution scale factors and pt resolution
42  auto jerinfo = analysis.jetResolutionInfo(jerpt_,jersf_);
43 
44  // Physics Objects Collections
45  analysis.addTree<Jet> ("Jets",jetsCol_);
46  analysis.addTree<GenParticle> ("GenParticles",genParticleCol_);
47  analysis.addTree<GenJet> ("GenJets",genjetsCol_);
48 
49  // Analysis of events
50  std::cout << "This analysis has " << analysis.size() << " events" << std::endl;
51  int nevts = analysis.size();
52  if ( nevts > 0 ) nevts = nevtmax_;
53  for ( int i = 0 ; i < nevts ; ++i )
54  {
55  analysis.event(i);
56 
57  std::cout << "++++++ ENTRY " << i;
58 
59  // EventInfo
60  std::cout << ", Run = " << analysis.run();
61  std::cout << ", Event = " << analysis.event();
62  std::cout << ", LumiSection = " << analysis.lumiSection();
63  std::cout << "\n" << std::endl;
64 
65  // GenParticles
66  auto particles = analysis.collection<GenParticle>("GenParticles"); // of course one can also loop over the particles
67  // GenJets
68  auto genjets = analysis.collection<GenJet>("GenJets"); // of course one can also loop over the genjets
69  // Jets
70  auto jets = analysis.collection<Jet>("Jets");
71 
72  // associate partons to jets
73  jets->associatePartons(particles,0.4,1);
74  // associate genjets to jets
75  jets->addGenJets(genjets); // if not defined -> jerMatch = false, smearing will be done using the stochastic method only
76 
77  for ( int j = 0 ; j < jets->size() ; ++j )
78  {
79  Jet jet = jets->at(j);
80 
81  // BTAGGING
82  float btag = -10000.;
83  if ( btagalgo_ == "deepcsv" )
84  btag = jet.btag("btag_deepb") + jet.btag("btag_deepbb");
85  if ( btagalgo_ == "deepflavour" )
86  btag = jet.btag("btag_dfb") + jet.btag("btag_dfbb") + jet.btag("btag_dflepb");
87 
88  // b-tag scale factors central, up and down
89  double jet_bscalefactor = jet.btagSF(bsf_reader); // OR jet.btagSF(analysis.btagCalibration());
90  double jet_bscalefactorup = jet.btagSFup(bsf_reader,2);
91  double jet_bscalefactordown = jet.btagSFdown(bsf_reader,2);
92 
93  // JER
94  jet.jerInfo(*jerinfo,0.2); // this also performs matching to the added gen jets above, with delta R < 0.2 which is default and can be omitted
95 
96 
97  // PRINTOUT
98 
99  std::cout << "Jet #" << j << ": ";
100  std::cout << "pT = " << jet.pt() << ", ";
101  std::cout << "eta = " << jet.eta() << ", ";
102  std::cout << "phi = " << jet.phi() << ", ";
103  std::cout << "flavour = " << jet.flavour() << " (extended = " << jet.extendedFlavour() << "), ";
104  std::cout << "btag = " << btag << " with scale factor = " << jet_bscalefactor;
105  std::cout << " up = " << jet_bscalefactorup << " down = " << jet_bscalefactordown << std::endl;
106  std::cout << " JER SF = " << jet.jerSF() << ", match = " << jet.jerMatch() << " jer corr = " << jet.jerCorrection() << " + ";
107  std::cout << jet.jerCorrection("up") << " - " << jet.jerCorrection("down") << std::endl;
108  std::cout << " quark-gluon likelihood = " << jet.qgLikelihood() << std::endl;
109  std::cout << " pileup jet id full discriminant = " << jet.pileupJetIdFullDiscriminant() << std::endl;
110  std::cout << " pileup jet id full id = " << jet.pileupJetIdFullId() << std::endl;
111 
112  }
113  std::cout << "===================" << std::endl;
114  }
115  std::cout << "Btag WP = " << btagwp_ << std::endl;
116 
117 //
118 }
float eta() const
returns the pseudorapidity
Definition: Candidate.cc:134
float jerSF() const
returns jet energy resolution SF
Definition: Jet.cc:72
double btagSFdown(std::shared_ptr< BTagCalibrationReader > reader, const float &nsig=1, const std::string &flavalgo="Hadron") const
Definition: Jet.cc:320
std::string jetsCol_
Definition: macro_config.h:132
int pileupJetIdFullId() const
Definition: Jet.cc:233
std::string jersf_
Definition: macro_config.h:71
std::string extendedFlavour() const
returns the extended flavour definition
Definition: Jet.cc:66
int macro_config(int argc, char *argv[])
Definition: macro_config.h:145
float jerCorrection(const std::string &var="nominal", const float &nsig=1) const
Definition: Jet.cc:179
float phi() const
returns the azimuthal angle
Definition: Candidate.cc:135
int flavour() const
returns the flavour with the Hadron definition (=0 for data)
Definition: Jet.cc:58
std::string jerpt_
Definition: macro_config.h:70
std::string inputlist_
Definition: macro_config.h:21
float pt() const
returns the transverse momentum
Definition: Candidate.cc:133
double btagSF(std::shared_ptr< BTagCalibrationReader > reader, const std::string &flavalgo="Hadron") const
btag SF
Definition: Jet.cc:307
void jerInfo(const JetResolutionInfo &, const std::string &)
Definition: Jet.cc:199
std::string btagwp_
Definition: macro_config.h:107
std::string genParticleCol_
Definition: macro_config.h:137
std::string genjetsCol_
Definition: macro_config.h:138
int nevtmax_
Definition: macro_config.h:15
float btag(const std::string &) const
returns the btag value of btag_csvivf
Definition: Jet.cc:57
std::string btagalgo_
Definition: macro_config.h:106
double btagSFup(std::shared_ptr< BTagCalibrationReader > reader, const float &nsig=1, const std::string &flavalgo="Hadron") const
Definition: Jet.cc:311
float qgLikelihood() const
quark-gluon separation
Definition: Jet.cc:231
float pileupJetIdFullDiscriminant() const
pile up jet id
Definition: Jet.cc:232
bool jerMatch(const std::string &)
JER matching.
Definition: Jet.cc:101
std::string btagsf_
Definition: macro_config.h:67