00001
00002
00003
00004
00005
00006
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 #include <TH3F.h>
00029
00030 class RecoTauPlotDiscriminator : public edm::EDAnalyzer {
00031 public:
00032 RecoTauPlotDiscriminator(const edm::ParameterSet& pset);
00033 virtual ~RecoTauPlotDiscriminator() {}
00034 void analyze(const edm::Event &evt, const edm::EventSetup &es);
00035 private:
00036 edm::InputTag src_;
00037 typedef std::map<std::string, TH1*> HistoMap;
00038 typedef std::map<std::string, HistoMap> DiscMap;
00039 typedef std::vector<edm::InputTag> VInputTag;
00040 VInputTag discs_;
00041 DiscMap histos_;
00042 };
00043
00044 RecoTauPlotDiscriminator::RecoTauPlotDiscriminator(const edm::ParameterSet &pset)
00045 :src_(pset.getParameter<edm::InputTag>("src")) {
00046 uint32_t nbins = pset.getParameter<uint32_t>("nbins");
00047 double min = pset.getParameter<double>("min");
00048 double max = pset.getParameter<double>("max");
00049 edm::Service<TFileService> fs;
00050
00051 discs_ = pset.getParameter<std::vector<edm::InputTag> >("discriminators");
00052
00053 BOOST_FOREACH(const edm::InputTag &tag, discs_) {
00054 HistoMap discMap;
00055 discMap["plain"] =
00056 fs->make<TH1F>(tag.label().c_str(), tag.label().c_str(),
00057 nbins, min, max);
00058
00059
00060 std::string vs_pt_name = tag.label()+"_pt";
00061 discMap["vs_pt"] =
00062 fs->make<TH2F>(vs_pt_name.c_str(), vs_pt_name.c_str(),
00063 nbins, min, max, 100, 0, 100);
00064
00065
00066 std::string vs_jetpt_name = tag.label()+"_jetPt";
00067 discMap["vs_jetPt"] =
00068 fs->make<TH2F>(vs_jetpt_name.c_str(), vs_jetpt_name.c_str(),
00069 nbins, min, max, 100, 0, 100);
00070
00071
00072 std::string vs_embedpt_name = tag.label()+"_embedPt";
00073 discMap["vs_embedPt"] =
00074 fs->make<TH2F>(vs_embedpt_name.c_str(), vs_embedpt_name.c_str(),
00075 nbins, min, max, 100, 0, 100);
00076
00077
00078 std::string vs_pt_jetPt_name = tag.label()+"_pt_jetPt";
00079 discMap["vs_pt_jetPt"] =
00080 fs->make<TH3F>(vs_pt_jetPt_name.c_str(), vs_pt_jetPt_name.c_str(),
00081 nbins, min, max, 100, 0, 100, 100, 0, 100);
00082
00083 std::string vs_pt_embedPt_name = tag.label()+"_pt_embedPt";
00084 discMap["vs_pt_embedPt"] =
00085 fs->make<TH3F>(vs_pt_embedPt_name.c_str(), vs_pt_embedPt_name.c_str(),
00086 nbins, min, max, 100, 0, 100, 100, 0, 100);
00087
00088
00089 std::string vs_eta_name = tag.label()+"_eta";
00090 discMap["vs_eta"] =
00091 fs->make<TH2F>(vs_eta_name.c_str(), vs_eta_name.c_str(),
00092 nbins, min, max, 100, -2.5, 2.5);
00093
00094 std::string vs_dm_name = tag.label()+"_dm";
00095 discMap["vs_dm"] =
00096 fs->make<TH2F>(vs_dm_name.c_str(), vs_dm_name.c_str(),
00097 nbins, min, max, 15, -0.5, 14.5);
00098 histos_[tag.label()] = discMap;
00099 }
00100 }
00101
00102 void
00103 RecoTauPlotDiscriminator::analyze(const edm::Event &evt,
00104 const edm::EventSetup &es) {
00105
00106 edm::Handle<reco::CandidateView> input;
00107 evt.getByLabel(src_, input);
00108
00109
00110 reco::PFTauRefVector inputRefs =
00111 reco::tau::castView<reco::PFTauRefVector>(input);
00112
00113
00114 BOOST_FOREACH(const reco::PFTauRef& tau, inputRefs) {
00115
00116 BOOST_FOREACH(const edm::InputTag &tag, discs_) {
00117 edm::Handle<reco::PFTauDiscriminator> discHandle;
00118 evt.getByLabel(tag, discHandle);
00119
00120 double result = (*discHandle)[tau];
00121 HistoMap& mymap = histos_[tag.label()];
00122 mymap["plain"]->Fill(result);
00123 mymap["vs_pt"]->Fill(result, tau->pt());
00124 mymap["vs_jetPt"]->Fill(result, tau->jetRef()->pt());
00125 mymap["vs_embedPt"]->Fill(result, tau->alternatLorentzVect().pt());
00126 dynamic_cast<TH3F*>(mymap["vs_pt_jetPt"])->Fill(
00127 result, tau->pt(), tau->jetRef()->pt());
00128 dynamic_cast<TH3F*>(mymap["vs_pt_embedPt"])->Fill(
00129 result, tau->pt(), tau->alternatLorentzVect().pt());
00130 mymap["vs_eta"]->Fill(result, tau->eta());
00131 mymap["vs_dm"]->Fill(result, tau->decayMode());
00132 }
00133 }
00134 }
00135
00136 #include "FWCore/Framework/interface/MakerMacros.h"
00137 DEFINE_FWK_MODULE(RecoTauPlotDiscriminator);