Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:02

0001 #include "DataFormats/Common/interface/Handle.h"
0002 #include "DataFormats/JetReco/interface/Jet.h"
0003 #include "PhysicsTools/PatExamples/interface/AnalysisTasksAnalyzerBTag.h"
0004 
0005 /// default constructor
0006 AnalysisTasksAnalyzerBTag::AnalysisTasksAnalyzerBTag(const edm::ParameterSet& cfg, TFileDirectory& fs)
0007     : edm::BasicAnalyzer::BasicAnalyzer(cfg, fs),
0008       Jets_(cfg.getParameter<edm::InputTag>("Jets")),
0009       bTagAlgo_(cfg.getParameter<std::string>("bTagAlgo")),
0010       bins_(cfg.getParameter<unsigned int>("bins")),
0011       lowerbin_(cfg.getParameter<double>("lowerbin")),
0012       upperbin_(cfg.getParameter<double>("upperbin")) {
0013   hists_["BTag_b"] = fs.make<TH1F>("BTag_b", "BTag_b", bins_, lowerbin_, upperbin_);
0014   hists_["BTag_g"] = fs.make<TH1F>("BTag_g", "BTag_g", bins_, lowerbin_, upperbin_);
0015   hists_["BTag_c"] = fs.make<TH1F>("BTag_c", "BTag_c", bins_, lowerbin_, upperbin_);
0016   hists_["BTag_uds"] = fs.make<TH1F>("BTag_uds", "BTag_uds", bins_, lowerbin_, upperbin_);
0017   hists_["BTag_other"] = fs.make<TH1F>("BTag_other", "BTag_other", bins_, lowerbin_, upperbin_);
0018   hists_["effBTag_b"] = fs.make<TH1F>("effBTag_b", "effBTag_b", bins_, lowerbin_, upperbin_);
0019   hists_["effBTag_g"] = fs.make<TH1F>("effBTag_g", "effBTag_g", bins_, lowerbin_, upperbin_);
0020   hists_["effBTag_c"] = fs.make<TH1F>("effBTag_c", "effBTag_c", bins_, lowerbin_, upperbin_);
0021   hists_["effBTag_uds"] = fs.make<TH1F>("effBTag_uds", "effBTag_uds", bins_, lowerbin_, upperbin_);
0022   hists_["effBTag_other"] = fs.make<TH1F>("effBTag_other", "effBTag_other", bins_, lowerbin_, upperbin_);
0023 }
0024 AnalysisTasksAnalyzerBTag::AnalysisTasksAnalyzerBTag(const edm::ParameterSet& cfg,
0025                                                      TFileDirectory& fs,
0026                                                      edm::ConsumesCollector&& iC)
0027     : edm::BasicAnalyzer::BasicAnalyzer(cfg, fs),
0028       Jets_(cfg.getParameter<edm::InputTag>("Jets")),
0029       JetsToken_(iC.consumes<std::vector<pat::Jet> >(Jets_)),
0030       bTagAlgo_(cfg.getParameter<std::string>("bTagAlgo")),
0031       bins_(cfg.getParameter<unsigned int>("bins")),
0032       lowerbin_(cfg.getParameter<double>("lowerbin")),
0033       upperbin_(cfg.getParameter<double>("upperbin")) {
0034   hists_["BTag_b"] = fs.make<TH1F>("BTag_b", "BTag_b", bins_, lowerbin_, upperbin_);
0035   hists_["BTag_g"] = fs.make<TH1F>("BTag_g", "BTag_g", bins_, lowerbin_, upperbin_);
0036   hists_["BTag_c"] = fs.make<TH1F>("BTag_c", "BTag_c", bins_, lowerbin_, upperbin_);
0037   hists_["BTag_uds"] = fs.make<TH1F>("BTag_uds", "BTag_uds", bins_, lowerbin_, upperbin_);
0038   hists_["BTag_other"] = fs.make<TH1F>("BTag_other", "BTag_other", bins_, lowerbin_, upperbin_);
0039   hists_["effBTag_b"] = fs.make<TH1F>("effBTag_b", "effBTag_b", bins_, lowerbin_, upperbin_);
0040   hists_["effBTag_g"] = fs.make<TH1F>("effBTag_g", "effBTag_g", bins_, lowerbin_, upperbin_);
0041   hists_["effBTag_c"] = fs.make<TH1F>("effBTag_c", "effBTag_c", bins_, lowerbin_, upperbin_);
0042   hists_["effBTag_uds"] = fs.make<TH1F>("effBTag_uds", "effBTag_uds", bins_, lowerbin_, upperbin_);
0043   hists_["effBTag_other"] = fs.make<TH1F>("effBTag_other", "effBTag_other", bins_, lowerbin_, upperbin_);
0044 }
0045 AnalysisTasksAnalyzerBTag::~AnalysisTasksAnalyzerBTag() {
0046   for (unsigned int i = 0; i < bins_; ++i) {
0047     hists_["effBTag_b"]->SetBinContent(i,
0048                                        hists_["BTag_b"]->Integral(i, hists_["BTag_b"]->GetNbinsX() + 1) /
0049                                            hists_["BTag_b"]->Integral(0, hists_["BTag_b"]->GetNbinsX() + 1));
0050     hists_["effBTag_g"]->SetBinContent(i,
0051                                        hists_["BTag_g"]->Integral(i, hists_["BTag_g"]->GetNbinsX() + 1) /
0052                                            hists_["BTag_g"]->Integral(0, hists_["BTag_g"]->GetNbinsX() + 1));
0053     hists_["effBTag_c"]->SetBinContent(i,
0054                                        hists_["BTag_c"]->Integral(i, hists_["BTag_c"]->GetNbinsX() + 1) /
0055                                            hists_["BTag_c"]->Integral(0, hists_["BTag_c"]->GetNbinsX() + 1));
0056     hists_["effBTag_uds"]->SetBinContent(i,
0057                                          hists_["BTag_uds"]->Integral(i, hists_["BTag_uds"]->GetNbinsX() + 1) /
0058                                              hists_["BTag_uds"]->Integral(0, hists_["BTag_uds"]->GetNbinsX() + 1));
0059     hists_["effBTag_other"]->SetBinContent(
0060         i,
0061         hists_["BTag_other"]->Integral(i, hists_["BTag_other"]->GetNbinsX() + 1) /
0062             hists_["BTag_other"]->Integral(0, hists_["BTag_other"]->GetNbinsX() + 1));
0063   }
0064 }
0065 /// everything that needs to be done during the event loop
0066 void AnalysisTasksAnalyzerBTag::analyze(const edm::EventBase& event) {
0067   // define what Jet you are using; this is necessary as FWLite is not
0068   // capable of reading edm::Views
0069   using pat::Jet;
0070 
0071   // Handle to the Jet collection
0072   edm::Handle<std::vector<Jet> > Jets;
0073   event.getByLabel(Jets_, Jets);
0074 
0075   // loop Jet collection and fill histograms
0076   for (std::vector<Jet>::const_iterator Jet_it = Jets->begin(); Jet_it != Jets->end(); ++Jet_it) {
0077     const pat::Jet& Jet(*Jet_it);
0078 
0079     //Categorize the Jets
0080     if (abs(Jet.partonFlavour()) == 5) {
0081       hists_["BTag_b"]->Fill(Jet.bDiscriminator(bTagAlgo_));
0082     } else {
0083       if (abs(Jet.partonFlavour()) == 21 || abs(Jet.partonFlavour()) == 9) {
0084         hists_["BTag_g"]->Fill(Jet.bDiscriminator(bTagAlgo_));
0085       } else {
0086         if (abs(Jet.partonFlavour()) == 4) {
0087           hists_["BTag_c"]->Fill(Jet.bDiscriminator(bTagAlgo_));
0088         } else {
0089           if (abs(Jet.partonFlavour()) == 1 || abs(Jet.partonFlavour()) == 2 || abs(Jet.partonFlavour()) == 3) {
0090             hists_["BTag_uds"]->Fill(Jet.bDiscriminator(bTagAlgo_));
0091           } else {
0092             hists_["BTag_other"]->Fill(Jet.bDiscriminator(bTagAlgo_));
0093           }
0094         }
0095       }
0096     }
0097   }
0098 }