DESY Hbb Analysis Framework
Candidate.cc
Go to the documentation of this file.
1 // system include files
2 //
3 #include <iostream>
4 // user include files
6 
7 
8 //
9 // class declaration
10 //
11 
12 using namespace analysis;
13 using namespace analysis::tools;
14 
15 //
16 // constructors and destructor
17 //
19 {
20  q_ = 0;
21  p4_.SetPtEtaPhiE(0.,0.,0.,0.);
22 }
23 
24 Candidate::Candidate(const float & pt, const float & eta, const float & phi, const float & e, const float & q)
25 {
26  q_ = q;
27  p4_.SetPtEtaPhiE(pt,eta,phi,e);
28 }
29 
30 Candidate::Candidate(const float & px, const float & py, const float & pz)
31 {
32  q_ = 0;
33  p4_.SetXYZM(px,py,pz,0.);
34 }
35 
36 Candidate::Candidate(const float & px, const float & py, const float & pz, const float & q)
37 {
38  q_ = q;
39  p4_.SetXYZM(px,py,pz,0.);
40 }
41 
42 
44 {
45 }
46 
47 
48 //
49 // member functions
50 //
51 bool Candidate::matchTo(const std::vector<Candidate> * cands, const std::string & name, const float & deltaR)
52 {
53  bool status = false;
54 
55 
56  if ( ! cands )
57  {
58  this -> matched_[name] = nullptr;
59  return status;
60  }
61 
62  const Candidate * cand = nullptr;
63  const Candidate * nearcand = nullptr;
64  float minDeltaR = 100.;
65  for ( size_t i = 0; i < cands->size() ; ++i )
66  {
67  cand = &(cands->at(i));
68  if(this->deltaR(*cand) < minDeltaR)
69  {
70  minDeltaR = this->deltaR(*cand);
71  nearcand = cand;
72  }
73  }
74 
75  if(minDeltaR < deltaR)
76  {
77  this->matched_[name]=nearcand;
78  status = true;
79  }
80 
81  else {
82  this -> matched_[name] = nullptr;
83  }
84 
85  return status;
86 }
87 
88 bool Candidate::matchTo(const std::vector<Candidate> * cands, const std::string & name, const float & delta_pT, const float & deltaR)
89 {
90  bool status = false;
91 
92 
93  if ( ! cands )
94  {
95  this -> matched_[name] = nullptr;
96  return status;
97  }
98 
99  const Candidate * cand = nullptr;
100  const Candidate * nearcand = nullptr;
101  float minDeltaR = deltaR + 1; // Assign more real value;
102  float dpT = 0.;
103  float dpTmin = delta_pT + 1;
104  for ( size_t i = 0; i < cands->size() ; ++i )
105  {
106  cand = &(cands->at(i));
107  dpT = std::abs(this->pt() - cand->pt());
108  if(this->deltaR(*cand) < minDeltaR && dpT < dpTmin)
109  {
110  minDeltaR = this->deltaR(*cand);
111  dpTmin = dpT;
112  nearcand = cand;
113  }
114  }
115 
116  if(minDeltaR < deltaR && dpTmin < delta_pT)
117  {
118  this->matched_[name]=nearcand;
119  status = true;
120  }
121 
122  else {
123  this -> matched_[name] = nullptr;
124  }
125 
126  return status;
127 }
128 
129 // Gets
130 float Candidate::px() const { return p4_.Px() ; }
131 float Candidate::py() const { return p4_.Py() ; }
132 float Candidate::pz() const { return p4_.Pz() ; }
133 float Candidate::pt() const { return p4_.Pt() ; }
134 float Candidate::eta() const { return p4_.Eta(); }
135 float Candidate::phi() const { return p4_.Phi(); }
136 float Candidate::e() const { return p4_.E() ; }
137 float Candidate::m() const { return p4_.M() ; }
138 float Candidate::mass() const { return p4_.M() ; }
139 int Candidate::q() const { return q_; }
140 float Candidate::deltaR(const Candidate &cand) const { return p4_.DeltaR(cand.p4()) ;}
141 float Candidate::deltaPhi(const Candidate &cand) const { return p4_.DeltaPhi(cand.p4()) ;}
142 
143 TLorentzVector Candidate::p4() const { return p4_; }
144 TVector3 Candidate::p3() const { return p4_.Vect(); }
145 
146 
147 const Candidate * Candidate::matched(const std::string & name) { return matched_[name]; }
148 const Candidate * Candidate::matched(const std::string & name) const { return matched_.find(name) != matched_.end() ? matched_.find(name)->second : 0; }
149 
150 void Candidate::unmatch(const std::string & name)
151 {
152  matched_[name] = nullptr;
153 }
154 
155 // Sets
156 void Candidate::p4(const TLorentzVector & p4) { p4_ = p4; }
157 void Candidate::px(const float & px) { p4_.SetPx(px); }
158 void Candidate::py(const float & py) { p4_.SetPy(py); }
159 void Candidate::pz(const float & pz) { p4_.SetPy(pz); }
160 void Candidate::e (const float & e ) { p4_.SetE(e); }
161 void Candidate::q (const float & q) { q_ = q; }
TLorentzVector p4() const
returns the 4-momentum (TLorentzVector)
Definition: Candidate.cc:143
Candidate()
default constructor
Definition: Candidate.cc:18
float eta() const
returns the pseudorapidity
Definition: Candidate.cc:134
float mass() const
returns the mass
Definition: Candidate.cc:138
def status(submission_dir, failed_only=False)
Definition: naf_submit.py:362
float m() const
returns the mass
Definition: Candidate.cc:137
float deltaPhi(const Candidate &) const
returns the deltaPhi between this and another candidate
Definition: Candidate.cc:141
float e() const
returns the energy
Definition: Candidate.cc:136
std::map< std::string, const Candidate * > matched_
map of matched candidates
Definition: Candidate.h:116
float py() const
returns the y component of the momentum
Definition: Candidate.cc:131
float px() const
returns the x component of the momentum
Definition: Candidate.cc:130
float phi() const
returns the azimuthal angle
Definition: Candidate.cc:135
float pz() const
returns the z component of the momentum
Definition: Candidate.cc:132
const Candidate * matched(const std::string &name)
returns the pointer to the matched candidate object
Definition: Candidate.cc:147
float q_
the charge
Definition: Candidate.h:112
void unmatch(const std::string &)
unmatch a matched candidate object, i.e. set it to nullptr, useful to remove possible duplicated matc...
Definition: Candidate.cc:150
float pt() const
returns the transverse momentum
Definition: Candidate.cc:133
TLorentzVector p4_
the 4-momentum
Definition: Candidate.h:114
TVector3 p3() const
returns the 4-momentum (TVector3)
Definition: Candidate.cc:144
float deltaR(const Candidate &) const
returns the deltaR between this and another candidate
Definition: Candidate.cc:140
int q() const
returns the charge
Definition: Candidate.cc:139
virtual ~Candidate()
destructor
Definition: Candidate.cc:43
virtual bool matchTo(const std::vector< Candidate > *cands, const std::string &name, const float &deltaR=0.5)
function to match this candidate to another object from a list of pointers with a name ...
Definition: Candidate.cc:51