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 // -*- C++ -*-
2 
3 // Package: METReco
4 // Class: MET
5 //
6 // Original authors: Michael Schmitt, Richard Cavanaugh The University of Florida
7 // changes by: Freya Blekman, Cornell University
8 //
9 
10 //____________________________________________________________________________||
12 #include "TVector.h"
13 
14 //____________________________________________________________________________||
15 using namespace std;
16 using namespace reco;
17 
18 //____________________________________________________________________________||
19 MET::MET()
20 {
21  sumet = 0.0;
22  elongit = 0.0;
23  signif_dxx=signif_dyy=signif_dyx=signif_dxy=0.;
24 }
25 
26 // Constructer for the case when only p4_ = (mEx, mEy, 0, mEt) is known.
27 // The vertex information is currently not used (but may be in the future)
28 // and is required by the RecoCandidate constructer.
29 //____________________________________________________________________________||
30 MET::MET( const LorentzVector& p4_, const Point& vtx_ ) :
31  RecoCandidate( 0, p4_, vtx_ )
32 {
33  sumet = 0.0;
34  elongit = 0.0;
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 // Constructor for the case when the SumET, the corrections which
51 // were applied to the MET, as well the MET itself p4_ = (mEx, mEy, 0, mEt)
52 // are all known. See above concerning the vertex information.
53 //____________________________________________________________________________||
54 MET::MET( double sumet_, const std::vector<CorrMETData>& corr_,
55  const LorentzVector& p4_, const Point& vtx_ ) :
56  RecoCandidate( 0, p4_, vtx_ )
57 {
58  sumet = sumet_;
59  elongit = 0.0;
61  //-----------------------------------
62  // Fill the vector containing the corrections (corr) with vector of
63  // known corrections (corr_) passed in via the constructor.
64  std::vector<CorrMETData>::const_iterator i;
65  for( i = corr_.begin(); i != corr_.end(); i++ )
66  {
67  corr.push_back( *i );
68  }
69 }
70 
71 //____________________________________________________________________________||
72 MET * MET::clone() const {
73  return new MET( * this );
74 }
75 
76 // function that calculates the MET significance from the vector information.
77 //____________________________________________________________________________||
78 double MET::significance() const {
79  if(signif_dxx==0 && signif_dyy==0 && signif_dxy==0 && signif_dyx==0)
80  return -1;
81  TMatrixD metmat = getSignificanceMatrix();
82  TVectorD metvec(2);
83  metvec(0)=this->px();
84  metvec(1)=this->py();
85  double signif = -1;
86  if(std::fabs(metmat.Determinant())>0.000001){
87  metmat.Invert();
88  signif = metvec * (metmat * metvec);
89  }
90  return signif;
91 }
92 
93 // Returns the vector of all corrections applied to the x component of the
94 // missing transverse momentum, mEx
95 //____________________________________________________________________________||
96 std::vector<double> MET::dmEx() const
97 {
98  std::vector<double> deltas;
99  std::vector<CorrMETData>::const_iterator i;
100  for( i = corr.begin(); i != corr.end(); i++ )
101  {
102  deltas.push_back( i->mex );
103  }
104  return deltas;
105 }
106 
107 // Returns the vector of all corrections applied to the y component of the
108 // missing transverse momentum, mEy
109 //____________________________________________________________________________||
110 std::vector<double> MET::dmEy() const
111 {
112  std::vector<double> deltas;
113  std::vector<CorrMETData>::const_iterator i;
114  for( i = corr.begin(); i != corr.end(); i++ )
115  {
116  deltas.push_back( i->mey );
117  }
118  return deltas;
119 }
120 
121 // Returns the vector of all corrections applied to the scalar sum of the
122 // transverse energy (over all objects)
123 //____________________________________________________________________________||
124 std::vector<double> MET::dsumEt() const
125 {
126  std::vector<double> deltas;
127  std::vector<CorrMETData>::const_iterator i;
128  for( i = corr.begin(); i != corr.end(); i++ )
129  {
130  deltas.push_back( i->sumet );
131  }
132  return deltas;
133 }
134 
135 // returns the significance matrix
136 //____________________________________________________________________________||
137 TMatrixD MET::getSignificanceMatrix(void) const
138 {
139  TMatrixD result(2,2);
140  result(0,0)=signif_dxx;
141  result(0,1)=signif_dxy;
142  result(1,0)=signif_dyx;
143  result(1,1)=signif_dyy;
144  return result;
145 }
146 
147 // Required RecoCandidate polymorphism
148 //____________________________________________________________________________||
149 bool MET::overlap( const Candidate & ) const
150 {
151  return false;
152 }
153 
154 //____________________________________________________________________________||
155 void MET::setSignificanceMatrix(const TMatrixD &inmatrix)
156 {
157  signif_dxx=inmatrix(0,0);
158  signif_dxy=inmatrix(0,1);
159  signif_dyx=inmatrix(1,0);
160  signif_dyy=inmatrix(1,1);
161  return;
162 }
163 
164 //____________________________________________________________________________||
std::vector< double > dsumEt() const
Definition: MET.cc:124
int i
Definition: DBlmapReader.cc:9
double signif_dyx
Definition: MET.h:80
MET * clone() const
returns a clone of the Candidate object
Definition: MET.cc:72
double signif_dxy
Definition: MET.h:81
double signif_dxx
Definition: MET.h:78
std::vector< double > dmEx() const
Definition: MET.cc:96
MET()
Definition: MET.cc:19
TMatrixD getSignificanceMatrix(void) const
Definition: MET.cc:137
Definition: MET.h:39
tuple result
Definition: query.py:137
std::vector< CorrMETData > corr
Definition: MET.h:82
void setSignificanceMatrix(const TMatrixD &matrix)
Definition: MET.cc:155
std::vector< double > dmEy() const
Definition: MET.cc:110
double signif_dyy
Definition: MET.h:79
virtual double px() const
x coordinate of momentum vector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
virtual bool overlap(const Candidate &) const
check overlap with another candidate
Definition: MET.cc:149
math::XYZPoint Point
point in the space
Definition: Candidate.h:45
double elongit
Definition: MET.h:76
double significance() const
Definition: MET.cc:78
double sumet
Definition: MET.h:75
virtual double py() const
y coordinate of momentum vector