CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/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   resolutionSmearFactor_(cfg.getParameter<double>("resolutionSmearFactor"))
00026 {
00027   if(cfg.exists("udscResolutions") && cfg.exists("bResolutions")){
00028     udscResolutions_ = cfg.getParameter <std::vector<edm::ParameterSet> >("udscResolutions");
00029     bResolutions_    = cfg.getParameter <std::vector<edm::ParameterSet> >("bResolutions");
00030   }
00031   else if(cfg.exists("udscResolutions") || cfg.exists("bResolutions")){
00032     if(cfg.exists("udscResolutions")) throw cms::Exception("WrongConfig") << "Parameter 'bResolutions' is needed if parameter 'udscResolutions' is defined!\n";
00033     else                              throw cms::Exception("WrongConfig") << "Parameter 'udscResolutions' is needed if parameter 'bResolutions' is defined!\n";
00034   }
00035 
00036   // define kinematic fit interface
00037   kinFitter = new TtFullHadKinFitter::KinFit(useBTagging_, bTags_, bTagAlgo_, minBTagValueBJet_, maxBTagValueNonBJet_,
00038                                              udscResolutions_, bResolutions_, resolutionSmearFactor_ ,
00039                                              jetCorrectionLevel_, maxNJets_, maxNComb_,
00040                                              maxNrIter_, maxDeltaS_, maxF_, jetParam_, constraints_, mW_, mTop_);
00041 
00042   // produces the following collections
00043   produces< std::vector<pat::Particle> >("PartonsB");
00044   produces< std::vector<pat::Particle> >("PartonsBBar");
00045   produces< std::vector<pat::Particle> >("PartonsLightQ");
00046   produces< std::vector<pat::Particle> >("PartonsLightQBar");
00047   produces< std::vector<pat::Particle> >("PartonsLightP");
00048   produces< std::vector<pat::Particle> >("PartonsLightPBar");
00049 
00050   produces< std::vector<std::vector<int> > >();
00051   produces< std::vector<double> >("Chi2");
00052   produces< std::vector<double> >("Prob");
00053   produces< std::vector<int> >("Status");
00054 }
00055 
00057 TtFullHadKinFitProducer::~TtFullHadKinFitProducer()
00058 {
00059   delete kinFitter;
00060 }
00061 
00063 void 
00064 TtFullHadKinFitProducer::produce(edm::Event& event, const edm::EventSetup& setup)
00065 {
00066   // get jet collection
00067   edm::Handle<std::vector<pat::Jet> > jets;
00068   event.getByLabel(jets_, jets);
00069 
00070   // get match in case that useOnlyMatch_ is true
00071   std::vector<int> match;
00072   bool invalidMatch=false;
00073   if(useOnlyMatch_) {
00074     kinFitter->setUseOnlyMatch(true);
00075     // in case that only a ceratin match should be used, get match here
00076     edm::Handle<std::vector<std::vector<int> > > matches;
00077     event.getByLabel(match_, matches);
00078     match = *(matches->begin());
00079     // check if match is valid
00080     if( match.size()!=nPartons ){
00081       invalidMatch=true;
00082     }
00083     else {
00084       for(unsigned int idx=0; idx<match.size(); ++idx) {
00085         if(match[idx]<0 || match[idx]>=(int)jets->size()) {
00086           invalidMatch=true;
00087           break;
00088         }
00089       }
00090     }
00092     kinFitter->setMatch(match);
00093   }
00094 
00096   kinFitter->setMatchInvalidity(invalidMatch);
00097 
00098   std::list<TtFullHadKinFitter::KinFitResult> fitResults = kinFitter->fit(*jets);
00099 
00100   // pointer for output collections
00101   std::auto_ptr< std::vector<pat::Particle> > pPartonsB( new std::vector<pat::Particle> );
00102   std::auto_ptr< std::vector<pat::Particle> > pPartonsBBar( new std::vector<pat::Particle> );
00103   std::auto_ptr< std::vector<pat::Particle> > pPartonsLightQ   ( new std::vector<pat::Particle> );
00104   std::auto_ptr< std::vector<pat::Particle> > pPartonsLightQBar( new std::vector<pat::Particle> );
00105   std::auto_ptr< std::vector<pat::Particle> > pPartonsLightP   ( new std::vector<pat::Particle> );
00106   std::auto_ptr< std::vector<pat::Particle> > pPartonsLightPBar( new std::vector<pat::Particle> );
00107   // pointer for meta information
00108   std::auto_ptr< std::vector<std::vector<int> > > pCombi ( new std::vector<std::vector<int> > );
00109   std::auto_ptr< std::vector<double> > pChi2  ( new std::vector<double> );
00110   std::auto_ptr< std::vector<double> > pProb  ( new std::vector<double> );
00111   std::auto_ptr< std::vector<int> > pStatus( new std::vector<int> );
00112 
00113   unsigned int iComb = 0;
00114   for(std::list<TtFullHadKinFitter::KinFitResult>::const_iterator res = fitResults.begin(); res != fitResults.end(); ++res){
00115     if(maxNComb_>=1 && iComb==(unsigned int)maxNComb_){ 
00116       break;
00117     }
00118     ++iComb;
00119 
00120     pPartonsB        ->push_back( res->B         );
00121     pPartonsBBar     ->push_back( res->BBar      );
00122     pPartonsLightQ   ->push_back( res->LightQ    );
00123     pPartonsLightQBar->push_back( res->LightQBar );
00124     pPartonsLightP   ->push_back( res->LightP    );
00125     pPartonsLightPBar->push_back( res->LightPBar );
00126 
00127     pCombi ->push_back( res->JetCombi );
00128     pChi2  ->push_back( res->Chi2     );
00129     pProb  ->push_back( res->Prob     );
00130     pStatus->push_back( res->Status   );
00131 
00132   }
00133 
00134   event.put(pCombi);
00135   event.put(pPartonsB        , "PartonsB"        );
00136   event.put(pPartonsBBar     , "PartonsBBar"     );
00137   event.put(pPartonsLightQ   , "PartonsLightQ"   );
00138   event.put(pPartonsLightQBar, "PartonsLightQBar");
00139   event.put(pPartonsLightP   , "PartonsLightP"   );
00140   event.put(pPartonsLightPBar, "PartonsLightPBar");
00141   event.put(pChi2   , "Chi2"   );
00142   event.put(pProb   , "Prob"   );
00143   event.put(pStatus , "Status" );
00144 }
00145 
00146 #include "FWCore/Framework/interface/MakerMacros.h"
00147 DEFINE_FWK_MODULE(TtFullHadKinFitProducer);