CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/TopQuarkAnalysis/TopTools/plugins/TtJetPartonMatch.h

Go to the documentation of this file.
00001 #ifndef TtJetPartonMatch_h
00002 #define TtJetPartonMatch_h
00003 
00004 #include <memory>
00005 #include <vector>
00006 
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "FWCore/Framework/interface/EDProducer.h"
00009 #include "FWCore/Framework/interface/Frameworkfwd.h"
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011 
00012 #include "DataFormats/JetReco/interface/Jet.h"
00013 #include "AnalysisDataFormats/TopObjects/interface/TtGenEvent.h"
00014 #include "TopQuarkAnalysis/TopTools/interface/JetPartonMatching.h"
00015 
00016 // ----------------------------------------------------------------------
00017 // template class for: 
00018 //
00019 //  * TtFullHadJetPartonMatch
00020 //  * TtSemiLepJetPartonMatch
00021 //  * TtFullLepJetPartonMatch
00022 //
00023 //  the class provides plugins for jet-parton matching corresponding 
00024 //  to the JetPartonMatching class; expected input is one of the 
00025 //  classes in:
00026 //
00027 //  AnalysisDataFormats/TopObjects/interface/TtFullHadEvtPartons.h
00028 //  AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h
00029 //  AnalysisDataFormats/TopObjects/interface/TtFullLepEvtPartons.h
00030 //
00031 //  output is:
00032 //  * a vector of vectors containing the indices of the jets in the
00033 //    input collection matched to the partons in the order defined in
00034 //    the corresponding Tt*EvtPartons class
00035 //  * a set of vectors with quality parameters of the matching 
00036 //
00037 //  the matching can be performed on an arbitrary number of jet 
00038 //  combinations; per default the combination which matches best 
00039 //  according to the quality parameters will be stored; the length
00040 //  of the vectors will be 1 then 
00041 // ----------------------------------------------------------------------
00042 
00043 template <typename C>
00044 class TtJetPartonMatch : public edm::EDProducer {
00045 
00046  public:
00047 
00049   explicit TtJetPartonMatch(const edm::ParameterSet&);
00051   ~TtJetPartonMatch();
00053   virtual void produce(edm::Event&, const edm::EventSetup&);
00054 
00055  private:
00056 
00058   JetPartonMatching::algorithms readAlgorithm(const std::string& str);
00059 
00061   C partons_;
00063   edm::InputTag jets_;
00066   int maxNJets_;
00069   int maxNComb_;
00071   JetPartonMatching::algorithms algorithm_;
00073   bool useDeltaR_;
00076   bool useMaxDist_;
00079   double maxDist_;
00081   int verbosity_;
00082 };
00083 
00084 template<typename C>
00085 TtJetPartonMatch<C>::TtJetPartonMatch(const edm::ParameterSet& cfg):
00086   partons_   (cfg.getParameter<std::vector<std::string> >("partonsToIgnore")),
00087   jets_      (cfg.getParameter<edm::InputTag>            ("jets"           )),
00088   maxNJets_  (cfg.getParameter<int>                      ("maxNJets"       )),
00089   maxNComb_  (cfg.getParameter<int>                      ("maxNComb"       )),
00090   algorithm_ (readAlgorithm(cfg.getParameter<std::string>("algorithm"      ))),
00091   useDeltaR_ (cfg.getParameter<bool>                     ("useDeltaR"      )),
00092   useMaxDist_(cfg.getParameter<bool>                     ("useMaxDist"     )),
00093   maxDist_   (cfg.getParameter<double>                   ("maxDist"        )),
00094   verbosity_ (cfg.getParameter<int>                      ("verbosity"      ))
00095 {
00096   // produces a vector of jet/lepton indices in the order of
00097   //  * TtSemiLepEvtPartons
00098   //  * TtFullHadEvtPartons
00099   //  * TtFullLepEvtPartons
00100   // and vectors of the corresponding quality parameters
00101   produces< std::vector<std::vector<int> > >();
00102   produces< std::vector<double> >("SumPt");
00103   produces< std::vector<double> >("SumDR");
00104 }
00105 
00106 template<typename C>
00107 TtJetPartonMatch<C>::~TtJetPartonMatch()
00108 {
00109 }
00110 
00111 template<typename C>
00112 void
00113 TtJetPartonMatch<C>::produce(edm::Event& evt, const edm::EventSetup& setup)
00114 {
00115   edm::Handle<TtGenEvent> genEvt;
00116   evt.getByLabel("genEvt", genEvt);
00117   
00118   edm::Handle<edm::View<reco::Jet> > topJets;
00119   evt.getByLabel(jets_, topJets);
00120 
00121   // fill vector of partons in the order of
00122   //  * TtFullLepEvtPartons
00123   //  * TtSemiLepEvtPartons 
00124   //  * TtFullHadEvtPartons
00125   std::vector<const reco::Candidate*> partons = partons_.vec(*genEvt);
00126 
00127   // prepare vector of jets
00128   std::vector<const reco::Candidate*> jets;
00129   for(unsigned int ij=0; ij<topJets->size(); ++ij) {
00130     // take all jets if maxNJets_ == -1; otherwise use
00131     // maxNJets_ if maxNJets_ is big enough or use same
00132     // number of jets as partons if maxNJets_ < number 
00133     // of partons
00134     if(maxNJets_!=-1) {
00135       if(maxNJets_>=(int)partons.size()) {
00136         if((int)ij==maxNJets_) break;
00137       }
00138       else {
00139         if(ij==partons.size()) break;
00140       }
00141     }
00142     jets.push_back( (const reco::Candidate*) &(*topJets)[ij] );
00143   }
00144 
00145   // do the matching with specified parameters
00146   JetPartonMatching jetPartonMatch(partons, jets, algorithm_, useMaxDist_, useDeltaR_, maxDist_);
00147 
00148   // print some info for each event
00149   // if corresponding verbosity level set
00150   if(verbosity_>0)
00151     jetPartonMatch.print();
00152 
00153   // write 
00154   // * parton match 
00155   // * sumPt 
00156   // * sumDR
00157   // to the event
00158   std::auto_ptr< std::vector<std::vector<int> > > match(new std::vector<std::vector<int> >);
00159   std::auto_ptr< std::vector<double> > sumPt(new std::vector<double>);
00160   std::auto_ptr< std::vector<double> > sumDR(new std::vector<double>);
00161 
00162   for(unsigned int ic=0; ic<jetPartonMatch.getNumberOfAvailableCombinations(); ++ic) {
00163     if((int)ic>=maxNComb_ && maxNComb_>=0) break;
00164     std::vector<int> matches = jetPartonMatch.getMatchesForPartons(ic);
00165     partons_.expand(matches); // insert dummy indices for partons that were chosen to be ignored
00166     match->push_back( matches );
00167     sumPt->push_back( jetPartonMatch.getSumDeltaPt(ic) );
00168     sumDR->push_back( jetPartonMatch.getSumDeltaR (ic) );
00169   }
00170   evt.put(match);
00171   evt.put(sumPt, "SumPt");
00172   evt.put(sumDR, "SumDR");
00173 }
00174 
00175 template<typename C>
00176 JetPartonMatching::algorithms 
00177 TtJetPartonMatch<C>::readAlgorithm(const std::string& str)
00178 {
00179   if     (str == "totalMinDist"    ) return JetPartonMatching::totalMinDist;
00180   else if(str == "minSumDist"      ) return JetPartonMatching::minSumDist;
00181   else if(str == "ptOrderedMinDist") return JetPartonMatching::ptOrderedMinDist;
00182   else if(str == "unambiguousOnly" ) return JetPartonMatching::unambiguousOnly;
00183   else throw cms::Exception("Configuration")
00184     << "Chosen algorithm is not supported: " << str << "\n";
00185 }
00186 
00187 #endif