CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepJetCombMVAComputer.cc

Go to the documentation of this file.
00001 #include "PhysicsTools/JetMCUtils/interface/combination.h"
00002 
00003 #include "DataFormats/PatCandidates/interface/Jet.h"
00004 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00005 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h"
00006 #include "TopQuarkAnalysis/TopJetCombination/interface/TtSemiLepJetCombEval.h"
00007 #include "TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepJetCombMVAComputer.h"
00008 
00009 TtSemiLepJetCombMVAComputer::TtSemiLepJetCombMVAComputer(const edm::ParameterSet& cfg):
00010   leps_    (cfg.getParameter<edm::InputTag>("leps")),
00011   jets_    (cfg.getParameter<edm::InputTag>("jets")),
00012   mets_    (cfg.getParameter<edm::InputTag>("mets")),
00013   maxNJets_(cfg.getParameter<int>("maxNJets")),
00014   maxNComb_(cfg.getParameter<int>("maxNComb"))
00015 {
00016   produces< std::vector<std::vector<int> > >();
00017   produces< std::vector<double>            >("Discriminators");
00018   produces< std::string                    >("Method");
00019 }
00020 
00021 TtSemiLepJetCombMVAComputer::~TtSemiLepJetCombMVAComputer()
00022 {
00023 }
00024 
00025 void
00026 TtSemiLepJetCombMVAComputer::produce(edm::Event& evt, const edm::EventSetup& setup)
00027 {
00028   std::auto_ptr< std::vector<std::vector<int> > >pOut    (new std::vector<std::vector<int> >);
00029   std::auto_ptr< std::vector<double>            >pOutDisc(new std::vector<double>);
00030   std::auto_ptr< std::string                    >pOutMeth(new std::string);
00031 
00032   mvaComputer.update<TtSemiLepJetCombMVARcd>(setup, "ttSemiLepJetCombMVA");
00033 
00034   // read name of the processor that provides the MVA discriminator
00035   // (to be used as meta information)
00036   edm::ESHandle<PhysicsTools::Calibration::MVAComputerContainer> calibContainer;
00037   setup.get<TtSemiLepJetCombMVARcd>().get( calibContainer );
00038   std::vector<PhysicsTools::Calibration::VarProcessor*> processors
00039     = (calibContainer->find("ttSemiLepJetCombMVA")).getProcessors();
00040   *pOutMeth = ( processors[ processors.size()-3 ] )->getInstanceName();
00041   evt.put(pOutMeth, "Method");
00042 
00043   // get lepton, jets and mets
00044   edm::Handle< edm::View<reco::RecoCandidate> > leptons; 
00045   evt.getByLabel(leps_, leptons);
00046 
00047   edm::Handle< std::vector<pat::Jet> > jets;
00048   evt.getByLabel(jets_, jets);
00049 
00050   edm::Handle< std::vector<pat::MET> > mets;
00051   evt.getByLabel(mets_, mets);
00052 
00053   unsigned int nPartons = 4;
00054 
00055   // skip events with no appropriate lepton candidate,
00056   // empty METs vector or less jets than partons
00057   if( leptons->empty() || mets->empty() || jets->size() < nPartons ) {
00058     std::vector<int> invalidCombi;
00059     for(unsigned int i = 0; i < nPartons; ++i) 
00060       invalidCombi.push_back( -1 );
00061     pOut->push_back( invalidCombi );
00062     evt.put(pOut);
00063     pOutDisc->push_back( 0. );
00064     evt.put(pOutDisc, "Discriminators");
00065     return;
00066   }
00067 
00068   const math::XYZTLorentzVector lepton = leptons->begin()->p4();
00069 
00070   const pat::MET *met = &(*mets)[0];
00071 
00072   // analyze jet combinations
00073   std::vector<int> jetIndices;
00074   for(unsigned int i=0; i<jets->size(); ++i){
00075     if(maxNJets_ >= (int) nPartons && maxNJets_ == (int) i) break;
00076     jetIndices.push_back(i);
00077   }
00078   
00079   std::vector<int> combi;
00080   for(unsigned int i=0; i<nPartons; ++i) 
00081     combi.push_back(i);
00082 
00083   typedef std::pair<double, std::vector<int> > discCombPair;
00084   std::list<discCombPair> discCombList;
00085 
00086   do{
00087     for(int cnt = 0; cnt < TMath::Factorial( combi.size() ); ++cnt){
00088       // take into account indistinguishability of the two jets from the hadr. W decay,
00089       // reduces combinatorics by a factor of 2
00090       if(combi[TtSemiLepEvtPartons::LightQ] < combi[TtSemiLepEvtPartons::LightQBar]) {
00091 
00092         TtSemiLepJetComb jetComb(*jets, combi, lepton, *met);
00093 
00094         // feed MVA input variables into a ValueList
00095         PhysicsTools::Variable::ValueList values;
00096         evaluateTtSemiLepJetComb(values, jetComb);
00097 
00098         // get discriminator from the MVAComputer
00099         double discrim = mvaComputer->eval( values );
00100 
00101         discCombList.push_back( std::make_pair(discrim, combi) );
00102 
00103       }
00104       next_permutation( combi.begin() , combi.end() );
00105     }
00106   }
00107   while(stdcomb::next_combination( jetIndices.begin(), jetIndices.end(), combi.begin(), combi.end() ));
00108 
00109   // sort results w.r.t. discriminator values
00110   discCombList.sort();
00111 
00112   // write result into the event
00113   // (starting with the JetComb having the highest discriminator value -> reverse iterator)
00114   unsigned int iDiscComb = 0;
00115   typedef std::list<discCombPair>::reverse_iterator discCombIterator;
00116   for(discCombIterator discCombPair = discCombList.rbegin(); discCombPair != discCombList.rend(); ++discCombPair) {
00117     if(maxNComb_ >= 1 && iDiscComb == (unsigned int) maxNComb_) break;
00118     pOut    ->push_back( discCombPair->second );
00119     pOutDisc->push_back( discCombPair->first  );
00120     iDiscComb++;
00121   }
00122   evt.put(pOut);
00123   evt.put(pOutDisc, "Discriminators");
00124 }
00125 
00126 void 
00127 TtSemiLepJetCombMVAComputer::beginJob()
00128 {
00129 }
00130 
00131 void 
00132 TtSemiLepJetCombMVAComputer::endJob()
00133 {
00134 }
00135 
00136 // implement the plugins for the computer container
00137 // -> register TtSemiLepJetCombMVARcd
00138 // -> define TtSemiLepJetCombMVAFileSource
00139 MVA_COMPUTER_CONTAINER_IMPLEMENT(TtSemiLepJetCombMVA);