CMS 3D CMS Logo

TtFullHadSignalSelMVATrainer.cc
Go to the documentation of this file.
1 #include "TMath.h"
2 #include <algorithm>
3 
6 
9 
11 
12 
14  jetsToken_ (consumes< std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
15  genEvtToken_ (consumes<TtGenEvent>(edm::InputTag("genEvt"))),
16  whatData_ (cfg.getParameter<int>("whatData")),
17  maxEv_ (cfg.getParameter<int>("maxEv")),
18  weight_ (cfg.getParameter<double>("weight"))
19 {
20 }
21 
23 {
24 }
25 
26 void
28 {
29  //communication with CMSSW CondDB
30  mvaComputer.update<TtFullHadSignalSelMVARcd>("trainer", setup, "traintreeSaver");
31 
32  // can occur in the last iteration when the
33  // MVATrainer is about to save the result
34  if(!mvaComputer) return;
35 
36  // get the jets out of the event
38  evt.getByToken(jetsToken_, jets);
39 
40  //count the number of selected events
41  selEv++;
42  //skip event if enough events are already selected
43  if(selEv>maxEv_ && maxEv_!=-1) return;
44 
45  //calculation of InputVariables
46  //see TopQuarkAnalysis/TopTools/interface/TtFullHadSignalSel.h
47  // /src/TtFullHadSignalSel.cc
48  //all objects, jets, which are needed for the calculation
49  //of the input-variables have to be passed to this class
50 
52 
53  if(whatData_==-1) {
54  //your training-file contains both, signal and background events
55 
57  evt.getByToken(genEvtToken_, genEvt);
58 
59  bool isSignal = false;
60  if(genEvt->isTtBar()){
61  if(genEvt->isFullHadronic()) isSignal = true;
62  }
63  evaluateTtFullHadSignalSel(mvaComputer, selection, weight_, true, isSignal);
64  }
65  else {
66 
67  if(whatData_==1){
68  //your tree contains only signal events
69 
70  evaluateTtFullHadSignalSel(mvaComputer, selection, weight_, true, true);
71  }
72  else if(whatData_==0){
73  //your tree contains only background events
74 
75  evaluateTtFullHadSignalSel(mvaComputer, selection, weight_, true, false);
76  }
77  else std::cout<<"Config File Error!! Please check <whatData> in TtFullHadSignalSelMVATrainer_cfi";
78  }
79 }
80 
82  selEv = 0;
83  if(whatData_!=-1 && whatData_!=0 && whatData_!=1){
84  std::cout<<"Config File Error!! Please check <whatData> in TtFullHadSignalSelMVATrainer_cfi"<<std::endl;;
85  return;
86  }
87 }
88 
89 // implement the plugins for the trainer
90 // -> defines TtFullHadSignalSelMVAContainerSaveCondDB
91 // -> defines TtFullHadSignalSelMVASaveFile
92 // -> defines TtFullHadSignalSelMVATrainerLooper
93 MVA_TRAINER_IMPLEMENT(TtFullHadSignalSelMVA);
void analyze(const edm::Event &evt, const edm::EventSetup &setup) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
#define MVA_TRAINER_IMPLEMENT(N)
Definition: HelperMacros.h:40
selection
main part
Definition: corrVsCorr.py:98
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: HeavyIon.h:7
bool isFullHadronic(bool excludeTauLeptons=false) const
check if the event can be classified as full hadronic
Definition: TtGenEvent.h:36
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
Definition: Jet.py:1
vector< PseudoJet > jets
double evaluateTtFullHadSignalSel(PhysicsTools::MVAComputerCache &mvaComputer, const TtFullHadSignalSel &sigsel, double weight=1.0, const bool training=false, const bool isSignal=false)
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
edm::EDGetTokenT< TtGenEvent > genEvtToken_
bool update(const Calibration::MVAComputer *computer)
HLT enums.
TtFullHadSignalSelMVATrainer(const edm::ParameterSet &)
PhysicsTools::MVAComputerCache mvaComputer
bool isTtBar() const
check if the event can be classified as ttbar
Definition: TtGenEvent.h:30