CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MET.cc
Go to the documentation of this file.
1 // File: MET.cc
2 // Description: see MET.h
3 // Author: Michael Schmitt, R. Cavanaugh University of Florida
4 // Creation Date: MHS MAY 30, 2005 initial version
5 
7 #include "TVector.h"
8 
9 using namespace std;
10 using namespace reco;
11 
12 //---------------------------------------------------------------------------
13 // Default Constructor;
14 //-----------------------------------
15 MET::MET()
16 {
17  sumet = 0.0;
18  elongit = 0.0;
19  signif_dxx=signif_dyy=signif_dyx=signif_dxy=0.;
20 }
21 //---------------------------------------------------------------------------
22 
23 //---------------------------------------------------------------------------
24 // Constructer for the case when only p4_ = (mEx, mEy, 0, mEt) is known.
25 // The vertex information is currently not used (but may be in the future)
26 // and is required by the RecoCandidate constructer.
27 //-----------------------------------
28 MET::MET( const LorentzVector& p4_, const Point& vtx_ ) :
29  RecoCandidate( 0, p4_, vtx_ )
30 {
31  sumet = 0.0;
32  elongit = 0.0;
34 }
35 //---------------------------------------------------------------------------
36 
37 //---------------------------------------------------------------------------
38 // Constructer for the case when the SumET is known in addition to
39 // p4_ = (mEx, mEy, 0, mEt). The vertex information is currently not used
40 // (but see above).
41 //-----------------------------------
42 MET::MET( double sumet_, const LorentzVector& p4_, const Point& vtx_ ) :
43  RecoCandidate( 0, p4_, vtx_ )
44 {
45  sumet = sumet_;
46  elongit = 0.0;
48 }
49 //---------------------------------------------------------------------------
50 
51 //---------------------------------------------------------------------------
52 // Constructor for the case when the SumET, the corrections which
53 // were applied to the MET, as well the MET itself p4_ = (mEx, mEy, 0, mEt)
54 // are all known. See above concerning the vertex information.
55 //-----------------------------------
56 MET::MET( double sumet_, const std::vector<CorrMETData>& corr_,
57  const LorentzVector& p4_, const Point& vtx_ ) :
58  RecoCandidate( 0, p4_, vtx_ )
59 {
60  sumet = sumet_;
61  elongit = 0.0;
63  //-----------------------------------
64  // Fill the vector containing the corrections (corr) with vector of
65  // known corrections (corr_) passed in via the constructor.
66  std::vector<CorrMETData>::const_iterator i;
67  for( i = corr_.begin(); i != corr_.end(); i++ )
68  {
69  corr.push_back( *i );
70  }
71 }
72 //---------------------------------------------------------------------------
73 
74 
75 //-----------------------------------------------------------------
76 //explicit clone function
77 MET * MET::clone() const {
78  return new MET( * this );
79 }
80 //-----------------------------------------------------------------
81 
82 // function that calculates the MET significance from the vector information.
83 double MET::significance() const {
84  if(signif_dxx==0 && signif_dyy==0 && signif_dxy==0 && signif_dyx==0)
85  return -1;
86  TMatrixD metmat = getSignificanceMatrix();
87  TVectorD metvec(2);
88  metvec(0)=this->px();
89  metvec(1)=this->py();
90  double signif = -1;
91  if(std::fabs(metmat.Determinant())>0.000001){
92  metmat.Invert();
93  signif = metvec * (metmat * metvec);
94  }
95  return signif;
96 }
97 
98 //---------------------------------------------------------------------------
99 // Returns the vector of all corrections applied to the x component of the
100 // missing transverse momentum, mEx
101 //-----------------------------------
102 std::vector<double> MET::dmEx() const
103 {
104  std::vector<double> deltas;
105  std::vector<CorrMETData>::const_iterator i;
106  for( i = corr.begin(); i != corr.end(); i++ )
107  {
108  deltas.push_back( i->mex );
109  }
110  return deltas;
111 }
112 //---------------------------------------------------------------------------
113 
114 //---------------------------------------------------------------------------
115 // Returns the vector of all corrections applied to the y component of the
116 // missing transverse momentum, mEy
117 //-----------------------------------
118 std::vector<double> MET::dmEy() const
119 {
120  std::vector<double> deltas;
121  std::vector<CorrMETData>::const_iterator i;
122  for( i = corr.begin(); i != corr.end(); i++ )
123  {
124  deltas.push_back( i->mey );
125  }
126  return deltas;
127 }
128 //---------------------------------------------------------------------------
129 
130 //---------------------------------------------------------------------------
131 // Returns the vector of all corrections applied to the scalar sum of the
132 // transverse energy (over all objects)
133 //-----------------------------------
134 std::vector<double> MET::dsumEt() const
135 {
136  std::vector<double> deltas;
137  std::vector<CorrMETData>::const_iterator i;
138  for( i = corr.begin(); i != corr.end(); i++ )
139  {
140  deltas.push_back( i->sumet );
141  }
142  return deltas;
143 }
144 
145 //---------------------------------------------------------------------------
146 // returns the significance matrix
147 //---------------------------------------------------------------------------
148 
149 TMatrixD MET::getSignificanceMatrix(void) const
150 {
151  TMatrixD result(2,2);
152  result(0,0)=signif_dxx;
153  result(0,1)=signif_dxy;
154  result(1,0)=signif_dyx;
155  result(1,1)=signif_dyy;
156  return result;
157 }
158 
159 //---------------------------------------------------------------------------
160 // Returns the vector of all corrections applied to the significance of the
161 // transverse energy (over all objects)
162 //-----------------------------------
163 std::vector<double> MET::dSignificance() const
164 {
165  std::vector<double> deltas;
166  std::vector<CorrMETData>::const_iterator i;
167  for( i = corr.begin(); i != corr.end(); i++ )
168  {
169  deltas.push_back( i->significance );
170  }
171  return deltas;
172 }
173 //---------------------------------------------------------------------------
174 
175 //---------------------------------------------------------------------------
176 // Required RecoCandidate polymorphism
177 //-----------------------------------
178 bool MET::overlap( const Candidate & ) const
179 {
180  return false;
181 }
182 //---------------------------------------------------------------------------
183 
184 void
185 MET::setSignificanceMatrix(const TMatrixD &inmatrix){
186  signif_dxx=inmatrix(0,0);
187  signif_dxy=inmatrix(0,1);
188  signif_dyx=inmatrix(1,0);
189  signif_dyy=inmatrix(1,1);
190  return;
191 }
std::vector< double > dsumEt() const
Definition: MET.cc:134
int i
Definition: DBlmapReader.cc:9
double signif_dyx
Definition: MET.h:74
MET * clone() const
returns a clone of the Candidate object
Definition: MET.cc:77
double signif_dxy
Definition: MET.h:75
double signif_dxx
Definition: MET.h:72
std::vector< double > dmEx() const
Definition: MET.cc:102
MET()
Definition: MET.cc:15
virtual double py() const GCC11_FINAL
y coordinate of momentum vector
TMatrixD getSignificanceMatrix(void) const
Definition: MET.cc:149
Definition: MET.h:32
tuple result
Definition: query.py:137
std::vector< CorrMETData > corr
Definition: MET.h:76
virtual double px() const GCC11_FINAL
x coordinate of momentum vector
void setSignificanceMatrix(const TMatrixD &matrix)
Definition: MET.cc:185
std::vector< double > dmEy() const
Definition: MET.cc:118
double signif_dyy
Definition: MET.h:73
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
virtual bool overlap(const Candidate &) const
check overlap with another candidate
Definition: MET.cc:178
std::vector< double > dSignificance() const
Definition: MET.cc:163
math::XYZPoint Point
point in the space
Definition: Candidate.h:45
double elongit
Definition: MET.h:70
double significance() const
Definition: MET.cc:83
double sumet
Definition: MET.h:69