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 "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
00024 #include "DataFormats/VertexReco/interface/Vertex.h"
00025 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00026
00027 #include "FWCore/ServiceRegistry/interface/Service.h"
00028 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00029
00030 #include <TH1F.h>
00031 #include <TH2F.h>
00032 #include <TH3F.h>
00033
00034 class RecoTauPlotDiscriminator : public edm::EDAnalyzer {
00035 public:
00036 RecoTauPlotDiscriminator(const edm::ParameterSet& pset);
00037 virtual ~RecoTauPlotDiscriminator() {}
00038 void analyze(const edm::Event &evt, const edm::EventSetup &es);
00039 private:
00040 edm::InputTag src_;
00041 typedef std::map<std::string, TH1*> HistoMap;
00042 typedef std::map<std::string, HistoMap> DiscMap;
00043 typedef std::vector<edm::InputTag> VInputTag;
00044 VInputTag discs_;
00045 bool plotPU_;
00046 double pileupPtCut_;
00047 edm::InputTag pileupInfoSrc_;
00048 edm::InputTag pileupVerticesSrc_;
00049 DiscMap histos_;
00050 };
00051
00052 RecoTauPlotDiscriminator::RecoTauPlotDiscriminator(const edm::ParameterSet &pset)
00053 :src_(pset.getParameter<edm::InputTag>("src")) {
00054 uint32_t nbins = pset.getParameter<uint32_t>("nbins");
00055 double min = pset.getParameter<double>("min");
00056 double max = pset.getParameter<double>("max");
00057 edm::Service<TFileService> fs;
00058
00059 discs_ = pset.getParameter<std::vector<edm::InputTag> >("discriminators");
00060
00061 plotPU_ = pset.getParameter<bool>("plotPU");
00062 if (plotPU_) {
00063 pileupPtCut_ = pset.getParameter<double>("pileupTauPtCut");
00064 pileupInfoSrc_ = pset.getParameter<edm::InputTag>("pileupInfo");
00065 pileupVerticesSrc_ = pset.getParameter<edm::InputTag>("pileupVertices");
00066 }
00067
00068 BOOST_FOREACH(const edm::InputTag &tag, discs_) {
00069 HistoMap discMap;
00070 discMap["plain"] =
00071 fs->make<TH1F>(tag.label().c_str(), tag.label().c_str(),
00072 nbins, min, max);
00073
00074
00075 std::string vs_pt_name = tag.label()+"_pt";
00076 discMap["vs_pt"] =
00077 fs->make<TH2F>(vs_pt_name.c_str(), vs_pt_name.c_str(),
00078 nbins, min, max, 100, 0, 200);
00079
00080
00081 std::string vs_jetpt_name = tag.label()+"_jetPt";
00082 discMap["vs_jetPt"] =
00083 fs->make<TH2F>(vs_jetpt_name.c_str(), vs_jetpt_name.c_str(),
00084 nbins, min, max, 100, 0, 200);
00085
00086
00087 std::string vs_embedpt_name = tag.label()+"_embedPt";
00088 discMap["vs_embedPt"] =
00089 fs->make<TH2F>(vs_embedpt_name.c_str(), vs_embedpt_name.c_str(),
00090 nbins, min, max, 100, 0, 200);
00091
00092
00093 std::string vs_pt_jetPt_name = tag.label()+"_pt_jetPt";
00094 discMap["vs_pt_jetPt"] =
00095 fs->make<TH3F>(vs_pt_jetPt_name.c_str(), vs_pt_jetPt_name.c_str(),
00096 nbins, min, max, 100, 0, 200, 100, 0, 200);
00097
00098 std::string vs_pt_embedPt_name = tag.label()+"_pt_embedPt";
00099 discMap["vs_pt_embedPt"] =
00100 fs->make<TH3F>(vs_pt_embedPt_name.c_str(), vs_pt_embedPt_name.c_str(),
00101 nbins, min, max, 100, 0, 200, 100, 0, 200);
00102
00103
00104 std::string vs_eta_name = tag.label()+"_eta";
00105 discMap["vs_eta"] =
00106 fs->make<TH2F>(vs_eta_name.c_str(), vs_eta_name.c_str(),
00107 nbins, min, max, 100, -2.5, 2.5);
00108
00109 std::string vs_dm_name = tag.label()+"_dm";
00110 discMap["vs_dm"] =
00111 fs->make<TH2F>(vs_dm_name.c_str(), vs_dm_name.c_str(),
00112 nbins, min, max, 15, -0.5, 14.5);
00113
00114 if (plotPU_) {
00115 std::string vs_truePU_name = tag.label()+"_truePU";
00116 discMap["vs_truePU"] = fs->make<TH2F>(vs_truePU_name.c_str(),
00117 vs_truePU_name.c_str(), nbins, min, max, 15, -0.5, 14.5);
00118 std::string vs_recoPU_name = tag.label()+"_recoPU";
00119 discMap["vs_recoPU"] = fs->make<TH2F>(vs_recoPU_name.c_str(),
00120 vs_recoPU_name.c_str(), nbins, min, max, 15, -0.5, 14.5);
00121 }
00122
00123 histos_[tag.label()] = discMap;
00124 }
00125 }
00126
00127 void
00128 RecoTauPlotDiscriminator::analyze(const edm::Event &evt,
00129 const edm::EventSetup &es) {
00130
00131 edm::Handle<reco::CandidateView> input;
00132 evt.getByLabel(src_, input);
00133
00134
00135 reco::PFTauRefVector inputRefs =
00136 reco::tau::castView<reco::PFTauRefVector>(input);
00137
00138 edm::Handle<PileupSummaryInfo> puInfo;
00139 edm::Handle<reco::VertexCollection> puVertices;
00140 if (plotPU_) {
00141 evt.getByLabel(pileupInfoSrc_, puInfo);
00142 evt.getByLabel(pileupVerticesSrc_, puVertices);
00143 }
00144
00145
00146 BOOST_FOREACH(const reco::PFTauRef& tau, inputRefs) {
00147
00148 BOOST_FOREACH(const edm::InputTag &tag, discs_) {
00149 edm::Handle<reco::PFTauDiscriminator> discHandle;
00150 evt.getByLabel(tag, discHandle);
00151
00152 double result = (*discHandle)[tau];
00153 HistoMap& mymap = histos_[tag.label()];
00154 mymap["plain"]->Fill(result);
00155 mymap["vs_pt"]->Fill(result, tau->pt());
00156 mymap["vs_jetPt"]->Fill(result, tau->jetRef()->pt());
00157 mymap["vs_embedPt"]->Fill(result, tau->alternatLorentzVect().pt());
00158 dynamic_cast<TH3F*>(mymap["vs_pt_jetPt"])->Fill(
00159 result, tau->pt(), tau->jetRef()->pt());
00160 dynamic_cast<TH3F*>(mymap["vs_pt_embedPt"])->Fill(
00161 result, tau->pt(), tau->alternatLorentzVect().pt());
00162 mymap["vs_eta"]->Fill(result, tau->eta());
00163 mymap["vs_dm"]->Fill(result, tau->decayMode());
00164 if (plotPU_ && tau->pt() > pileupPtCut_) {
00165 if (puInfo.isValid())
00166 mymap["vs_truePU"]->Fill(result, puInfo->getPU_NumInteractions());
00167 mymap["vs_recoPU"]->Fill(result, puVertices->size());
00168 }
00169 }
00170 }
00171 }
00172
00173 #include "FWCore/Framework/interface/MakerMacros.h"
00174 DEFINE_FWK_MODULE(RecoTauPlotDiscriminator);