CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/TopQuarkAnalysis/TopEventProducers/interface/TtEvtBuilder.h

Go to the documentation of this file.
00001 #ifndef TtEvtBuilder_h
00002 #define TtEvtBuilder_h
00003 
00004 #include <vector>
00005 
00006 #include "FWCore/Framework/interface/Event.h"
00007 #include "FWCore/Framework/interface/EDProducer.h"
00008 #include "FWCore/Framework/interface/Frameworkfwd.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 #include "AnalysisDataFormats/TopObjects/interface/TtEvent.h"
00011 
00028 template <typename C>
00029 class TtEvtBuilder : public edm::EDProducer {
00030 
00031  public:
00032 
00034   explicit TtEvtBuilder(const edm::ParameterSet&);
00036   ~TtEvtBuilder(){};
00037   
00038  private:
00039 
00042   virtual void produce(edm::Event&, const edm::EventSetup&);
00044   virtual void fillSpecific(C&, const edm::Event&);
00045 
00046  private:
00047 
00049   int verbosity_;
00051   std::vector<std::string> hyps_;
00053   edm::InputTag genEvt_;
00056   int decayChnTop1_;
00057   int decayChnTop2_;
00058 
00061   edm::ParameterSet kinFit_;
00062   edm::InputTag fitChi2_;
00063   edm::InputTag fitProb_;
00066   edm::ParameterSet kinSolution_;
00067   edm::InputTag solWeight_;  
00068   edm::InputTag wrongCharge_;   
00071   edm::ParameterSet genMatch_;
00072   edm::InputTag sumPt_;
00073   edm::InputTag sumDR_;
00076   edm::ParameterSet mvaDisc_;
00077   edm::InputTag meth_;
00078   edm::InputTag disc_;
00079 };
00080 
00081 template <typename C>
00082 TtEvtBuilder<C>::TtEvtBuilder(const edm::ParameterSet& cfg) :
00083   verbosity_   (cfg.getParameter<int>                      ("verbosity"    )),
00084   hyps_        (cfg.getParameter<std::vector<std::string> >("hypotheses"   )),
00085   genEvt_      (cfg.getParameter<edm::InputTag>            ("genEvent"     )),
00086   decayChnTop1_(cfg.getParameter<int>                      ("decayChannel1")),
00087   decayChnTop2_(cfg.getParameter<int>                      ("decayChannel2"))
00088 {
00089   // parameter subsets for kKinFit
00090   if( cfg.exists("kinFit") ) {
00091     kinFit_  = cfg.getParameter<edm::ParameterSet>("kinFit");
00092     fitChi2_ = kinFit_.getParameter<edm::InputTag>("chi2");
00093     fitProb_ = kinFit_.getParameter<edm::InputTag>("prob");
00094   }
00095   // parameter subsets for kKinSolution
00096   if( cfg.exists("kinSolution") ) {
00097     kinSolution_  = cfg.getParameter<edm::ParameterSet>("kinSolution");
00098     solWeight_    = kinSolution_.getParameter<edm::InputTag>("solWeight");
00099     wrongCharge_  = kinSolution_.getParameter<edm::InputTag>("wrongCharge");
00100   }  
00101   // parameter subsets for kGenMatch
00102   if( cfg.exists("genMatch") ) {
00103     genMatch_ = cfg.getParameter<edm::ParameterSet>("genMatch");
00104     sumPt_    = genMatch_.getParameter<edm::InputTag>("sumPt");
00105     sumDR_    = genMatch_.getParameter<edm::InputTag>("sumDR");
00106   }
00107   // parameter subsets for kMvaDisc
00108   if( cfg.exists("mvaDisc") ) {
00109     mvaDisc_ = cfg.getParameter<edm::ParameterSet>("mvaDisc");
00110     meth_    = mvaDisc_.getParameter<edm::InputTag>("meth");
00111     disc_    = mvaDisc_.getParameter<edm::InputTag>("disc");
00112   }
00113   // produces a TtEventEvent for:
00114   //  * TtSemiLeptonicEvent 
00115   //  * TtFullLeptonicEvent
00116   //  * TtFullHadronicEvent
00117   // from hypotheses and associated extra information
00118   produces<C>();
00119 }
00120 
00121 template <typename C>
00122 void
00123 TtEvtBuilder<C>::produce(edm::Event& evt, const edm::EventSetup& setup)
00124 {
00125   C ttEvent;
00126 
00127   // set leptonic decay channels
00128   ttEvent.setLepDecays( WDecay::LeptonType(decayChnTop1_), WDecay::LeptonType(decayChnTop2_) );
00129 
00130   // set genEvent (if available)
00131   edm::Handle<TtGenEvent> genEvt;
00132   if( evt.getByLabel(genEvt_, genEvt) )
00133     ttEvent.setGenEvent(genEvt);
00134 
00135   // add event hypotheses for all given 
00136   // hypothesis classes to the TtEvent
00137   typedef std::vector<std::string>::const_iterator EventHypo;
00138   for(EventHypo h=hyps_.begin(); h!=hyps_.end(); ++h){
00139     edm::Handle<int> key; 
00140     evt.getByLabel(*h, "Key", key);
00141 
00142     edm::Handle<std::vector<TtEvent::HypoCombPair> > hypMatchVec; 
00143     evt.getByLabel(*h, hypMatchVec);
00144 
00145     typedef std::vector<TtEvent::HypoCombPair>::const_iterator HypMatch;
00146     for(HypMatch hm=hypMatchVec->begin(); hm != hypMatchVec->end(); ++hm){
00147       ttEvent.addEventHypo((TtEvent::HypoClassKey&)*key, *hm);
00148     }
00149   }
00150 
00151   // set kKinFit extras
00152   if( ttEvent.isHypoAvailable(TtEvent::kKinFit) ) {
00153     edm::Handle<std::vector<double> > fitChi2;
00154     evt.getByLabel(fitChi2_, fitChi2);
00155     ttEvent.setFitChi2( *fitChi2 );
00156     
00157     edm::Handle<std::vector<double> > fitProb;
00158     evt.getByLabel(fitProb_, fitProb);
00159     ttEvent.setFitProb( *fitProb );
00160   }
00161 
00162   // set kGenMatch extras
00163   if( ttEvent.isHypoAvailable(TtEvent::kGenMatch) ) {
00164     edm::Handle<std::vector<double> > sumPt;
00165     evt.getByLabel(sumPt_, sumPt);
00166     ttEvent.setGenMatchSumPt( *sumPt );
00167 
00168     edm::Handle<std::vector<double> > sumDR;
00169     evt.getByLabel(sumDR_, sumDR);
00170     ttEvent.setGenMatchSumDR( *sumDR );
00171   }
00172 
00173   // set kMvaDisc extras
00174   if( ttEvent.isHypoAvailable(TtEvent::kMVADisc) ) {
00175     edm::Handle<std::string> meth;
00176     evt.getByLabel(meth_, meth);
00177     ttEvent.setMvaMethod( *meth );
00178 
00179     edm::Handle<std::vector<double> > disc;
00180     evt.getByLabel(disc_, disc);
00181     ttEvent.setMvaDiscriminators( *disc );
00182   }
00183 
00184   // fill data members that are decay-channel specific
00185   fillSpecific(ttEvent, evt);
00186 
00187   // print summary via MessageLogger if verbosity_>0
00188   ttEvent.print(verbosity_);
00189 
00190   // write object into the edm::Event
00191   std::auto_ptr<C> pOut(new C);
00192   *pOut=ttEvent;
00193   evt.put(pOut);
00194 }
00195 
00196 template <>
00197 void TtEvtBuilder<TtFullHadronicEvent>::fillSpecific(TtFullHadronicEvent& ttEvent, const edm::Event& evt)
00198 {
00199 }
00200 
00201 template <>
00202 void TtEvtBuilder<TtFullLeptonicEvent>::fillSpecific(TtFullLeptonicEvent& ttEvent, const edm::Event& evt)
00203 {
00204 
00205   // set kKinSolution extras  
00206   if( ttEvent.isHypoAvailable(TtEvent::kKinSolution) ) {
00207     edm::Handle<std::vector<double> > solWeight;
00208     evt.getByLabel(solWeight_, solWeight);
00209     ttEvent.setSolWeight( *solWeight );
00210     
00211     edm::Handle<bool> wrongCharge;
00212     evt.getByLabel(wrongCharge_, wrongCharge);
00213     ttEvent.setWrongCharge( *wrongCharge );   
00214   }
00215 
00216 }
00217 
00218 template <>
00219 void TtEvtBuilder<TtSemiLeptonicEvent>::fillSpecific(TtSemiLeptonicEvent& ttEvent, const edm::Event& evt)
00220 {
00221 
00222   // set number of real neutrino solutions for all hypotheses
00223   typedef std::vector<std::string>::const_iterator EventHypo;
00224   for(EventHypo h=hyps_.begin(); h!=hyps_.end(); ++h){
00225     edm::Handle<int> key; 
00226     evt.getByLabel(*h, "Key", key);
00227 
00228     edm::Handle<int> numberOfRealNeutrinoSolutions;
00229     evt.getByLabel(*h, "NumberOfRealNeutrinoSolutions", numberOfRealNeutrinoSolutions);
00230 
00231     ttEvent.setNumberOfRealNeutrinoSolutions((TtEvent::HypoClassKey&)*key, *numberOfRealNeutrinoSolutions);
00232   }
00233 
00234 }
00235 
00236 #endif