CMS 3D CMS Logo

AnalysisTasksAnalyzerBTag.cc
Go to the documentation of this file.
4 
5 
9  Jets_(cfg.getParameter<edm::InputTag>("Jets")),
10  bTagAlgo_(cfg.getParameter<std::string>("bTagAlgo")),
11  bins_(cfg.getParameter<unsigned int>("bins")),
12  lowerbin_(cfg.getParameter<double>("lowerbin")),
13  upperbin_(cfg.getParameter<double>("upperbin"))
14 {
15  hists_["BTag_b"] = fs.make<TH1F>("BTag_b" , "BTag_b" , bins_, lowerbin_, upperbin_);
16  hists_["BTag_g"] = fs.make<TH1F>("BTag_g" , "BTag_g" , bins_, lowerbin_, upperbin_);
17  hists_["BTag_c"] = fs.make<TH1F>("BTag_c" , "BTag_c" , bins_, lowerbin_, upperbin_);
18  hists_["BTag_uds"] = fs.make<TH1F>("BTag_uds", "BTag_uds", bins_, lowerbin_, upperbin_);
19  hists_["BTag_other"] = fs.make<TH1F>("BTag_other", "BTag_other", bins_, lowerbin_, upperbin_);
20  hists_["effBTag_b"] = fs.make<TH1F>("effBTag_b" , "effBTag_b" , bins_, lowerbin_, upperbin_);
21  hists_["effBTag_g"] = fs.make<TH1F>("effBTag_g" , "effBTag_g" , bins_, lowerbin_, upperbin_);
22  hists_["effBTag_c"] = fs.make<TH1F>("effBTag_c" , "effBTag_c" , bins_, lowerbin_, upperbin_);
23  hists_["effBTag_uds"] = fs.make<TH1F>("effBTag_uds", "effBTag_uds", bins_, lowerbin_, upperbin_);
24  hists_["effBTag_other"] = fs.make<TH1F>("effBTag_other", "effBTag_other", bins_, lowerbin_, upperbin_);
25 }
27  edm::BasicAnalyzer::BasicAnalyzer(cfg, fs),
28  Jets_(cfg.getParameter<edm::InputTag>("Jets")),
29  JetsToken_(iC.consumes<std::vector<pat::Jet> >(Jets_)),
30  bTagAlgo_(cfg.getParameter<std::string>("bTagAlgo")),
31  bins_(cfg.getParameter<unsigned int>("bins")),
32  lowerbin_(cfg.getParameter<double>("lowerbin")),
33  upperbin_(cfg.getParameter<double>("upperbin"))
34 {
35  hists_["BTag_b"] = fs.make<TH1F>("BTag_b" , "BTag_b" , bins_, lowerbin_, upperbin_);
36  hists_["BTag_g"] = fs.make<TH1F>("BTag_g" , "BTag_g" , bins_, lowerbin_, upperbin_);
37  hists_["BTag_c"] = fs.make<TH1F>("BTag_c" , "BTag_c" , bins_, lowerbin_, upperbin_);
38  hists_["BTag_uds"] = fs.make<TH1F>("BTag_uds", "BTag_uds", bins_, lowerbin_, upperbin_);
39  hists_["BTag_other"] = fs.make<TH1F>("BTag_other", "BTag_other", bins_, lowerbin_, upperbin_);
40  hists_["effBTag_b"] = fs.make<TH1F>("effBTag_b" , "effBTag_b" , bins_, lowerbin_, upperbin_);
41  hists_["effBTag_g"] = fs.make<TH1F>("effBTag_g" , "effBTag_g" , bins_, lowerbin_, upperbin_);
42  hists_["effBTag_c"] = fs.make<TH1F>("effBTag_c" , "effBTag_c" , bins_, lowerbin_, upperbin_);
43  hists_["effBTag_uds"] = fs.make<TH1F>("effBTag_uds", "effBTag_uds", bins_, lowerbin_, upperbin_);
44  hists_["effBTag_other"] = fs.make<TH1F>("effBTag_other", "effBTag_other", bins_, lowerbin_, upperbin_);
45 }
47 {
48  for(unsigned int i=0; i< bins_; ++i){
49  hists_["effBTag_b"]->SetBinContent(i,hists_["BTag_b"]->Integral(i,hists_["BTag_b"]->GetNbinsX()+1)/hists_["BTag_b"]->Integral(0,hists_["BTag_b"]->GetNbinsX()+1) );
50  hists_["effBTag_g"]->SetBinContent(i,hists_["BTag_g"]->Integral(i,hists_["BTag_g"]->GetNbinsX()+1)/hists_["BTag_g"]->Integral(0,hists_["BTag_g"]->GetNbinsX()+1) );
51  hists_["effBTag_c"]->SetBinContent(i,hists_["BTag_c"]->Integral(i,hists_["BTag_c"]->GetNbinsX()+1)/hists_["BTag_c"]->Integral(0,hists_["BTag_c"]->GetNbinsX()+1) );
52  hists_["effBTag_uds"]->SetBinContent(i,hists_["BTag_uds"]->Integral(i,hists_["BTag_uds"]->GetNbinsX()+1)/hists_["BTag_uds"]->Integral(0,hists_["BTag_uds"]->GetNbinsX()+1) );
53  hists_["effBTag_other"]->SetBinContent(i,hists_["BTag_other"]->Integral(i,hists_["BTag_other"]->GetNbinsX()+1)/hists_["BTag_other"]->Integral(0,hists_["BTag_other"]->GetNbinsX()+1) );
54  }
55 }
57 void
59 {
60  // define what Jet you are using; this is necessary as FWLite is not
61  // capable of reading edm::Views
62  using pat::Jet;
63 
64  // Handle to the Jet collection
66  event.getByLabel(Jets_, Jets);
67 
68  // loop Jet collection and fill histograms
69  for(std::vector<Jet>::const_iterator Jet_it=Jets->begin(); Jet_it!=Jets->end(); ++Jet_it){
70 
71  pat::Jet Jet(*Jet_it);
72 
73  //Categorize the Jets
74  if( abs(Jet.partonFlavour())==5){
75  hists_["BTag_b"]->Fill(Jet.bDiscriminator(bTagAlgo_));
76  }
77  else{
78  if( abs(Jet.partonFlavour())==21 || abs(Jet.partonFlavour())==9 ){
79  hists_["BTag_g"]->Fill(Jet.bDiscriminator(bTagAlgo_));
80  }
81  else{
82  if( abs(Jet.partonFlavour())==4){
83  hists_["BTag_c"]->Fill(Jet.bDiscriminator(bTagAlgo_));
84  }
85  else{
86  if( abs(Jet.partonFlavour())==1 || abs(Jet.partonFlavour())==2 || abs(Jet.partonFlavour())==3){
87  hists_["BTag_uds"]->Fill(Jet.bDiscriminator(bTagAlgo_));
88  }
89  else{
90  hists_["BTag_other"]->Fill(Jet.bDiscriminator(bTagAlgo_));
91  }
92  }
93  }
94  }
95  }
96 }
void analyze(const edm::EventBase &event) override
everything that needs to be done during the event loop
~AnalysisTasksAnalyzerBTag() override
default destructor
edm::EDGetTokenT< std::vector< pat::Jet > > JetsToken_
float bDiscriminator(const std::string &theLabel) const
-— methods for accessing b-tagging info -—
Abstract base class for FWLite and EDM friendly analyzers.
Definition: HeavyIon.h:7
std::map< std::string, TH1 * > hists_
histograms
Definition: Jet.py:1
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T * make(const Args &...args) const
make new ROOT object
int partonFlavour() const
return the parton-based flavour of the jet
AnalysisTasksAnalyzerBTag(const edm::ParameterSet &cfg, TFileDirectory &fs)
default constructor
Analysis-level calorimeter jet class.
Definition: Jet.h:80
HLT enums.
edm::InputTag Jets_
input tag for mouns
Definition: event.py:1