CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtSemiLepJetCombMVAComputer.cc
Go to the documentation of this file.
2 
4 
6  lepsToken_ (consumes< edm::View<reco::RecoCandidate>>(cfg.getParameter<edm::InputTag>("leps"))),
7  jetsToken_ (consumes< std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
8  metsToken_ (consumes< std::vector<pat::MET> >(cfg.getParameter<edm::InputTag>("mets"))),
9  maxNJets_(cfg.getParameter<int>("maxNJets")),
10  maxNComb_(cfg.getParameter<int>("maxNComb"))
11 {
12  produces<std::vector<std::vector<int> > >();
13  produces<std::vector<double> >("Discriminators");
14  produces<std::string >("Method");
15  produces<int >("NumberOfConsideredJets");
16 }
17 
19 {
20 }
21 
22 void
24 {
25  std::auto_ptr<std::vector<std::vector<int> > > pOut (new std::vector<std::vector<int> >);
26  std::auto_ptr<std::vector<double> > pOutDisc(new std::vector<double>);
27  std::auto_ptr<std::string > pOutMeth(new std::string);
28  std::auto_ptr<int > pJetsConsidered(new int);
29 
30  mvaComputer.update<TtSemiLepJetCombMVARcd>(setup, "ttSemiLepJetCombMVA");
31 
32  // read name of the processor that provides the MVA discriminator
33  // (to be used as meta information)
35  setup.get<TtSemiLepJetCombMVARcd>().get( calibContainer );
36  std::vector<PhysicsTools::Calibration::VarProcessor*> processors
37  = (calibContainer->find("ttSemiLepJetCombMVA")).getProcessors();
38  *pOutMeth = ( processors[ processors.size()-3 ] )->getInstanceName();
39  evt.put(pOutMeth, "Method");
40 
41  // get lepton, jets and mets
43  evt.getByToken(lepsToken_, leptons);
44 
46  evt.getByToken(jetsToken_, jets);
47 
49  evt.getByToken(metsToken_, mets);
50 
51  const unsigned int nPartons = 4;
52 
53  // skip events with no appropriate lepton candidate,
54  // empty METs vector or less jets than partons
55  if( leptons->empty() || mets->empty() || jets->size() < nPartons ) {
56  std::vector<int> invalidCombi;
57  for(unsigned int i = 0; i < nPartons; ++i)
58  invalidCombi.push_back( -1 );
59  pOut->push_back( invalidCombi );
60  evt.put(pOut);
61  pOutDisc->push_back( 0. );
62  evt.put(pOutDisc, "Discriminators");
63  *pJetsConsidered = jets->size();
64  evt.put(pJetsConsidered, "NumberOfConsideredJets");
65  return;
66  }
67 
68  const math::XYZTLorentzVector lepton = leptons->begin()->p4();
69 
70  const pat::MET *met = &(*mets)[0];
71 
72  // analyze jet combinations
73  std::vector<int> jetIndices;
74  for(unsigned int i=0; i<jets->size(); ++i){
75  if(maxNJets_ >= (int) nPartons && maxNJets_ == (int) i) {
76  *pJetsConsidered = i;
77  break;
78  }
79  jetIndices.push_back(i);
80  }
81 
82  std::vector<int> combi;
83  for(unsigned int i=0; i<nPartons; ++i)
84  combi.push_back(i);
85 
86  typedef std::pair<double, std::vector<int> > discCombPair;
87  std::list<discCombPair> discCombList;
88 
89  do{
90  for(int cnt = 0; cnt < TMath::Factorial( combi.size() ); ++cnt){
91  // take into account indistinguishability of the two jets from the hadr. W decay,
92  // reduces combinatorics by a factor of 2
94 
95  TtSemiLepJetComb jetComb(*jets, combi, lepton, *met);
96 
97  // feed MVA input variables into a ValueList
99  evaluateTtSemiLepJetComb(values, jetComb);
100 
101  // get discriminator from the MVAComputer
102  double discrim = mvaComputer->eval( values );
103 
104  discCombList.push_back( std::make_pair(discrim, combi) );
105 
106  }
107  next_permutation( combi.begin() , combi.end() );
108  }
109  }
110  while(stdcomb::next_combination( jetIndices.begin(), jetIndices.end(), combi.begin(), combi.end() ));
111 
112  // sort results w.r.t. discriminator values
113  discCombList.sort();
114 
115  // write result into the event
116  // (starting with the JetComb having the highest discriminator value -> reverse iterator)
117  unsigned int iDiscComb = 0;
118  typedef std::list<discCombPair>::reverse_iterator discCombIterator;
119  for(discCombIterator discCombPair = discCombList.rbegin(); discCombPair != discCombList.rend(); ++discCombPair) {
120  if(maxNComb_ >= 1 && iDiscComb == (unsigned int) maxNComb_) break;
121  pOut ->push_back( discCombPair->second );
122  pOutDisc->push_back( discCombPair->first );
123  iDiscComb++;
124  }
125  evt.put(pOut);
126  evt.put(pOutDisc, "Discriminators");
127  evt.put(pJetsConsidered, "NumberOfConsideredJets");
128 }
129 
130 void
132 {
133 }
134 
135 void
137 {
138 }
139 
140 // implement the plugins for the computer container
141 // -> register TtSemiLepJetCombMVARcd
142 // -> define TtSemiLepJetCombMVAFileSource
143 MVA_COMPUTER_CONTAINER_IMPLEMENT(TtSemiLepJetCombMVA);
Analysis-level MET class.
Definition: MET.h:43
int i
Definition: DBlmapReader.cc:9
static const unsigned int nPartons
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
tuple cfg
Definition: looper.py:293
double eval(Iterator_t first, Iterator_t last) const
evaluate variables given by a range of iterators given by first and last
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
vector< PseudoJet > jets
void get(HolderT &iHolder) const
PhysicsTools::MVAComputerCache mvaComputer
Helper class that can contain an list of identifier-value pairs.
Definition: Variable.h:81
bool update(const Calibration::MVAComputer *computer)
virtual void produce(edm::Event &evt, const edm::EventSetup &setup)
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
void evaluateTtSemiLepJetComb(PhysicsTools::Variable::ValueList &values, const TtSemiLepJetComb &jetComb)
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
#define MVA_COMPUTER_CONTAINER_IMPLEMENT(N)
Definition: HelperMacros.h:46
TtSemiLepJetCombMVAComputer(const edm::ParameterSet &)
Common calculator class to keep multivariate analysis variables for jet combinations in semi-leptonic...
bool next_combination(BidIt n_begin, BidIt n_end, BidIt r_begin, BidIt r_end)
Definition: combination.h:22
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")