CMS 3D CMS Logo

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 //____________________________________________________________________________||
20  sumet = 0.0;
21  elongit = 0.0;
22  signif_dxx = signif_dyy = signif_dyx = signif_dxy = 0.;
23  mIsWeighted = false;
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_, bool isWeighted) : RecoCandidate(0, p4_, vtx_) {
31  sumet = 0.0;
32  elongit = 0.0;
35 }
36 
37 // Constructer for the case when the SumET is known in addition to
38 // p4_ = (mEx, mEy, 0, mEt). The vertex information is currently not used
39 // (but see above).
40 //____________________________________________________________________________||
41 MET::MET(double sumet_, const LorentzVector& p4_, const Point& vtx_, bool isWeighted) : RecoCandidate(0, p4_, vtx_) {
42  sumet = sumet_;
43  elongit = 0.0;
46 }
47 
48 // Constructor for the case when the SumET, the corrections which
49 // were applied to the MET, as well the MET itself p4_ = (mEx, mEy, 0, mEt)
50 // are all known. See above concerning the vertex information.
51 //____________________________________________________________________________||
53  double sumet_, const std::vector<CorrMETData>& corr_, const LorentzVector& p4_, const Point& vtx_, bool isWeighted)
54  : RecoCandidate(0, p4_, vtx_) {
55  sumet = sumet_;
56  elongit = 0.0;
59  //-----------------------------------
60  // Fill the vector containing the corrections (corr) with vector of
61  // known corrections (corr_) passed in via the constructor.
62  std::vector<CorrMETData>::const_iterator i;
63  for (i = corr_.begin(); i != corr_.end(); i++) {
64  corr.push_back(*i);
65  }
66 }
67 
68 //____________________________________________________________________________||
69 MET* MET::clone() const { return new MET(*this); }
70 
71 // function that calculates the MET significance from the vector information.
72 //____________________________________________________________________________||
73 double MET::significance() const {
74  if (signif_dxx == 0 && signif_dyy == 0 && signif_dxy == 0 && signif_dyx == 0)
75  return -1;
77  ROOT::Math::SVector<double, 2> metvec;
78  metvec(0) = this->px();
79  metvec(1) = this->py();
80  double signif = -1;
81  double det = 0;
82  metmat.Det2(det);
83  if (std::abs(det) > 0.000001) {
84  metmat.Invert();
85  signif = ROOT::Math::Dot(metvec, (metmat * metvec));
86  }
87  return signif;
88 }
89 
90 // Returns the vector of all corrections applied to the x component of the
91 // missing transverse momentum, mEx
92 //____________________________________________________________________________||
93 std::vector<double> MET::dmEx() const {
94  std::vector<double> deltas;
95  std::vector<CorrMETData>::const_iterator i;
96  for (i = corr.begin(); i != corr.end(); i++) {
97  deltas.push_back(i->mex);
98  }
99  return deltas;
100 }
101 
102 // Returns the vector of all corrections applied to the y component of the
103 // missing transverse momentum, mEy
104 //____________________________________________________________________________||
105 std::vector<double> MET::dmEy() const {
106  std::vector<double> deltas;
107  std::vector<CorrMETData>::const_iterator i;
108  for (i = corr.begin(); i != corr.end(); i++) {
109  deltas.push_back(i->mey);
110  }
111  return deltas;
112 }
113 
114 // Returns the vector of all corrections applied to the scalar sum of the
115 // transverse energy (over all objects)
116 //____________________________________________________________________________||
117 std::vector<double> MET::dsumEt() const {
118  std::vector<double> deltas;
119  std::vector<CorrMETData>::const_iterator i;
120  for (i = corr.begin(); i != corr.end(); i++) {
121  deltas.push_back(i->sumet);
122  }
123  return deltas;
124 }
125 
126 // returns the significance matrix
127 //____________________________________________________________________________||
130  result(0, 0) = signif_dxx;
131  result(0, 1) = signif_dxy;
132  result(1, 0) = signif_dyx;
133  result(1, 1) = signif_dyy;
134  return result;
135 }
136 
137 // Required RecoCandidate polymorphism
138 //____________________________________________________________________________||
139 bool MET::overlap(const Candidate&) const { return false; }
140 
141 //____________________________________________________________________________||
143  signif_dxx = inmatrix(0, 0);
144  signif_dxy = inmatrix(0, 1);
145  signif_dyx = inmatrix(1, 0);
146  signif_dyy = inmatrix(1, 1);
147  return;
148 }
149 
150 //____________________________________________________________________________||
std::vector< double > dsumEt() const
Definition: MET.cc:117
double signif_dyx
Definition: MET.h:88
double signif_dxy
Definition: MET.h:89
std::vector< double > dmEy() const
Definition: MET.cc:105
reco::METCovMatrix getSignificanceMatrix(void) const
Definition: MET.cc:128
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:39
void setSignificanceMatrix(const reco::METCovMatrix &matrix)
Definition: MET.cc:142
double signif_dxx
Definition: MET.h:86
MET()
Definition: MET.cc:19
MET * clone() const override
returns a clone of the Candidate object
Definition: MET.cc:69
int isWeighted() const
boolean if weights were applied by algorithm (e.g. PUPPI weights)
Definition: MET.h:79
double px() const final
x coordinate of momentum vector
Definition: MET.h:41
std::vector< CorrMETData > corr
Definition: MET.h:90
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool overlap(const Candidate &) const override
check overlap with another candidate
Definition: MET.cc:139
double py() const final
y coordinate of momentum vector
double signif_dyy
Definition: MET.h:87
double significance() const
Definition: MET.cc:73
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
fixed size matrix
Structure Point Contains parameters of Gaussian fits to DMRs.
double elongit
Definition: MET.h:84
bool mIsWeighted
Definition: MET.h:91
double sumet
Definition: MET.h:83
std::vector< double > dmEx() const
Definition: MET.cc:93