CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
TtSemiLepJetCombMVATrainer Class Reference

#include <TtSemiLepJetCombMVATrainer.h>

Inheritance diagram for TtSemiLepJetCombMVATrainer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 TtSemiLepJetCombMVATrainer (const edm::ParameterSet &)
 
 ~TtSemiLepJetCombMVATrainer ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

virtual void analyze (const edm::Event &evt, const edm::EventSetup &setup)
 
virtual void beginJob ()
 
virtual void endJob ()
 
WDecay::LeptonType readLeptonType (const std::string &str)
 

Private Attributes

edm::InputTag jets_
 
edm::InputTag leps_
 
WDecay::LeptonType leptonType_
 
edm::InputTag matching_
 
int maxNJets_
 
edm::InputTag mets_
 
PhysicsTools::MVAComputerCache mvaComputer
 
unsigned int nEvents [5]
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 23 of file TtSemiLepJetCombMVATrainer.h.

Constructor & Destructor Documentation

TtSemiLepJetCombMVATrainer::TtSemiLepJetCombMVATrainer ( const edm::ParameterSet cfg)
explicit

Definition at line 17 of file TtSemiLepJetCombMVATrainer.cc.

17  :
18  leps_ (cfg.getParameter<edm::InputTag>("leps" )),
19  jets_ (cfg.getParameter<edm::InputTag>("jets" )),
20  mets_ (cfg.getParameter<edm::InputTag>("mets" )),
21  matching_ (cfg.getParameter<edm::InputTag>("matching")),
22  maxNJets_ (cfg.getParameter<int> ("maxNJets")),
24 {
25 }
T getParameter(std::string const &) const
WDecay::LeptonType readLeptonType(const std::string &str)
TtSemiLepJetCombMVATrainer::~TtSemiLepJetCombMVATrainer ( )

Definition at line 27 of file TtSemiLepJetCombMVATrainer.cc.

28 {
29 }

Member Function Documentation

void TtSemiLepJetCombMVATrainer::analyze ( const edm::Event evt,
const edm::EventSetup setup 
)
privatevirtual

Implements edm::EDAnalyzer.

Definition at line 39 of file TtSemiLepJetCombMVATrainer.cc.

References PhysicsTools::Variable::ValueList::add(), PhysicsTools::MVAComputer::eval(), evaluateTtSemiLepJetComb(), TtGenEvtProducer_cfi::genEvt, edm::Event::getByLabel(), i, fwrapper::jets, jets_, PhysicsTools::MVATrainer::kTargetId, leps_, EgammaValidation_Wenu_cff::leptons, leptonType_, TtSemiLepEvtPartons::LightQ, TtSemiLepEvtPartons::LightQBar, matching_, maxNJets_, CaloMET_cfi::met, mets_, mvaComputer, nEvents, stdcomb::next_combination(), nPartons, HcalObjRepresent::setup(), PhysicsTools::MVAComputerCache::update(), and makeHLTPrescaleTable::values.

