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