CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoTauTag/TauTagTools/plugins/RecoTauMVATrainer.cc

Go to the documentation of this file.
00001 /*
00002  * RecoTauMVATrainer
00003  *
00004  * Pass either signal or background events (with an option weight) to the
00005  * MVATrainer interface.
00006  *
00007  * Author: Evan K. Friis
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     typedef std::vector<edm::InputTag> VInputTag;
00033     reco::tau::RecoTauMVAHelper mva_;
00034     edm::InputTag signalSrc_;
00035     edm::InputTag backgroundSrc_;
00036     bool applyWeights_;
00037     edm::InputTag signalWeightsSrc_;
00038     edm::InputTag backgroundWeightsSrc_;
00039     VInputTag eventWeights_;
00040 };
00041 
00042 RecoTauMVATrainer::RecoTauMVATrainer(const edm::ParameterSet &pset)
00043   : mva_(pset.getParameter<std::string>("computerName"),
00044          pset.getParameter<std::string>("dbLabel"),
00045          pset.getParameter<edm::ParameterSet>("discriminantOptions")),
00046     signalSrc_(pset.getParameter<edm::InputTag>("signalSrc")),
00047     backgroundSrc_(pset.getParameter<edm::InputTag>("backgroundSrc")) {
00048       // Check if we want to apply weights
00049       applyWeights_ = false;
00050       if (pset.exists("signalWeights")) {
00051         applyWeights_ = true;
00052         signalWeightsSrc_ = pset.getParameter<edm::InputTag>("signalWeights");
00053         backgroundWeightsSrc_ = pset.getParameter<edm::InputTag>("backgroundWeights");
00054       }
00055       if (pset.exists("eventWeights")) {
00056         eventWeights_ = pset.getParameter<VInputTag>("eventWeights");
00057       }
00058     }
00059 
00060 namespace {
00061 
00062 // Upload a view to the MVA
00063 void uploadTrainingData(reco::tau::RecoTauMVAHelper *helper,
00064                         const edm::Handle<reco::CandidateView>& taus,
00065                         const edm::Handle<reco::PFTauDiscriminator>& weights,
00066                         bool isSignal, double eventWeight) {
00067   // Convert to a vector of refs
00068   reco::PFTauRefVector tauRefs =
00069       reco::tau::castView<reco::PFTauRefVector>(taus);
00070   // Loop over our taus and pass each to the MVA interface
00071   BOOST_FOREACH(reco::PFTauRef tau, tauRefs) {
00072     // Lookup the weight if desired
00073     double weight = (weights.isValid()) ? (*weights)[tau] : 1.0;
00074     helper->train(tau, isSignal, weight*eventWeight);
00075   }
00076 }
00077 
00078 }
00079 
00080 
00081 void RecoTauMVATrainer::analyze(const edm::Event &evt,
00082                                 const edm::EventSetup &es) {
00083   // Get a view to our taus
00084   edm::Handle<reco::CandidateView> signal;
00085   edm::Handle<reco::CandidateView> background;
00086 
00087   bool signalExists = false;
00088   evt.getByLabel(signalSrc_, signal);
00089   if (signal.isValid() && signal->size())
00090     signalExists = true;
00091 
00092   bool backgroundExists = false;
00093   evt.getByLabel(backgroundSrc_, background);
00094   if (background.isValid() && background->size())
00095     backgroundExists = true;
00096 
00097   // Check if we have anything to do
00098   bool somethingToDo = signalExists || backgroundExists;
00099   if (!somethingToDo)
00100     return;
00101 
00102   // Make sure the MVA is up to date from the DB
00103   mva_.setEvent(evt, es);
00104 
00105   // Get event weights if specified
00106   double eventWeight = 1.0;
00107   BOOST_FOREACH(const edm::InputTag& weightTag, eventWeights_) {
00108     edm::Handle<double> weightHandle;
00109     evt.getByLabel(weightTag, weightHandle);
00110     if (weightHandle.isValid())
00111       eventWeight *= *weightHandle;
00112   }
00113 
00114   // Get weights if desired
00115   edm::Handle<reco::PFTauDiscriminator> signalWeights;
00116   edm::Handle<reco::PFTauDiscriminator> backgroundWeights;
00117   if (applyWeights_ && signalExists)
00118     evt.getByLabel(signalWeightsSrc_, signalWeights);
00119   if (applyWeights_ && backgroundExists)
00120     evt.getByLabel(backgroundWeightsSrc_, backgroundWeights);
00121 
00122   if (signalExists)
00123     uploadTrainingData(&mva_, signal, signalWeights, true, eventWeight);
00124   if (backgroundExists)
00125     uploadTrainingData(&mva_, background, backgroundWeights,
00126         false, eventWeight);
00127 }
00128 
00129 #include "FWCore/Framework/interface/MakerMacros.h"
00130 DEFINE_FWK_MODULE(RecoTauMVATrainer);