40 {
41  mvaComputer.update<TtSemiLepJetCombMVARcd>("trainer", setup, "ttSemiLepJetCombMVA");
42 
43  // can occur in the last iteration when the
44  // MVATrainer is about to save the result
45  if(!mvaComputer) return;
46 
47  nEvents[0]++;
48 
50  evt.getByLabel("genEvt", genEvt);
51 
53  evt.getByLabel(leps_, leptons);
54 
55  // skip events with no appropriate lepton candidate
56  if( leptons->empty() ) return;
57 
58  nEvents[1]++;
59 
60  const math::XYZTLorentzVector lepton = leptons->begin()->p4();
61 
63  evt.getByLabel(mets_, mets);
64 
65  // skip events with empty METs vector
66  if( mets->empty() ) return;
67 
68  nEvents[2]++;
69 
70  const pat::MET *met = &(*mets)[0];
71 
73  evt.getByLabel(jets_, jets);
74 
75  // skip events with less jets than partons
76  unsigned int nPartons = 4;
77  if( jets->size() < nPartons ) return;
78 
79  nEvents[3]++;
80 
82  std::vector<int> matching;
83  // get jet-parton matching if signal channel
84  if(genEvt->semiLeptonicChannel() == leptonType_) {
85  evt.getByLabel(matching_, matchingHandle);
86  if( matchingHandle->empty() )
87  throw cms::Exception("EmptyProduct")
88  << "Empty vector from jet-parton matching. This should not happen! \n";
89  matching = *(matchingHandle->begin());
90  if(matching.size() < nPartons) return;
91  // skip events that were affected by the outlier
92  // rejection in the jet-parton matching
93  for(unsigned int i = 0; i < matching.size(); ++i) {
94  if(matching[i] == -3) continue; // -3: parton was chosen to be excluded from jet-parton matching
95  if(matching[i] < 0 || matching[i] >= (int)jets->size())
96  return;
97  }
98  }
99  // use dummy for matching if not signal channel
100  else
101  for(unsigned int i = 0; i < nPartons; i++)
102  matching.push_back( -1 );
103 
104  nEvents[4]++;
105 
106  // take into account indistinguishability of the two jets from the hadr. W decay,
107  // reduces combinatorics by a factor of 2
109  int iTemp = matching[TtSemiLepEvtPartons::LightQ];
111  matching[TtSemiLepEvtPartons::LightQBar] = iTemp;
112  }
113 
114  std::vector<int> jetIndices;
115  for(unsigned int i=0; i<jets->size(); ++i) {
116  if(maxNJets_ >= (int) nPartons && i == (unsigned int) maxNJets_) break;
117  jetIndices.push_back(i);
118  }
119 
120  std::vector<int> combi;
121  for(unsigned int i = 0; i < nPartons; ++i)
122  combi.push_back(i);
123 
124  // use a ValueList to collect all variables for the MVAComputer
126 
127  // set target variable to false like for background
128  // this is necessary when passing all combinations at once
130 
131  do {
132  // number of possible combinations from number of partons: e.g. 4! = 24
133  for(unsigned int cnt = 0; cnt < TMath::Factorial( combi.size() ); ++cnt) {
134 
135  // take into account indistinguishability of the two jets from the hadr. W decay,
136  // reduces combinatorics by a factor of 2
137  if(combi[TtSemiLepEvtPartons::LightQ] < combi[TtSemiLepEvtPartons::LightQBar]) {
138 
139  TtSemiLepJetComb jetComb(*jets, combi, lepton, *met);
140 
141  bool trueCombi = false;
142  // true combination only if signal channel
143  // and in agreement with matching
144  if(genEvt->semiLeptonicChannel()==leptonType_) {
145  trueCombi = true;
146  for(unsigned int i = 0; i < matching.size(); ++i) {
147  if(combi[i]!=matching[i] && matching[i]!=-3) {
148  trueCombi = false;
149  break;
150  }
151  }
152  }
153 
154  // feed MVA input variables for this jetComb into the ValueList
155  values.add("target", trueCombi);
156  evaluateTtSemiLepJetComb(values, jetComb);
157 
158  }
159 
160  next_permutation( combi.begin() , combi.end() );
161  }
162  }
163  while(stdcomb::next_combination( jetIndices.begin(), jetIndices.end(), combi.begin(), combi.end() ));
164 
165  // pass MVA input variables for all jet combinations in this event
166  // to the MVAComputer for training
167  mvaComputer->eval( values );
168 }
PhysicsTools::MVAComputerCache mvaComputer
Analysis-level MET class.
Definition: MET.h:42
int i
Definition: DBlmapReader.cc:9
static const unsigned int nPartons
double eval(Iterator_t first, Iterator_t last) const
evaluate variables given by a range of iterators given by first and last
static const AtomicId kTargetId
Definition: MVATrainer.h:59
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
vector< PseudoJet > jets
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
Helper class that can contain an list of identifier-value pairs.
Definition: Variable.h:82
bool update(const Calibration::MVAComputer *computer)
Helper class that can contain an identifier-value pair.
Definition: Variable.h:52
void evaluateTtSemiLepJetComb(PhysicsTools::Variable::ValueList &values, const TtSemiLepJetComb &jetComb)
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="")
void add(AtomicId id, double value)
Definition: Variable.h:107
void TtSemiLepJetCombMVATrainer::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 32 of file TtSemiLepJetCombMVATrainer.cc.

