Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:04:45

0001 // -*- C++ -*-
0002 
0003 // Package:    METReco
0004 // Class:      MET
0005 //
0006 // Original authors: Michael Schmitt, Richard Cavanaugh The University of Florida
0007 // changes by: Freya Blekman, Cornell University
0008 //
0009 
0010 //____________________________________________________________________________||
0011 #include "DataFormats/METReco/interface/MET.h"
0012 #include "TVector.h"
0013 
0014 //____________________________________________________________________________||
0015 using namespace std;
0016 using namespace reco;
0017 
0018 //____________________________________________________________________________||
0019 MET::MET() {
0020   sumet = 0.0;
0021   elongit = 0.0;
0022   signif_dxx = signif_dyy = signif_dyx = signif_dxy = 0.;
0023   mIsWeighted = false;
0024 }
0025 
0026 // Constructer for the case when only p4_ =  (mEx, mEy, 0, mEt) is known.
0027 // The vertex information is currently not used (but may be in the future)
0028 // and is required by the RecoCandidate constructer.
0029 //____________________________________________________________________________||
0030 MET::MET(const LorentzVector& p4_, const Point& vtx_, bool isWeighted) : RecoCandidate(0, p4_, vtx_) {
0031   sumet = 0.0;
0032   elongit = 0.0;
0033   signif_dxx = signif_dyy = signif_dyx = signif_dxy = 0.;
0034   mIsWeighted = isWeighted;
0035 }
0036 
0037 // Constructer for the case when the SumET is known in addition to
0038 // p4_ = (mEx, mEy, 0, mEt).  The vertex information is currently not used
0039 // (but see above).
0040 //____________________________________________________________________________||
0041 MET::MET(double sumet_, const LorentzVector& p4_, const Point& vtx_, bool isWeighted) : RecoCandidate(0, p4_, vtx_) {
0042   sumet = sumet_;
0043   elongit = 0.0;
0044   signif_dxx = signif_dyy = signif_dyx = signif_dxy = 0.;
0045   mIsWeighted = isWeighted;
0046 }
0047 
0048 // Constructor for the case when the SumET, the corrections which
0049 // were applied to the MET, as well the MET itself p4_ = (mEx, mEy, 0, mEt)
0050 // are all known.  See above concerning the vertex information.
0051 //____________________________________________________________________________||
0052 MET::MET(
0053     double sumet_, const std::vector<CorrMETData>& corr_, const LorentzVector& p4_, const Point& vtx_, bool isWeighted)
0054     : RecoCandidate(0, p4_, vtx_) {
0055   sumet = sumet_;
0056   elongit = 0.0;
0057   signif_dxx = signif_dyy = signif_dyx = signif_dxy = 0.;
0058   mIsWeighted = isWeighted;
0059   //-----------------------------------
0060   // Fill the vector containing the corrections (corr) with vector of
0061   // known corrections (corr_) passed in via the constructor.
0062   std::vector<CorrMETData>::const_iterator i;
0063   for (i = corr_.begin(); i != corr_.end(); i++) {
0064     corr.push_back(*i);
0065   }
0066 }
0067 
0068 //____________________________________________________________________________||
0069 MET* MET::clone() const { return new MET(*this); }
0070 
0071 // function that calculates the MET significance from the vector information.
0072 //____________________________________________________________________________||
0073 double MET::significance() const {
0074   if (signif_dxx == 0 && signif_dyy == 0 && signif_dxy == 0 && signif_dyx == 0)
0075     return -1;
0076   METCovMatrix metmat = getSignificanceMatrix();
0077   ROOT::Math::SVector<double, 2> metvec;
0078   metvec(0) = this->px();
0079   metvec(1) = this->py();
0080   double signif = -1;
0081   double det = 0;
0082   metmat.Det2(det);
0083   if (std::abs(det) > 0.000001) {
0084     metmat.Invert();
0085     signif = ROOT::Math::Dot(metvec, (metmat * metvec));
0086   }
0087   return signif;
0088 }
0089 
0090 // Returns the vector of all corrections applied to the x component of the
0091 // missing transverse momentum, mEx
0092 //____________________________________________________________________________||
0093 std::vector<double> MET::dmEx() const {
0094   std::vector<double> deltas;
0095   std::vector<CorrMETData>::const_iterator i;
0096   for (i = corr.begin(); i != corr.end(); i++) {
0097     deltas.push_back(i->mex);
0098   }
0099   return deltas;
0100 }
0101 
0102 // Returns the vector of all corrections applied to the y component of the
0103 // missing transverse momentum, mEy
0104 //____________________________________________________________________________||
0105 std::vector<double> MET::dmEy() const {
0106   std::vector<double> deltas;
0107   std::vector<CorrMETData>::const_iterator i;
0108   for (i = corr.begin(); i != corr.end(); i++) {
0109     deltas.push_back(i->mey);
0110   }
0111   return deltas;
0112 }
0113 
0114 // Returns the vector of all corrections applied to the scalar sum of the
0115 // transverse energy (over all objects)
0116 //____________________________________________________________________________||
0117 std::vector<double> MET::dsumEt() const {
0118   std::vector<double> deltas;
0119   std::vector<CorrMETData>::const_iterator i;
0120   for (i = corr.begin(); i != corr.end(); i++) {
0121     deltas.push_back(i->sumet);
0122   }
0123   return deltas;
0124 }
0125 
0126 // returns the significance matrix
0127 //____________________________________________________________________________||
0128 METCovMatrix MET::getSignificanceMatrix(void) const {
0129   METCovMatrix result;
0130   result(0, 0) = signif_dxx;
0131   result(0, 1) = signif_dxy;
0132   result(1, 0) = signif_dyx;
0133   result(1, 1) = signif_dyy;
0134   return result;
0135 }
0136 
0137 // Required RecoCandidate polymorphism
0138 //____________________________________________________________________________||
0139 bool MET::overlap(const Candidate&) const { return false; }
0140 
0141 //____________________________________________________________________________||
0142 void MET::setSignificanceMatrix(const METCovMatrix& inmatrix) {
0143   signif_dxx = inmatrix(0, 0);
0144   signif_dxy = inmatrix(0, 1);
0145   signif_dyx = inmatrix(1, 0);
0146   signif_dyy = inmatrix(1, 1);
0147   return;
0148 }
0149 
0150 //____________________________________________________________________________||