CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/TopQuarkAnalysis/TopKinFitter/plugins/TtFullHadKinFitProducer.cc

Go to the documentation of this file.
00001 #include "TopQuarkAnalysis/TopKinFitter/plugins/TtFullHadKinFitProducer.h"
00002 
00003 static const unsigned int nPartons=6;
00004 
00006 TtFullHadKinFitProducer::TtFullHadKinFitProducer(const edm::ParameterSet& cfg):
00007   jets_                       (cfg.getParameter<edm::InputTag>("jets")),
00008   match_                      (cfg.getParameter<edm::InputTag>("match")),
00009   useOnlyMatch_               (cfg.getParameter<bool>("useOnlyMatch")),
00010   bTagAlgo_                   (cfg.getParameter<std::string>("bTagAlgo")),
00011   minBTagValueBJet_           (cfg.getParameter<double>("minBTagValueBJet")),
00012   maxBTagValueNonBJet_        (cfg.getParameter<double>("maxBTagValueNonBJet")),
00013   useBTagging_                (cfg.getParameter<bool>("useBTagging")),
00014   bTags_                      (cfg.getParameter<unsigned int>("bTags")),
00015   jetCorrectionLevel_         (cfg.getParameter<std::string>("jetCorrectionLevel")),
00016   maxNJets_                   (cfg.getParameter<int>("maxNJets")),
00017   maxNComb_                   (cfg.getParameter<int>("maxNComb")),
00018   maxNrIter_                  (cfg.getParameter<unsigned int>("maxNrIter")),
00019   maxDeltaS_                  (cfg.getParameter<double>("maxDeltaS")),
00020   maxF_                       (cfg.getParameter<double>("maxF")),
00021   jetParam_                   (cfg.getParameter<unsigned>("jetParametrisation")),
00022   constraints_                (cfg.getParameter<std::vector<unsigned> >("constraints")),
00023   mW_                         (cfg.getParameter<double>("mW"  )),
00024   mTop_                       (cfg.getParameter<double>("mTop")),
00025   jetEnergyResolutionScaleFactors_(cfg.getParameter<std::vector<double> >("jetEnergyResolutionScaleFactors")),
00026   jetEnergyResolutionEtaBinning_  (cfg.getParameter<std::vector<double> >("jetEnergyResolutionEtaBinning"))
00027 {
00028   if(cfg.exists("udscResolutions") && cfg.exists("bResolutions")){
00029     udscResolutions_ = cfg.getParameter <std::vector<edm::ParameterSet> >("udscResolutions");
00030     bResolutions_    = cfg.getParameter <std::vector<edm::ParameterSet> >("bResolutions");
00031   }
00032   else if(cfg.exists("udscResolutions") || cfg.exists("bResolutions")){
00033     if(cfg.exists("udscResolutions")) throw cms::Exception("Configuration") << "Parameter 'bResolutions' is needed if parameter 'udscResolutions' is defined!\n";
00034     else                              throw cms::Exception("Configuration") << "Parameter 'udscResolutions' is needed if parameter 'bResolutions' is defined!\n";
00035   }
00036 
00037   // define kinematic fit interface
00038   kinFitter = new TtFullHadKinFitter::KinFit(useBTagging_, bTags_, bTagAlgo_, minBTagValueBJet_, maxBTagValueNonBJet_,
00039                                              udscResolutions_, bResolutions_, jetEnergyResolutionScaleFactors_, 
00040                                              jetEnergyResolutionEtaBinning_, jetCorrectionLevel_, maxNJets_, maxNComb_,
00041                                              maxNrIter_, maxDeltaS_, maxF_, jetParam_, constraints_, mW_, mTop_);
00042 
00043   // produces the following collections
00044   produces< std::vector<pat::Particle> >("PartonsB");
00045   produces< std::vector<pat::Particle> >("PartonsBBar");
00046   produces< std::vector<pat::Particle> >("PartonsLightQ");
00047   produces< std::vector<pat::Particle> >("PartonsLightQBar");
00048   produces< std::vector<pat::Particle> >("PartonsLightP");
00049   produces< std::vector<pat::Particle> >("PartonsLightPBar");
00050 
00051   produces< std::vector<std::vector<int> > >();
00052   produces< std::vector<double> >("Chi2");
00053   produces< std::vector<double> >("Prob");
00054   produces< std::vector<int> >("Status");
00055 }
00056 
00058 TtFullHadKinFitProducer::~TtFullHadKinFitProducer()
00059 {
00060   delete kinFitter;
00061 }
00062 
00064 void 
00065 TtFullHadKinFitProducer::produce(edm::Event& event, const edm::EventSetup& setup)
00066 {
00067   // get jet collection
00068   edm::Handle<std::vector<pat::Jet> > jets;
00069   event.getByLabel(jets_, jets);
00070 
00071   // get match in case that useOnlyMatch_ is true
00072   std::vector<int> match;
00073   bool invalidMatch=false;
00074   if(useOnlyMatch_) {
00075     kinFitter->setUseOnlyMatch(true);
00076     // in case that only a ceratin match should be used, get match here
00077     edm::Handle<std::vector<std::vector<int> > > matches;
00078     event.getByLabel(match_, matches);
00079     match = *(matches->begin());
00080     // check if match is valid
00081     if( match.size()!=nPartons ){
00082       invalidMatch=true;
00083     }
00084     else {
00085       for(unsigned int idx=0; idx<match.size(); ++idx) {
00086         if(match[idx]<0 || match[idx]>=(int)jets->size()) {
00087           invalidMatch=true;
00088           break;
00089         }
00090       }
00091     }
00093     kinFitter->setMatch(match);
00094   }
00095 
00097   kinFitter->setMatchInvalidity(invalidMatch);
00098 
00099   std::list<TtFullHadKinFitter::KinFitResult> fitResults = kinFitter->fit(*jets);
00100 
00101   // pointer for output collections
00102   std::auto_ptr< std::vector<pat::Particle> > pPartonsB( new std::vector<pat::Particle> );
00103   std::auto_ptr< std::vector<pat::Particle> > pPartonsBBar( new std::vector<pat::Particle> );
00104   std::auto_ptr< std::vector<pat::Particle> > pPartonsLightQ   ( new std::vector<pat::Particle> );
00105   std::auto_ptr< std::vector<pat::Particle> > pPartonsLightQBar( new std::vector<pat::Particle> );
00106   std::auto_ptr< std::vector<pat::Particle> > pPartonsLightP   ( new std::vector<pat::Particle> );
00107   std::auto_ptr< std::vector<pat::Particle> > pPartonsLightPBar( new std::vector<pat::Particle> );
00108   // pointer for meta information
00109   std::auto_ptr< std::vector<std::vector<int> > > pCombi ( new std::vector<std::vector<int> > );
00110   std::auto_ptr< std::vector<double> > pChi2  ( new std::vector<double> );
00111   std::auto_ptr< std::vector<double> > pProb  ( new std::vector<double> );
00112   std::auto_ptr< std::vector<int> > pStatus( new std::vector<int> );
00113 
00114   unsigned int iComb = 0;
00115   for(std::list<TtFullHadKinFitter::KinFitResult>::const_iterator res = fitResults.begin(); res != fitResults.end(); ++res){
00116     if(maxNComb_>=1 && iComb==(unsigned int)maxNComb_){ 
00117       break;
00118     }
00119     ++iComb;
00120 
00121     pPartonsB        ->push_back( res->B         );
00122     pPartonsBBar     ->push_back( res->BBar      );
00123     pPartonsLightQ   ->push_back( res->LightQ    );
00124     pPartonsLightQBar->push_back( res->LightQBar );
00125     pPartonsLightP   ->push_back( res->LightP    );
00126     pPartonsLightPBar->push_back( res->LightPBar );
00127 
00128     pCombi ->push_back( res->JetCombi );
00129     pChi2  ->push_back( res->Chi2     );
00130     pProb  ->push_back( res->Prob     );
00131     pStatus->push_back( res->Status   );
00132 
00133   }
00134 
00135   event.put(pCombi);
00136   event.put(pPartonsB        , "PartonsB"        );
00137   event.put(pPartonsBBar     , "PartonsBBar"     );
00138   event.put(pPartonsLightQ   , "PartonsLightQ"   );
00139   event.put(pPartonsLightQBar, "PartonsLightQBar");
00140   event.put(pPartonsLightP   , "PartonsLightP"   );
00141   event.put(pPartonsLightPBar, "PartonsLightPBar");
00142   event.put(pChi2   , "Chi2"   );
00143   event.put(pProb   , "Prob"   );
00144   event.put(pStatus , "Status" );
00145 }
00146 
00147 #include "FWCore/Framework/interface/MakerMacros.h"
00148 DEFINE_FWK_MODULE(TtFullHadKinFitProducer);