CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoTauTag/TauTagTools/plugins/RecoTauPlotDiscriminator.cc

Go to the documentation of this file.
00001 /*
00002  * RecoTauPlotDiscriminator
00003  *
00004  * Plot the output of a PFTauDiscriminator using TFileService.
00005  *
00006  * Author: Evan K. Friis (UC Davis)
00007  */
00008 
00009 #include <boost/foreach.hpp>
00010 
00011 #include "FWCore/Framework/interface/EDAnalyzer.h"
00012 #include "FWCore/Framework/interface/ESHandle.h"
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015 #include "FWCore/Utilities/interface/InputTag.h"
00016 
00017 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00018 
00019 #include "DataFormats/TauReco/interface/PFTau.h"
00020 #include "DataFormats/TauReco/interface/PFTauFwd.h"
00021 #include "DataFormats/TauReco/interface/PFTauDiscriminator.h"
00022 
00023 #include "FWCore/ServiceRegistry/interface/Service.h"
00024 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00025 
00026 #include <TH1F.h>
00027 #include <TH2F.h>
00028 
00029 class RecoTauPlotDiscriminator : public edm::EDAnalyzer {
00030   public:
00031     RecoTauPlotDiscriminator(const edm::ParameterSet& pset);
00032     virtual ~RecoTauPlotDiscriminator() {}
00033     void analyze(const edm::Event &evt, const edm::EventSetup &es);
00034   private:
00035     edm::InputTag src_;
00036     typedef std::map<std::string, TH1*> HistoMap;
00037     typedef std::map<std::string, HistoMap> DiscMap;
00038     typedef std::vector<edm::InputTag> VInputTag;
00039     VInputTag discs_;
00040     DiscMap histos_;
00041 };
00042 
00043 RecoTauPlotDiscriminator::RecoTauPlotDiscriminator(const edm::ParameterSet &pset)
00044   :src_(pset.getParameter<edm::InputTag>("src")) {
00045   uint32_t nbins = pset.getParameter<uint32_t>("nbins");
00046   double min = pset.getParameter<double>("min");
00047   double max = pset.getParameter<double>("max");
00048   edm::Service<TFileService> fs;
00049   // Get the discriminators
00050   discs_ = pset.getParameter<std::vector<edm::InputTag> >("discriminators");
00051 
00052   BOOST_FOREACH(const edm::InputTag &tag, discs_) {
00053     HistoMap discMap;
00054     discMap["plain"] =
00055         fs->make<TH1F>(tag.label().c_str(), tag.label().c_str(),
00056                        nbins, min, max);
00057 
00058     // Make correlation plots
00059     std::string vs_pt_name = tag.label()+"_pt";
00060     discMap["vs_pt"] =
00061         fs->make<TH2F>(vs_pt_name.c_str(), vs_pt_name.c_str(),
00062                        nbins, min, max, 100, 0, 100);
00063 
00064     std::string vs_eta_name = tag.label()+"_eta";
00065     discMap["vs_eta"] =
00066         fs->make<TH2F>(vs_eta_name.c_str(), vs_eta_name.c_str(),
00067                        nbins, min, max, 100, -2.5, 2.5);
00068 
00069     std::string vs_dm_name = tag.label()+"_dm";
00070     discMap["vs_dm"] =
00071         fs->make<TH2F>(vs_dm_name.c_str(), vs_dm_name.c_str(),
00072                        nbins, min, max, 15, -0.5, 14.5);
00073     histos_[tag.label()] = discMap;
00074   }
00075 }
00076 
00077 void
00078 RecoTauPlotDiscriminator::analyze(const edm::Event &evt,
00079                                   const edm::EventSetup &es) {
00080   // Get the input collection to clean
00081   edm::Handle<reco::CandidateView> input;
00082   evt.getByLabel(src_, input);
00083 
00084   // Cast the input candidates to Refs to real taus
00085   reco::PFTauRefVector inputRefs =
00086       reco::tau::castView<reco::PFTauRefVector>(input);
00087 
00088   // Plot the discriminator output for each of our taus
00089   BOOST_FOREACH(const reco::PFTauRef& tau, inputRefs) {
00090     // Plot each discriminator
00091     BOOST_FOREACH(const edm::InputTag &tag, discs_) {
00092       edm::Handle<reco::PFTauDiscriminator> discHandle;
00093       evt.getByLabel(tag, discHandle);
00094       //const HistoMap &discHistos = disc.second;
00095       double result = (*discHandle)[tau];
00096       HistoMap& mymap = histos_[tag.label()];
00097       mymap["plain"]->Fill(result);
00098       mymap["vs_pt"]->Fill(result, tau->pt());
00099       mymap["vs_eta"]->Fill(result, tau->eta());
00100       mymap["vs_dm"]->Fill(result, tau->decayMode());
00101     }
00102   }
00103 }
00104 
00105 #include "FWCore/Framework/interface/MakerMacros.h"
00106 DEFINE_FWK_MODULE(RecoTauPlotDiscriminator);