References i, and nEvents.

33 {
34  for(unsigned int i = 0; i < 5; i++)
35  nEvents[i] = 0;
36 }
int i
Definition: DBlmapReader.cc:9
void TtSemiLepJetCombMVATrainer::endJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 171 of file TtSemiLepJetCombMVATrainer.cc.

References create_public_lumi_plots::log, and nEvents.

172 {
173  edm::LogInfo log("TtSemiLepJetCombMVATrainer");
174  log << "Number of events... \n"
175  << "...passed to the trainer : " << std::setw(7) << nEvents[0] << "\n"
176  << "...rejected since no lepton candidate : " << std::setw(7) << nEvents[0]-nEvents[1] << "\n"
177  << "...rejected since no MET object : " << std::setw(7) << nEvents[1]-nEvents[2] << "\n"
178  << "...rejected since not enough jets : " << std::setw(7) << nEvents[2]-nEvents[3] << "\n"
179  << "...rejected due to bad jet-parton matching : " << std::setw(7) << nEvents[3]-nEvents[4] << "\n"
180  << "...accepted for training : " << std::setw(7) << nEvents[4] << "\n";
181 }
WDecay::LeptonType TtSemiLepJetCombMVATrainer::readLeptonType ( const std::string &  str)
private

Definition at line 184 of file TtSemiLepJetCombMVATrainer.cc.

References edm::hlt::Exception, WDecay::kElec, WDecay::kMuon, and WDecay::kTau.

185 {
186  if (str == "kElec") return WDecay::kElec;
187  else if(str == "kMuon") return WDecay::kMuon;
188  else if(str == "kTau" ) return WDecay::kTau;
189  else throw cms::Exception("Configuration")
190  << "Chosen leptonType is not supported: " << str << "\n";
191 }

Member Data Documentation

edm::InputTag TtSemiLepJetCombMVATrainer::jets_
private

Definition at line 39 of file TtSemiLepJetCombMVATrainer.h.

Referenced by analyze().

edm::InputTag TtSemiLepJetCombMVATrainer::leps_
private

Definition at line 38 of file TtSemiLepJetCombMVATrainer.h.

Referenced by analyze().

WDecay::LeptonType TtSemiLepJetCombMVATrainer::leptonType_
private

Definition at line 45 of file TtSemiLepJetCombMVATrainer.h.

Referenced by analyze().

edm::InputTag TtSemiLepJetCombMVATrainer::matching_
private

Definition at line 41 of file TtSemiLepJetCombMVATrainer.h.

Referenced by analyze().

int TtSemiLepJetCombMVATrainer::maxNJets_
private

Definition at line 43 of file TtSemiLepJetCombMVATrainer.h.

Referenced by analyze().

edm::InputTag TtSemiLepJetCombMVATrainer::mets_
private

Definition at line 40 of file TtSemiLepJetCombMVATrainer.h.

Referenced by analyze().

PhysicsTools::MVAComputerCache TtSemiLepJetCombMVATrainer::mvaComputer
private

Definition at line 47 of file TtSemiLepJetCombMVATrainer.h.

Referenced by analyze().

unsigned int TtSemiLepJetCombMVATrainer::nEvents[5]
private

Definition at line 49 of file TtSemiLepJetCombMVATrainer.h.

Referenced by analyze(), beginJob(), and endJob().