Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include <boost/foreach.hpp>
00011 #include <string>
00012
00013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "FWCore/Framework/interface/EDAnalyzer.h"
00017
00018 #include "RecoTauTag/RecoTau/interface/RecoTauMVAHelper.h"
00019 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00020
00021 #include "DataFormats/TauReco/interface/PFTauFwd.h"
00022 #include "DataFormats/TauReco/interface/PFTau.h"
00023 #include "DataFormats/TauReco/interface/PFTauDiscriminator.h"
00024 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00025
00026 class RecoTauMVATrainer : public edm::EDAnalyzer {
00027 public:
00028 explicit RecoTauMVATrainer(const edm::ParameterSet &pset);
00029 virtual ~RecoTauMVATrainer() {};
00030 virtual void analyze(const edm::Event &evt, const edm::EventSetup &es);
00031 private:
00032 reco::tau::RecoTauMVAHelper mva_;
00033 edm::InputTag signalSrc_;
00034 edm::InputTag backgroundSrc_;
00035 bool applyWeights_;
00036 edm::InputTag signalWeightsSrc_;
00037 edm::InputTag backgroundWeightsSrc_;
00038 };
00039
00040 RecoTauMVATrainer::RecoTauMVATrainer(const edm::ParameterSet &pset)
00041 : mva_(pset.getParameter<std::string>("computerName"),
00042 pset.getParameter<std::string>("dbLabel"),
00043 pset.getParameter<edm::ParameterSet>("discriminantOptions")),
00044 signalSrc_(pset.getParameter<edm::InputTag>("signalSrc")),
00045 backgroundSrc_(pset.getParameter<edm::InputTag>("backgroundSrc")) {
00046
00047 applyWeights_ = false;
00048 if (pset.exists("signalWeights")) {
00049 applyWeights_ = true;
00050 signalWeightsSrc_ = pset.getParameter<edm::InputTag>("signalWeights");
00051 backgroundWeightsSrc_ = pset.getParameter<edm::InputTag>("backgroundWeights");
00052 }
00053 }
00054
00055 namespace {
00056
00057
00058 void uploadTrainingData(reco::tau::RecoTauMVAHelper *helper,
00059 const edm::Handle<reco::CandidateView>& taus,
00060 const edm::Handle<reco::PFTauDiscriminator>& weights,
00061 bool isSignal) {
00062
00063 reco::PFTauRefVector tauRefs =
00064 reco::tau::castView<reco::PFTauRefVector>(taus);
00065
00066 BOOST_FOREACH(reco::PFTauRef tau, tauRefs) {
00067
00068 double weight = (weights.isValid()) ? (*weights)[tau] : 1.0;
00069 helper->train(tau, isSignal, weight);
00070 }
00071 }
00072
00073 }
00074
00075
00076 void RecoTauMVATrainer::analyze(const edm::Event &evt,
00077 const edm::EventSetup &es) {
00078
00079 mva_.setEvent(evt, es);
00080
00081
00082 edm::Handle<reco::CandidateView> signal;
00083 edm::Handle<reco::CandidateView> background;
00084
00085 bool signalExists = true;
00086 try {
00087 evt.getByLabel(signalSrc_, signal);
00088 if (!signal.isValid())
00089 signalExists = false;
00090 } catch(...) {
00091 signalExists = false;
00092 }
00093
00094 bool backgroundExists = true;
00095 try {
00096 evt.getByLabel(backgroundSrc_, background);
00097 if (!background.isValid())
00098 backgroundExists = false;
00099 } catch(...) {
00100 backgroundExists = false;
00101 }
00102
00103
00104 edm::Handle<reco::PFTauDiscriminator> signalWeights;
00105 edm::Handle<reco::PFTauDiscriminator> backgroundWeights;
00106 if (applyWeights_ && signalExists)
00107 evt.getByLabel(signalWeightsSrc_, signalWeights);
00108 if (applyWeights_ && backgroundExists)
00109 evt.getByLabel(backgroundWeightsSrc_, backgroundWeights);
00110
00111 if (signalExists)
00112 uploadTrainingData(&mva_, signal, signalWeights, true);
00113 if (backgroundExists)
00114 uploadTrainingData(&mva_, background, backgroundWeights, false);
00115 }
00116
00117 #include "FWCore/Framework/interface/MakerMacros.h"
00118 DEFINE_FWK_MODULE(RecoTauMVATrainer);