CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/TopQuarkAnalysis/TopKinFitter/plugins/TtSemiLepKinFitProducer.h

Go to the documentation of this file.
00001 #ifndef TtSemiLepKinFitProducer_h
00002 #define TtSemiLepKinFitProducer_h
00003 
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "FWCore/Framework/interface/EDProducer.h"
00006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00007 
00008 #include "PhysicsTools/JetMCUtils/interface/combination.h"
00009 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h"
00010 #include "TopQuarkAnalysis/TopKinFitter/interface/TtSemiLepKinFitter.h"
00011 
00012 template <typename LeptonCollection>
00013 class TtSemiLepKinFitProducer : public edm::EDProducer {
00014   
00015  public:
00016   
00017   explicit TtSemiLepKinFitProducer(const edm::ParameterSet&);
00018   ~TtSemiLepKinFitProducer();
00019   
00020  private:
00021   // produce
00022   virtual void produce(edm::Event&, const edm::EventSetup&);
00023 
00024   // convert unsigned to Param
00025   TtSemiLepKinFitter::Param param(unsigned);
00026   // convert unsigned to Param
00027   TtSemiLepKinFitter::Constraint constraint(unsigned);
00028   // convert unsigned to Param
00029   std::vector<TtSemiLepKinFitter::Constraint> constraints(std::vector<unsigned>&);
00030   // helper function for b-tagging
00031   bool doBTagging(bool& useBTag_, edm::Handle<std::vector<pat::Jet> >& jets, std::vector<int>& combi,
00032                   std::string& bTagAlgo_, double& minBTagValueBJets_, double& maxBTagValueNonBJets_);
00033 
00034   edm::InputTag jets_;
00035   edm::InputTag leps_;
00036   edm::InputTag mets_;
00037   
00038   edm::InputTag match_;
00040   bool useOnlyMatch_;
00042   std::string bTagAlgo_;
00044   double minBTagValueBJet_;
00046   double maxBTagValueNonBJet_;
00048   bool useBTag_;
00050   int maxNJets_;
00052   int maxNComb_;
00053 
00055   unsigned int maxNrIter_;
00057   double maxDeltaS_;
00059   double maxF_;
00060   unsigned int jetParam_;
00061   unsigned int lepParam_;
00062   unsigned int metParam_;
00064   std::vector<unsigned> constraints_;
00065   double mW_;
00066   double mTop_;
00068   std::vector<double> jetEnergyResolutionScaleFactors_;
00069   std::vector<double> jetEnergyResolutionEtaBinning_;
00071   std::vector<edm::ParameterSet> udscResolutions_;
00072   std::vector<edm::ParameterSet> bResolutions_;
00073   std::vector<edm::ParameterSet> lepResolutions_;
00074   std::vector<edm::ParameterSet> metResolutions_;
00075 
00076   TtSemiLepKinFitter* fitter;
00077 
00078   struct KinFitResult {
00079     int Status;
00080     double Chi2;
00081     double Prob;
00082     pat::Particle HadB;
00083     pat::Particle HadP;
00084     pat::Particle HadQ;
00085     pat::Particle LepB;
00086     pat::Particle LepL;
00087     pat::Particle LepN;
00088     std::vector<int> JetCombi;
00089     bool operator< (const KinFitResult& rhs) { return Chi2 < rhs.Chi2; };
00090   };
00091 };
00092 
00093 template<typename LeptonCollection>
00094 TtSemiLepKinFitProducer<LeptonCollection>::TtSemiLepKinFitProducer(const edm::ParameterSet& cfg):
00095   jets_                    (cfg.getParameter<edm::InputTag>("jets")),
00096   leps_                    (cfg.getParameter<edm::InputTag>("leps")),
00097   mets_                    (cfg.getParameter<edm::InputTag>("mets")),
00098   match_                   (cfg.getParameter<edm::InputTag>("match")),
00099   useOnlyMatch_            (cfg.getParameter<bool>         ("useOnlyMatch"        )),
00100   bTagAlgo_                (cfg.getParameter<std::string>  ("bTagAlgo"            )),
00101   minBTagValueBJet_        (cfg.getParameter<double>       ("minBDiscBJets"       )),
00102   maxBTagValueNonBJet_     (cfg.getParameter<double>       ("maxBDiscLightJets"   )),
00103   useBTag_                 (cfg.getParameter<bool>         ("useBTagging"         )),
00104   maxNJets_                (cfg.getParameter<int>          ("maxNJets"            )),
00105   maxNComb_                (cfg.getParameter<int>          ("maxNComb"            )),
00106   maxNrIter_               (cfg.getParameter<unsigned>     ("maxNrIter"           )),
00107   maxDeltaS_               (cfg.getParameter<double>       ("maxDeltaS"           )),
00108   maxF_                    (cfg.getParameter<double>       ("maxF"                )),
00109   jetParam_                (cfg.getParameter<unsigned>     ("jetParametrisation"  )),
00110   lepParam_                (cfg.getParameter<unsigned>     ("lepParametrisation"  )),
00111   metParam_                (cfg.getParameter<unsigned>     ("metParametrisation"  )),
00112   constraints_             (cfg.getParameter<std::vector<unsigned> >("constraints")),
00113   mW_                      (cfg.getParameter<double>       ("mW"                  )),
00114   mTop_                    (cfg.getParameter<double>       ("mTop"                )),
00115   jetEnergyResolutionScaleFactors_(cfg.getParameter<std::vector<double> >("jetEnergyResolutionScaleFactors")),
00116   jetEnergyResolutionEtaBinning_  (cfg.getParameter<std::vector<double> >("jetEnergyResolutionEtaBinning")),
00117   udscResolutions_(0), bResolutions_(0), lepResolutions_(0), metResolutions_(0)
00118 {
00119   if(cfg.exists("udscResolutions") && cfg.exists("bResolutions") && cfg.exists("lepResolutions") && cfg.exists("metResolutions")){
00120     udscResolutions_ = cfg.getParameter<std::vector<edm::ParameterSet> >("udscResolutions");
00121     bResolutions_    = cfg.getParameter<std::vector<edm::ParameterSet> >("bResolutions"   );
00122     lepResolutions_  = cfg.getParameter<std::vector<edm::ParameterSet> >("lepResolutions" );
00123     metResolutions_  = cfg.getParameter<std::vector<edm::ParameterSet> >("metResolutions" );
00124   }
00125   else if(cfg.exists("udscResolutions") || cfg.exists("bResolutions") || cfg.exists("lepResolutions") || cfg.exists("metResolutions") ){
00126     throw cms::Exception("Configuration") << "Parameters 'udscResolutions', 'bResolutions', 'lepResolutions', 'metResolutions' should be used together.\n";
00127   }
00128 
00129   fitter = new TtSemiLepKinFitter(param(jetParam_), param(lepParam_), param(metParam_), maxNrIter_, maxDeltaS_, maxF_,
00130                                   constraints(constraints_), mW_, mTop_, &udscResolutions_, &bResolutions_, &lepResolutions_, &metResolutions_,
00131                                   &jetEnergyResolutionScaleFactors_, &jetEnergyResolutionEtaBinning_);
00132 
00133   produces< std::vector<pat::Particle> >("PartonsHadP");
00134   produces< std::vector<pat::Particle> >("PartonsHadQ");
00135   produces< std::vector<pat::Particle> >("PartonsHadB");
00136   produces< std::vector<pat::Particle> >("PartonsLepB");
00137   produces< std::vector<pat::Particle> >("Leptons");
00138   produces< std::vector<pat::Particle> >("Neutrinos");
00139 
00140   produces< std::vector<std::vector<int> > >();
00141   produces< std::vector<double> >("Chi2");
00142   produces< std::vector<double> >("Prob");
00143   produces< std::vector<int> >("Status");
00144 
00145   produces<int>("NumberOfConsideredJets");
00146 }
00147 
00148 template<typename LeptonCollection>
00149 TtSemiLepKinFitProducer<LeptonCollection>::~TtSemiLepKinFitProducer()
00150 {
00151   delete fitter;
00152 }
00153 
00154 template<typename LeptonCollection>
00155 bool TtSemiLepKinFitProducer<LeptonCollection>::doBTagging(bool& useBTag_, edm::Handle<std::vector<pat::Jet> >& jets, std::vector<int>& combi,
00156                                                            std::string& bTagAlgo_, double& minBTagValueBJet_, double& maxBTagValueNonBJet_){
00157   
00158   if( !useBTag_ ) {
00159     return true;
00160   }
00161   if( useBTag_ &&
00162       (*jets)[combi[TtSemiLepEvtPartons::HadB     ]].bDiscriminator(bTagAlgo_) >= minBTagValueBJet_ &&
00163       (*jets)[combi[TtSemiLepEvtPartons::LepB     ]].bDiscriminator(bTagAlgo_) >= minBTagValueBJet_ &&
00164       (*jets)[combi[TtSemiLepEvtPartons::LightQ   ]].bDiscriminator(bTagAlgo_) <  maxBTagValueNonBJet_ &&
00165       (*jets)[combi[TtSemiLepEvtPartons::LightQBar]].bDiscriminator(bTagAlgo_) <  maxBTagValueNonBJet_ ) {
00166     return true;
00167   }
00168   else{
00169     return false;
00170   }
00171 }
00172 
00173 template<typename LeptonCollection>
00174 void TtSemiLepKinFitProducer<LeptonCollection>::produce(edm::Event& evt, const edm::EventSetup& setup)
00175 {
00176   std::auto_ptr< std::vector<pat::Particle> > pPartonsHadP( new std::vector<pat::Particle> );
00177   std::auto_ptr< std::vector<pat::Particle> > pPartonsHadQ( new std::vector<pat::Particle> );
00178   std::auto_ptr< std::vector<pat::Particle> > pPartonsHadB( new std::vector<pat::Particle> );
00179   std::auto_ptr< std::vector<pat::Particle> > pPartonsLepB( new std::vector<pat::Particle> );
00180   std::auto_ptr< std::vector<pat::Particle> > pLeptons    ( new std::vector<pat::Particle> );
00181   std::auto_ptr< std::vector<pat::Particle> > pNeutrinos  ( new std::vector<pat::Particle> );
00182 
00183   std::auto_ptr< std::vector<std::vector<int> > > pCombi ( new std::vector<std::vector<int> > );
00184   std::auto_ptr< std::vector<double>            > pChi2  ( new std::vector<double> );
00185   std::auto_ptr< std::vector<double>            > pProb  ( new std::vector<double> );
00186   std::auto_ptr< std::vector<int>               > pStatus( new std::vector<int> );
00187 
00188   std::auto_ptr<int> pJetsConsidered(new int);
00189 
00190   edm::Handle<std::vector<pat::Jet> > jets;
00191   evt.getByLabel(jets_, jets);
00192 
00193   edm::Handle<std::vector<pat::MET> > mets;
00194   evt.getByLabel(mets_, mets);
00195 
00196   edm::Handle<LeptonCollection> leps;
00197   evt.getByLabel(leps_, leps);
00198 
00199   const unsigned int nPartons = 4;
00200 
00201   std::vector<int> match;
00202   bool invalidMatch = false;
00203   if(useOnlyMatch_) {
00204     *pJetsConsidered = nPartons;
00205     edm::Handle<std::vector<std::vector<int> > > matchHandle;
00206     evt.getByLabel(match_, matchHandle);
00207     match = *(matchHandle->begin());
00208     // check if match is valid
00209     if(match.size()!=nPartons) invalidMatch=true;
00210     else {
00211       for(unsigned int idx=0; idx<match.size(); ++idx) {
00212         if(match[idx]<0 || match[idx]>=(int)jets->size()) {
00213           invalidMatch=true;
00214           break;
00215         }
00216       }
00217     }
00218   }
00219 
00220   // -----------------------------------------------------
00221   // skip events with no appropriate lepton candidate in
00222   // or empty MET or less jets than partons or invalid match
00223   // -----------------------------------------------------
00224 
00225   if( leps->empty() || mets->empty() || jets->size()<nPartons || invalidMatch ) {
00226     // the kinFit getters return empty objects here
00227     pPartonsHadP->push_back( fitter->fittedHadP()     );
00228     pPartonsHadQ->push_back( fitter->fittedHadQ()     );
00229     pPartonsHadB->push_back( fitter->fittedHadB()     );
00230     pPartonsLepB->push_back( fitter->fittedLepB()     );
00231     pLeptons    ->push_back( fitter->fittedLepton()   );
00232     pNeutrinos  ->push_back( fitter->fittedNeutrino() );
00233     // indices referring to the jet combination
00234     std::vector<int> invalidCombi;
00235     for(unsigned int i = 0; i < nPartons; ++i) 
00236       invalidCombi.push_back( -1 );
00237     pCombi->push_back( invalidCombi );
00238     // chi2
00239     pChi2->push_back( -1. );
00240     // chi2 probability
00241     pProb->push_back( -1. );
00242     // status of the fitter
00243     pStatus->push_back( -1 );
00244     // number of jets
00245     *pJetsConsidered = jets->size();
00246     // feed out all products
00247     evt.put(pCombi);
00248     evt.put(pPartonsHadP, "PartonsHadP");
00249     evt.put(pPartonsHadQ, "PartonsHadQ");
00250     evt.put(pPartonsHadB, "PartonsHadB");
00251     evt.put(pPartonsLepB, "PartonsLepB");
00252     evt.put(pLeptons    , "Leptons"    );
00253     evt.put(pNeutrinos  , "Neutrinos"  );
00254     evt.put(pChi2       , "Chi2"       );
00255     evt.put(pProb       , "Prob"       );
00256     evt.put(pStatus     , "Status"     );
00257     evt.put(pJetsConsidered, "NumberOfConsideredJets");
00258     return;
00259   }
00260 
00261   // -----------------------------------------------------
00262   // analyze different jet combinations using the KinFitter
00263   // (or only a given jet combination if useOnlyMatch=true)
00264   // -----------------------------------------------------
00265   
00266   std::vector<int> jetIndices;
00267   if(!useOnlyMatch_) {
00268     for(unsigned int i=0; i<jets->size(); ++i){
00269       if(maxNJets_ >= (int) nPartons && maxNJets_ == (int) i) {
00270         *pJetsConsidered = i;
00271         break;
00272       }
00273       jetIndices.push_back(i);
00274     }
00275   }
00276   
00277   std::vector<int> combi;
00278   for(unsigned int i=0; i<nPartons; ++i) {
00279     if(useOnlyMatch_) combi.push_back( match[i] );
00280     else combi.push_back(i);
00281   }
00282 
00283   std::list<KinFitResult> FitResultList;
00284 
00285   do{
00286     for(int cnt = 0; cnt < TMath::Factorial( combi.size() ); ++cnt){
00287       // take into account indistinguishability of the two jets from the hadr. W decay,
00288       // reduces combinatorics by a factor of 2
00289       if( (combi[TtSemiLepEvtPartons::LightQ] < combi[TtSemiLepEvtPartons::LightQBar]
00290          || useOnlyMatch_ ) && doBTagging(useBTag_, jets, combi, bTagAlgo_, minBTagValueBJet_, maxBTagValueNonBJet_) ){
00291         
00292         std::vector<pat::Jet> jetCombi;
00293         jetCombi.resize(nPartons);
00294         jetCombi[TtSemiLepEvtPartons::LightQ   ] = (*jets)[combi[TtSemiLepEvtPartons::LightQ   ]];
00295         jetCombi[TtSemiLepEvtPartons::LightQBar] = (*jets)[combi[TtSemiLepEvtPartons::LightQBar]];
00296         jetCombi[TtSemiLepEvtPartons::HadB     ] = (*jets)[combi[TtSemiLepEvtPartons::HadB     ]];
00297         jetCombi[TtSemiLepEvtPartons::LepB     ] = (*jets)[combi[TtSemiLepEvtPartons::LepB     ]];
00298 
00299         // do the kinematic fit
00300         const int status = fitter->fit(jetCombi, (*leps)[0], (*mets)[0]);
00301 
00302         if( status == 0 ) { // only take into account converged fits
00303           KinFitResult result;
00304           result.Status = status;
00305           result.Chi2 = fitter->fitS();
00306           result.Prob = fitter->fitProb();
00307           result.HadB = fitter->fittedHadB();
00308           result.HadP = fitter->fittedHadP();
00309           result.HadQ = fitter->fittedHadQ();
00310           result.LepB = fitter->fittedLepB();
00311           result.LepL = fitter->fittedLepton();
00312           result.LepN = fitter->fittedNeutrino();
00313           result.JetCombi = combi;
00314 
00315           FitResultList.push_back(result);
00316         }
00317 
00318       }
00319       if(useOnlyMatch_) break; // don't go through combinatorics if useOnlyMatch was chosen
00320       next_permutation( combi.begin(), combi.end() );
00321     }
00322     if(useOnlyMatch_) break; // don't go through combinatorics if useOnlyMatch was chosen
00323   }
00324   while(stdcomb::next_combination( jetIndices.begin(), jetIndices.end(), combi.begin(), combi.end() ));
00325 
00326   // sort results w.r.t. chi2 values
00327   FitResultList.sort();
00328   
00329   // -----------------------------------------------------
00330   // feed out result
00331   // starting with the JetComb having the smallest chi2
00332   // -----------------------------------------------------
00333 
00334   if( (unsigned)FitResultList.size() < 1 ) { // in case no fit results were stored in the list (all fits aborted)
00335     pPartonsHadP->push_back( fitter->fittedHadP()     );
00336     pPartonsHadQ->push_back( fitter->fittedHadQ()     );
00337     pPartonsHadB->push_back( fitter->fittedHadB()     );
00338     pPartonsLepB->push_back( fitter->fittedLepB()     );
00339     pLeptons    ->push_back( fitter->fittedLepton()   );
00340     pNeutrinos  ->push_back( fitter->fittedNeutrino() );
00341     // indices referring to the jet combination
00342     std::vector<int> invalidCombi;
00343     for(unsigned int i = 0; i < nPartons; ++i) 
00344       invalidCombi.push_back( -1 );
00345     pCombi->push_back( invalidCombi );
00346     // chi2
00347     pChi2->push_back( -1. );
00348     // chi2 probability
00349     pProb->push_back( -1. );
00350     // status of the fitter
00351     pStatus->push_back( -1 );
00352   }
00353   else {
00354     unsigned int iComb = 0;
00355     for(typename std::list<KinFitResult>::const_iterator result = FitResultList.begin(); result != FitResultList.end(); ++result) {
00356       if(maxNComb_ >= 1 && iComb == (unsigned int) maxNComb_) break;
00357       iComb++;
00358       // partons
00359       pPartonsHadP->push_back( result->HadP );
00360       pPartonsHadQ->push_back( result->HadQ );
00361       pPartonsHadB->push_back( result->HadB );
00362       pPartonsLepB->push_back( result->LepB );
00363       // lepton
00364       pLeptons->push_back( result->LepL );
00365       // neutrino
00366       pNeutrinos->push_back( result->LepN );
00367       // indices referring to the jet combination
00368       pCombi->push_back( result->JetCombi );
00369       // chi2
00370       pChi2->push_back( result->Chi2 );
00371       // chi2 probability
00372       pProb->push_back( result->Prob );
00373       // status of the fitter
00374       pStatus->push_back( result->Status );
00375     }
00376   }
00377   evt.put(pCombi);
00378   evt.put(pPartonsHadP, "PartonsHadP");
00379   evt.put(pPartonsHadQ, "PartonsHadQ");
00380   evt.put(pPartonsHadB, "PartonsHadB");
00381   evt.put(pPartonsLepB, "PartonsLepB");
00382   evt.put(pLeptons    , "Leptons"    );
00383   evt.put(pNeutrinos  , "Neutrinos"  );
00384   evt.put(pChi2       , "Chi2"       );
00385   evt.put(pProb       , "Prob"       );
00386   evt.put(pStatus     , "Status"     );
00387   evt.put(pJetsConsidered, "NumberOfConsideredJets");
00388 }
00389  
00390 template<typename LeptonCollection>
00391 TtSemiLepKinFitter::Param TtSemiLepKinFitProducer<LeptonCollection>::param(unsigned val) 
00392 {
00393   TtSemiLepKinFitter::Param result;
00394   switch(val){
00395   case TtSemiLepKinFitter::kEMom       : result=TtSemiLepKinFitter::kEMom;       break;
00396   case TtSemiLepKinFitter::kEtEtaPhi   : result=TtSemiLepKinFitter::kEtEtaPhi;   break;
00397   case TtSemiLepKinFitter::kEtThetaPhi : result=TtSemiLepKinFitter::kEtThetaPhi; break;
00398   default: 
00399     throw cms::Exception("Configuration") 
00400       << "Chosen jet parametrization is not supported: " << val << "\n";
00401     break;
00402   }
00403   return result;
00404 } 
00405 
00406 template<typename LeptonCollection>
00407 TtSemiLepKinFitter::Constraint TtSemiLepKinFitProducer<LeptonCollection>::constraint(unsigned val) 
00408 {
00409   TtSemiLepKinFitter::Constraint result;
00410   switch(val){
00411   case TtSemiLepKinFitter::kWHadMass       : result=TtSemiLepKinFitter::kWHadMass;       break;
00412   case TtSemiLepKinFitter::kWLepMass       : result=TtSemiLepKinFitter::kWLepMass;       break;
00413   case TtSemiLepKinFitter::kTopHadMass     : result=TtSemiLepKinFitter::kTopHadMass;     break;
00414   case TtSemiLepKinFitter::kTopLepMass     : result=TtSemiLepKinFitter::kTopLepMass;     break;
00415   case TtSemiLepKinFitter::kNeutrinoMass   : result=TtSemiLepKinFitter::kNeutrinoMass;   break;
00416   case TtSemiLepKinFitter::kEqualTopMasses : result=TtSemiLepKinFitter::kEqualTopMasses; break;
00417   case TtSemiLepKinFitter::kSumPt          : result=TtSemiLepKinFitter::kSumPt;          break;
00418   default: 
00419     throw cms::Exception("Configuration") 
00420       << "Chosen fit constraint is not supported: " << val << "\n";
00421     break;
00422   }
00423   return result;
00424 } 
00425 
00426 template<typename LeptonCollection>
00427 std::vector<TtSemiLepKinFitter::Constraint> TtSemiLepKinFitProducer<LeptonCollection>::constraints(std::vector<unsigned>& val)
00428 {
00429   std::vector<TtSemiLepKinFitter::Constraint> result;
00430   for(unsigned i=0; i<val.size(); ++i){
00431     result.push_back(constraint(val[i]));
00432   }
00433   return result; 
00434 }
00435 
00436 #endif