CMS 3D CMS Logo

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