CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/PhysicsTools/PatExamples/plugins/PatBTagAnalyzer.cc

Go to the documentation of this file.
00001 #include <map>
00002 #include <string>
00003 
00004 #include "TH1D.h"
00005 #include "TH2D.h"
00006 #include "TGraphErrors.h"
00007 
00008 #include "FWCore/Framework/interface/Event.h"
00009 #include "FWCore/Framework/interface/EDAnalyzer.h"
00010 #include "FWCore/Utilities/interface/InputTag.h"
00011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00012 #include "FWCore/ServiceRegistry/interface/Service.h"
00013 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00014 #include "FWCore/Framework/interface/MakerMacros.h"
00015 
00016 #include "DataFormats/PatCandidates/interface/Jet.h"
00017 #include "PhysicsTools/PatUtils/interface/bJetSelector.h"
00018 #include "PhysicsTools/PatExamples/interface/BTagPerformance.h"
00019 #include "PhysicsTools/PatExamples/interface/PatBTagCommonHistos.h"
00020 
00021 
00022 class PatBTagAnalyzer : public edm::EDAnalyzer {
00023 
00024 public:
00025 
00026   explicit PatBTagAnalyzer(const edm::ParameterSet&);
00027   ~PatBTagAnalyzer();
00028     
00029 private:
00030 
00031   virtual void beginJob() ;
00032   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00033   virtual void endJob() ;
00034   
00035   
00036   edm::InputTag jetLabel_;
00037   edm::ParameterSet PatBjet_;
00038 
00039   std::string  BTagpurity_;
00040   std::string  BTagmethod_;
00041   std::string  BTagdiscriminator_;
00042 
00043   bool    BTagverbose;
00044   double  BTagdisccut_;
00045   double  BTagdiscmax_;
00046 
00047 
00048   std::string  discname[10];   
00049   std::string  bname   [10];   
00050   std::string  cname   [10];   
00051   BTagPerformance BTagPerf[10];
00052   std::map<int,std::string> udsgname;   
00053 
00056   std::map<std::string,TH1D*> histocontainer_;   
00059   std::map<std::string,TH2D*> h2_; 
00062   std::map<std::string,TGraph*> graphcontainer_; 
00065   std::map<std::string,TGraphErrors*> grapherrorscontainer_; 
00066 
00067   bJetSelector BTagger;
00068   PatBTagCommonHistos BTagHistograms;
00069 };
00070 
00071 PatBTagAnalyzer::PatBTagAnalyzer(const edm::ParameterSet& iConfig):
00072   jetLabel_(iConfig.getUntrackedParameter<edm::InputTag>("jetTag")),
00073   PatBjet_(iConfig.getParameter< edm::ParameterSet >("BjetTag")),
00074   BTagpurity_(PatBjet_.getParameter<std::string>("purity")),
00075   BTagmethod_(PatBjet_.getUntrackedParameter<std::string>("tagger","TC2")),
00076   BTagdiscriminator_(PatBjet_.getParameter<std::string>("discriminator")),
00077   BTagverbose(PatBjet_.getUntrackedParameter<bool>("verbose",false)),
00078   BTagdisccut_(PatBjet_.getUntrackedParameter<double>("mindiscriminatorcut",5.0)),
00079   BTagdiscmax_(PatBjet_.getUntrackedParameter<double>("maxdiscriminatorcut",15.0)),
00080   BTagger(iConfig.getParameter< edm::ParameterSet >("BJetOperatingPoints")),
00081   BTagHistograms(iConfig)
00082 {
00083    //now do what ever initialization is needed
00084 }
00085 
00086 
00087 PatBTagAnalyzer::~PatBTagAnalyzer()
00088 {
00089  
00090    // do anything here that needs to be done at desctruction time
00091    // (e.g. close files, deallocate resources etc.)
00092 
00093 }
00094 
00095 // ------------ method called to for each event  ------------
00096 void
00097 PatBTagAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00098 {
00099    using namespace edm;
00100 
00101    // first: get all objects from the event.
00102 
00103    edm::Handle<edm::View<pat::Jet> > jetHandle;
00104    iEvent.getByLabel(jetLabel_,jetHandle);
00105    edm::View<pat::Jet> jets = *jetHandle; // get JETS
00106 
00107    // LOOP over all jets
00108    
00109    for(edm::View<pat::Jet>::const_iterator jet_iter = jets.begin(); jet_iter!=jets.end(); ++jet_iter){
00110 
00111        float bdiscriminator = jet_iter->bDiscriminator(BTagdiscriminator_);
00112        int flavor           = jet_iter->partonFlavour();
00113        //
00114        // Fill in for performance standard pt(uncorrected) >10 and abs(eta)<2.4 
00115        if( jet_iter->correctedJet("raw").pt() > 10  &&
00116                    fabs(jet_iter->eta()) < 2.4 ) {
00117                    
00118                    BTagPerf[0].Add(bdiscriminator,abs(flavor));
00119 
00120        }
00121        
00122            //Fill histograms
00123            BTagHistograms.Fill(jet_iter,"all");
00124            if (flavor ==  0  ) BTagHistograms.Fill(jet_iter,"no_flavor");
00125            if (flavor ==  5 || flavor ==  -5 ) BTagHistograms.Fill(jet_iter,"b");
00126            if (flavor ==  4 || flavor ==  -4 ) BTagHistograms.Fill(jet_iter,"c");
00127            if ((-4 < flavor && flavor < 4 && flavor != 0 )||(flavor == 21 || flavor == -21 ))
00128                    BTagHistograms.Fill(jet_iter,"udsg");
00129            
00130 
00131    }//end loop over jets 
00132 
00133 }
00134 // ------------ method called once each job just before starting event loop  ------------
00135 void 
00136 PatBTagAnalyzer::beginJob()
00137 {
00138   //
00139   // define some histograms using the framework tfileservice. Define the output file name in your .cfg.
00140   //
00141   edm::Service<TFileService> fs;
00142   
00143   TString suffix1="_test";
00144 
00145   //set performance variables collector
00146   for (int i=0; i < 10; i++){
00147     BTagPerf[i].Set(BTagmethod_);
00148     BTagPerf[i].SetMinDiscriminator(BTagdisccut_);
00149     BTagPerf[i].SetMaxDiscriminator(BTagdiscmax_);
00150   }
00151 
00152   histocontainer_["njets"]=fs->make<TH1D>("njets","jet multiplicity for jets with p_{T} > 50 GeV/c",10,0,10);
00153 // Std. 30 pt uncorr cut for performance
00154   discname[0]="disc"+BTagmethod_+"_udsg";   
00155   bname[0]   ="g"+BTagmethod_+"_b";   
00156   cname[0]   ="g"+BTagmethod_+"_c";   
00157   udsgname[0]="g"+BTagmethod_+"_udsg";   
00158 
00159 // 10 pt uncorr for performance + all,>0,>1,>2 tracks
00160   discname[1]="Uncor10_disc"+BTagmethod_+"_udsg";
00161   bname[1]   ="Uncor10_g"+BTagmethod_+"_b";   
00162   cname[1]   ="Uncor10_g"+BTagmethod_+"_c";   
00163   udsgname[1]="Uncor10_g"+BTagmethod_+"_udsg";
00164   discname[2]="Uncor10t0_disc"+BTagmethod_+"_udsg";
00165   bname[2]   ="Uncor10t0_g"+BTagmethod_+"_b";   
00166   cname[2]   ="Uncor10t0_g"+BTagmethod_+"_c";   
00167   udsgname[2]="Uncor10t0_g"+BTagmethod_+"_udsg";
00168   discname[3]="Uncor10t1_disc"+BTagmethod_+"_udsg";
00169   bname[3]   ="Uncor10t1_g"+BTagmethod_+"_b";   
00170   cname[3]   ="Uncor10t1_g"+BTagmethod_+"_c";   
00171   udsgname[3]="Uncor10t1_g"+BTagmethod_+"_udsg";
00172   discname[4]="Uncor10t2_disc"+BTagmethod_+"_udsg";
00173   bname[4]   ="Uncor10t2_g"+BTagmethod_+"_b";   
00174   cname[4]   ="Uncor10t2_g"+BTagmethod_+"_c";   
00175   udsgname[4]="Uncor10t2_g"+BTagmethod_+"_udsg";
00176 
00177 // 30 pt corr for performance + all,>0,>1,>2 tracks   
00178   discname[5]="Corr30_disc"+BTagmethod_+"_udsg";
00179   bname[5]   ="Corr30_g"+BTagmethod_+"_b";   
00180   cname[5]   ="Corr30_g"+BTagmethod_+"_c";   
00181   udsgname[5]="Corr30_g"+BTagmethod_+"_udsg";
00182   discname[6]="Corr30t0_disc"+BTagmethod_+"_udsg";
00183   bname[6]   ="Corr30t0_g"+BTagmethod_+"_b";   
00184   cname[6]   ="Corr30t0_g"+BTagmethod_+"_c";   
00185   udsgname[6]="Corr30t0_g"+BTagmethod_+"_udsg";
00186   discname[7]="Corr30t1_disc"+BTagmethod_+"_udsg";
00187   bname[7]   ="Corr30t1_g"+BTagmethod_+"_b";   
00188   cname[7]   ="Corr30t1_g"+BTagmethod_+"_c";   
00189   udsgname[7]="Corr30t1_g"+BTagmethod_+"_udsg";
00190   discname[8]="Corr30t2_disc"+BTagmethod_+"_udsg";
00191   bname[8]   ="Corr30t2_g"+BTagmethod_+"_b";   
00192   cname[8]   ="Corr30t2_g"+BTagmethod_+"_c";   
00193   udsgname[8]="Corr30t2_g"+BTagmethod_+"_udsg";
00194 
00195 // check filter
00196   discname[9]="check_disc"+BTagmethod_+"_udsg";   
00197   bname[9]   ="check_g"+BTagmethod_+"_b";   
00198   cname[9]   ="check_g"+BTagmethod_+"_c";   
00199   udsgname[9]="check_g"+BTagmethod_+"_udsg";   
00200 
00201   for(int i=1; i<10;i++){
00202     graphcontainer_[discname[i]]      =fs->make<TGraph>(BTagPerf[i].GetN());       graphcontainer_[discname[i]]->SetName(discname[i].c_str());
00203     grapherrorscontainer_[bname[i]]   =fs->make<TGraphErrors>(BTagPerf[i].GetN()); grapherrorscontainer_[bname[i]]   ->SetName(bname[i].c_str());
00204     grapherrorscontainer_[cname[i]]   =fs->make<TGraphErrors>(BTagPerf[i].GetN()); grapherrorscontainer_[cname[i]]   ->SetName(cname[i].c_str());
00205     grapherrorscontainer_[udsgname[i]]=fs->make<TGraphErrors>(BTagPerf[i].GetN()); grapherrorscontainer_[udsgname[i]]->SetName(udsgname[i].c_str());   
00206   }
00207   //Define histograms
00208   BTagHistograms.Set("all");  
00209   BTagHistograms.Set("no_flavor");  
00210   BTagHistograms.Set("b");  
00211   BTagHistograms.Set("c");  
00212   BTagHistograms.Set("udsg");  
00213 
00214   // Set to save histogram errors
00215   BTagHistograms.Sumw2();  
00216 
00217 }
00218 
00219 // ------------ method called once each job just after ending the event loop  ------------
00220 void 
00221 PatBTagAnalyzer::endJob() {
00222 //ed++
00223    edm::Service<TFileService> fs;
00224 
00225 // Save performance plots as Tgraphs
00226 
00227 
00228    for (int i=1;i<10;i++){
00229       BTagPerf[i].Eval();
00230       for (int n=0; n<BTagPerf[i].GetN();  n++ ){
00231          graphcontainer_[discname[i]]       ->SetPoint(n,BTagPerf[i].GetArray("udsg")[n],BTagPerf[i].GetArray("discriminator")[n]);
00232          grapherrorscontainer_[bname[i]]    ->SetPoint(n,BTagPerf[i].GetArray("b")[n],BTagPerf[i].GetArray("b")[n]);
00233          grapherrorscontainer_[cname[i]]    ->SetPoint(n,BTagPerf[i].GetArray("b")[n],BTagPerf[i].GetArray("c")[n]);
00234          grapherrorscontainer_[udsgname[i]] ->SetPoint(n,BTagPerf[i].GetArray("b")[n],BTagPerf[i].GetArray("udsg")[n]);
00235          grapherrorscontainer_[bname[i]]    ->SetPointError(n,BTagPerf[i].GetArray("bErr")[n],BTagPerf[i].GetArray("bErr")[n]);
00236          grapherrorscontainer_[cname[i]]    ->SetPointError(n,BTagPerf[i].GetArray("bErr")[n],BTagPerf[i].GetArray("cErr")[n]);
00237          grapherrorscontainer_[udsgname[i]] ->SetPointError(n,BTagPerf[i].GetArray("bErr")[n],BTagPerf[i].GetArray("udsgErr")[n]);
00238       }//end for over BTagPerf[i] elements
00239       graphcontainer_[discname[i]]     ->SetTitle("Jet udsg-mistagging");
00240       grapherrorscontainer_[bname[i]]   ->SetTitle("Jet b-efficiency");
00241       grapherrorscontainer_[cname[i]]   ->SetTitle("Jet c-mistagging");
00242       grapherrorscontainer_[udsgname[i]]->SetTitle("discriminator vs udsg-mistagging");
00243    }//end for over [i]
00244 
00245 
00246 // Save default cut performance plot
00247    BTagPerf[0].Eval();
00248 
00249 //  TFileDirectory TaggerDir = fs->mkdir(BTagmethod_);
00250 //  TGraphErrors *BTagger_b    = new TGraphErrors(BTagTool.GetN(),
00251   TGraphErrors *BTagger_b    = fs->mkdir(BTagmethod_).make<TGraphErrors>(BTagPerf[0].GetN(),
00252                                             BTagPerf[0].GetArray("b").GetArray(),BTagPerf[0].GetArray("b").GetArray(),
00253                                             BTagPerf[0].GetArray("bErr").GetArray(),BTagPerf[0].GetArray("bErr").GetArray());
00254         
00255   TGraphErrors *BTagger_c    = new TGraphErrors(BTagPerf[0].GetN(),
00256                                             BTagPerf[0].GetArray("b").GetArray(),BTagPerf[0].GetArray("c").GetArray(),
00257                                             BTagPerf[0].GetArray("bErr").GetArray(),BTagPerf[0].GetArray("cErr").GetArray());
00258                 
00259   TGraphErrors *BTagger_udsg = new TGraphErrors(BTagPerf[0].GetN(),
00260                                                BTagPerf[0].GetArray("b").GetArray(),BTagPerf[0].GetArray("udsg").GetArray(),
00261                                                BTagPerf[0].GetArray("bErr").GetArray(),BTagPerf[0].GetArray("udsgErr").GetArray());
00262   TGraph *discBTagger_udsg   = new TGraph(BTagPerf[0].GetN(),
00263                                    BTagPerf[0].GetArray("udsg").GetArray(),
00264                                    BTagPerf[0].GetArray("discriminator").GetArray());
00265  
00266   BTagger_b->SetName(bname[0].c_str());
00267   BTagger_c->SetName(cname[0].c_str());
00268   BTagger_udsg->SetName(udsgname[0].c_str());
00269   discBTagger_udsg->SetName(discname[0].c_str());
00270 
00271   BTagger_b->SetTitle("Jet b-efficiency");
00272   BTagger_c->SetTitle("Jet c-mistagging");
00273   BTagger_udsg->SetTitle("Jet udsg-mistagging");
00274   discBTagger_udsg->SetTitle("discriminator vs udsg-mistagging");
00275 
00276 
00277   for (int i=1;i<10;i++){
00278    graphcontainer_[discname[i]]      ->Write();
00279    grapherrorscontainer_[bname[i]]   ->Write();
00280    grapherrorscontainer_[cname[i]]   ->Write();
00281    grapherrorscontainer_[udsgname[i]]->Write();
00282   }
00283   
00284   BTagger_b->Write();
00285   BTagger_c->Write();
00286   BTagger_udsg->Write();
00287   discBTagger_udsg->Write();
00288     
00289 
00290 }
00291 
00292 //define this as a plug-in
00293 DEFINE_FWK_MODULE(PatBTagAnalyzer);