CMS 3D CMS Logo

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/PatCandidates/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/TtEventPartons.h
00028 //
00029 //  output is:
00030 //  * a vector of vectors in the orders defined in TtEventPartons 
00031 //    containing the index of the coresponding jet/lepton in the 
00032 //    input collection of jets/leptons
00033 //  * and a set of vectors with quality parameters of the matching 
00034 //
00035 //  the matching can be performed on an arbitary number of jet 
00036 //  combinations; per default the combination which matches best 
00037 //  according to the quality parameters will be stored; the length
00038 //  of the vectors will be 1 then 
00039 // ----------------------------------------------------------------------
00040 
00041 template <typename C>
00042 class TtJetPartonMatch : public edm::EDProducer {
00043 
00044  private:
00045   
00047   typedef std::vector<pat::Jet> TopJetCollection;
00048   
00049  public:
00050 
00052   explicit TtJetPartonMatch(const edm::ParameterSet&);
00054   ~TtJetPartonMatch();
00056   virtual void produce(edm::Event&, const edm::EventSetup&);
00057 
00058  private:
00059 
00062   int event_;
00064   edm::InputTag jets_;
00067   int maxNJets_;
00070   int maxNComb_;
00072   int algorithm_;
00074   bool useDeltaR_;
00077   bool useMaxDist_;
00080   double maxDist_;
00082   int verbosity_;
00083 };
00084 
00085 template<typename C>
00086 TtJetPartonMatch<C>::TtJetPartonMatch(const edm::ParameterSet& cfg): event_(0),
00087   jets_      (cfg.getParameter<edm::InputTag>("jets"      )),
00088   maxNJets_  (cfg.getParameter<int>          ("maxNJets"  )),
00089   maxNComb_  (cfg.getParameter<int>          ("maxNComb"  )),
00090   algorithm_ (cfg.getParameter<int>          ("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   //  * TtSemiLepEventPartons
00098   //  * TtFullHadEventPartons
00099   //  * TtFullLepEventPartons
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<TopJetCollection> topJets;
00119   evt.getByLabel(jets_, topJets);
00120 
00121   // fill vector of partons in the order of
00122   //  * TtFullLepEventPartons
00123   //  * TtSemiLepEventPartons 
00124   //  * TtFullHadEventPartons
00125   C parts;
00126   std::vector<const reco::Candidate*> partons = parts.vec(*genEvt);
00127 
00128   // prepare vector of jets
00129   std::vector<pat::Jet> jets;
00130   for(unsigned int ij=0; ij<topJets->size(); ++ij) {
00131     // take all jets if maxNJets_ == -1; otherwise use
00132     // maxNJets_ if maxNJets_ is big enough or use same
00133     // number of jets as partons if maxNJets_ < number 
00134     // of partons
00135     if(maxNJets_!=-1) {
00136       if(maxNJets_>=(int)partons.size()) {
00137         if((int)ij==maxNJets_) break;
00138       }
00139       else {
00140         if(ij==partons.size()) break;
00141       }
00142     }
00143     jets.push_back( (*topJets)[ij] );
00144   }
00145 
00146   // do the matching with specified parameters
00147   JetPartonMatching jetPartonMatch(partons, jets, algorithm_, useMaxDist_, useDeltaR_, maxDist_);
00148 
00149   // print some info for each event
00150   // if corresponding verbosity level set
00151   if(verbosity_>0 && event_<verbosity_){
00152     ++event_;
00153     jetPartonMatch.print();
00154   }
00155 
00156   // write 
00157   // * parton match 
00158   // * sumPt 
00159   // * sumDR
00160   // to the event
00161   std::auto_ptr< std::vector<std::vector<int> > > match(new std::vector<std::vector<int> >);
00162   std::auto_ptr< std::vector<double> > sumPt(new std::vector<double>);
00163   std::auto_ptr< std::vector<double> > sumDR(new std::vector<double>);
00164 
00165   for(unsigned int ic=0; ic<jetPartonMatch.getNumberOfAvailableCombinations(); ++ic) {
00166     if((int)ic>=maxNComb_ && maxNComb_>=0) break;
00167     match->push_back( jetPartonMatch.getMatchesForPartons(ic) );
00168     sumPt->push_back( jetPartonMatch.getSumDeltaPt       (ic) );
00169     sumDR->push_back( jetPartonMatch.getSumDeltaR        (ic) );
00170   }
00171   evt.put(match);
00172   evt.put(sumPt, "SumPt");
00173   evt.put(sumDR, "SumDR");
00174 }
00175 
00176 #endif

Generated on Tue Jun 9 17:48:15 2009 for CMSSW by  doxygen 1.5.4