CMS 3D CMS Logo

TtDilepEvtSolution.cc
Go to the documentation of this file.
1 //
2 //
3 
7 
9 {
10  jetCorrScheme_ = 0;
11  wpDecay_ = "NotDefined";
12  wmDecay_ = "NotDefined";
13  bestSol_ = false;
14  topmass_ = 0.;
15  weightmax_ = 0.;
16 }
17 
19 {
20 }
21 
22 //-------------------------------------------
23 // get calibrated base objects
24 //-------------------------------------------
26 {
27  // WARNING this is obsolete and only
28  // kept for backwards compatibility
29  if(jetCorrScheme_==1){
30  //jet calibrated according to MC truth
31  return jetB_->correctedJet("HAD", "B");
32  }
33  else if(jetCorrScheme_==2){
34  return jetB_->correctedJet("HAD", "B");
35  }
36  else{
37  return *jetB_;
38  }
39 }
40 
42 {
43  // WARNING this is obsolete and only
44  // kept for backwards compatibility
45  if(jetCorrScheme_==1){
46  //jet calibrated according to MC truth
47  return jetBbar_->correctedJet("HAD", "B");
48  }
49  else if(jetCorrScheme_==2){
50  return jetBbar_->correctedJet("HAD", "B");
51  }
52  else{
53  return *jetBbar_;
54  }
55 }
56 
57 //-------------------------------------------
58 // returns the 4-vector of the positive
59 // lepton, with the charge and the pdgId
60 //-------------------------------------------
62 {
64  if(wpDecay_ == "muon"){
65  p = reco::Particle(+1, getMuonp().p4() );
66  p.setPdgId(-11);
67  }
68  if(wpDecay_ == "electron"){
69  p = reco::Particle(+1, getElectronp().p4() );
70  p.setPdgId(-13);
71  }
72  if(wmDecay_ == "tau"){
73  p = reco::Particle(+1, getTaup().p4() );
74  p.setPdgId(-15);
75  }
76  return p;
77 }
78 
79 //-------------------------------------------
80 // miscellaneous
81 //-------------------------------------------
83 {
84  double distance = 0.;
85  if(!getGenB() || !getGenBbar()) return distance;
86  distance += reco::deltaR(getCalJetB(),*getGenB());
87  distance += reco::deltaR(getCalJetBbar(),*getGenBbar());
88  return distance;
89 }
90 
92 {
93  double distance = 0.;
94  if(!getGenLepp() || !getGenLepm()) return distance;
95  if(getWpDecay()=="electron")
96  distance += reco::deltaR(getElectronp(),*getGenLepp());
97  else if(getWpDecay()=="muon")
98  distance += reco::deltaR(getMuonp(),*getGenLepp());
99  else if(getWpDecay()=="tau")
100  distance += reco::deltaR(getTaup(),*getGenLepp());
101  if(getWmDecay()=="electron")
102  distance += reco::deltaR(getElectronm(),*getGenLepm());
103  else if(getWmDecay()=="muon")
104  distance += reco::deltaR(getMuonm(),*getGenLepm());
105  else if(getWmDecay()=="tau")
106  distance += reco::deltaR(getTaum(),*getGenLepm());
107  return distance;
108 }
109 
110 //-------------------------------------------
111 // returns the 4-vector of the negative
112 // lepton, with the charge and the pdgId
113 //-------------------------------------------
115 {
117  if(wmDecay_ == "electron"){
118  p = reco::Particle(-1, getElectronm().p4() );
119  p.setPdgId(11);
120  }
121  if(wmDecay_ == "muon"){
122  p = reco::Particle(-1, getMuonm().p4() );
123  p.setPdgId(13);
124  }
125  if(wmDecay_ == "tau"){
126  p = reco::Particle(-1, getTaum().p4() );
127  p.setPdgId(15);
128  }
129  return p;
130 }
131 
132 //-------------------------------------------
133 // get info on the outcome of the signal
134 //selection LR
135 //-------------------------------------------
136 double TtDilepEvtSolution::getLRSignalEvtObsVal(unsigned int selObs) const
137 {
138  double val = -999.;
139  for(size_t i=0; i<lrSignalEvtVarVal_.size(); i++){
140  if(lrSignalEvtVarVal_[i].first == selObs) val = lrSignalEvtVarVal_[i].second;
141  }
142  return val;
143 }
144 
145 //-------------------------------------------
146 // set the generated event
147 //-------------------------------------------
149  if( !aGenEvt->isFullLeptonic() ){
150  edm::LogInfo( "TtGenEventNotFilled" ) << "genEvt is not di-leptonic; TtGenEvent is not filled";
151  return;
152  }
154 }
155 
156 //-------------------------------------------
157 // set the outcome of the signal selection LR
158 //-------------------------------------------
159 void TtDilepEvtSolution::setLRSignalEvtObservables(const std::vector<std::pair<unsigned int, double> >& varval)
160 {
161  lrSignalEvtVarVal_.clear();
162  for(size_t ise = 0; ise<varval.size(); ise++) lrSignalEvtVarVal_.push_back(varval[ise]);
163 }
pat::Jet getCalJetBbar() const
const reco::GenParticle * getGenLepm() const
edm::Ref< std::vector< pat::Jet > > jetB_
const reco::GenParticle * getGenBbar() const
pat::Electron getElectronm() const
std::string getWmDecay() const
std::string getWpDecay() const
double getLRSignalEvtObsVal(unsigned int) const
bool isFullLeptonic(bool excludeTauLeptons=false) const
check if the event can be classified as full leptonic
Definition: TtGenEvent.h:40
void setLRSignalEvtObservables(const std::vector< std::pair< unsigned int, double > > &)
double p4[4]
Definition: TauolaWrapper.h:92
reco::Particle getLeptNeg() const
pat::Jet getJetBbar() const
pat::Electron getElectronp() const
edm::Ref< std::vector< pat::Jet > > jetBbar_
reco::Particle getLeptPos() const
pat::Tau getTaup() const
edm::RefProd< TtGenEvent > theGenEvt_
pat::Muon getMuonm() const
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
pat::Jet getJetB() const
double getLeptonResidual() const
void setGenEvt(const edm::Handle< TtGenEvent > &)
pat::Tau getTaum() const
pat::Jet getCalJetB() const
std::vector< std::pair< unsigned int, double > > lrSignalEvtVarVal_
void setPdgId(int pdgId)
Definition: Particle.h:136
const reco::GenParticle * getGenLepp() const
Analysis-level calorimeter jet class.
Definition: Jet.h:80
double getJetResidual() const
pat::Muon getMuonp() const
const reco::GenParticle * getGenB() const