CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/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   produces<int>("NumberOfConsideredJets");
00105 }
00106 
00107 template<typename C>
00108 TtJetPartonMatch<C>::~TtJetPartonMatch()
00109 {
00110 }
00111 
00112 template<typename C>
00113 void
00114 TtJetPartonMatch<C>::produce(edm::Event& evt, const edm::EventSetup& setup)
00115 {
00116   // will write 
00117   // * parton match 
00118   // * sumPt 
00119   // * sumDR
00120   // to the event
00121   std::auto_ptr<std::vector<std::vector<int> > > match(new std::vector<std::vector<int> >);
00122   std::auto_ptr<std::vector<double> > sumPt(new std::vector<double>);
00123   std::auto_ptr<std::vector<double> > sumDR(new std::vector<double>);
00124   std::auto_ptr<int> pJetsConsidered(new int);
00125 
00126   // get TtGenEvent and jet collection from the event
00127   edm::Handle<TtGenEvent> genEvt;
00128   evt.getByLabel("genEvt", genEvt);
00129   
00130   edm::Handle<edm::View<reco::Jet> > topJets;
00131   evt.getByLabel(jets_, topJets);
00132 
00133   // fill vector of partons in the order of
00134   //  * TtFullLepEvtPartons
00135   //  * TtSemiLepEvtPartons 
00136   //  * TtFullHadEvtPartons
00137   std::vector<const reco::Candidate*> partons = partons_.vec(*genEvt);
00138 
00139   // prepare vector of jets
00140   std::vector<const reco::Candidate*> jets;
00141   for(unsigned int ij=0; ij<topJets->size(); ++ij) {
00142     // take all jets if maxNJets_ == -1; otherwise use
00143     // maxNJets_ if maxNJets_ is big enough or use same
00144     // number of jets as partons if maxNJets_ < number 
00145     // of partons
00146     if(maxNJets_!=-1) {
00147       if(maxNJets_>=(int)partons.size()) {
00148         if((int)ij==maxNJets_) break;
00149       }
00150       else {
00151         if(ij==partons.size()) break;
00152       }
00153     }
00154     jets.push_back( (const reco::Candidate*) &(*topJets)[ij] );
00155   }
00156   *pJetsConsidered = jets.size();
00157 
00158   // do the matching with specified parameters
00159   JetPartonMatching jetPartonMatch(partons, jets, algorithm_, useMaxDist_, useDeltaR_, maxDist_);
00160 
00161   // print some info for each event
00162   // if corresponding verbosity level set
00163   if(verbosity_>0)
00164     jetPartonMatch.print();
00165 
00166   for(unsigned int ic=0; ic<jetPartonMatch.getNumberOfAvailableCombinations(); ++ic) {
00167     if((int)ic>=maxNComb_ && maxNComb_>=0) break;
00168     std::vector<int> matches = jetPartonMatch.getMatchesForPartons(ic);
00169     partons_.expand(matches); // insert dummy indices for partons that were chosen to be ignored
00170     match->push_back( matches );
00171     sumPt->push_back( jetPartonMatch.getSumDeltaPt(ic) );
00172     sumDR->push_back( jetPartonMatch.getSumDeltaR (ic) );
00173   }
00174   evt.put(match);
00175   evt.put(sumPt, "SumPt");
00176   evt.put(sumDR, "SumDR");
00177   evt.put(pJetsConsidered, "NumberOfConsideredJets");
00178 }
00179 
00180 template<typename C>
00181 JetPartonMatching::algorithms 
00182 TtJetPartonMatch<C>::readAlgorithm(const std::string& str)
00183 {
00184   if     (str == "totalMinDist"    ) return JetPartonMatching::totalMinDist;
00185   else if(str == "minSumDist"      ) return JetPartonMatching::minSumDist;
00186   else if(str == "ptOrderedMinDist") return JetPartonMatching::ptOrderedMinDist;
00187   else if(str == "unambiguousOnly" ) return JetPartonMatching::unambiguousOnly;
00188   else throw cms::Exception("Configuration")
00189     << "Chosen algorithm is not supported: " << str << "\n";
00190 }
00191 
00192 #endif