CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/DataFormats/METReco/src/MET.cc

Go to the documentation of this file.
00001 // File: MET.cc
00002 // Description: see MET.h
00003 // Author: Michael Schmitt, R. Cavanaugh University of Florida
00004 // Creation Date:  MHS MAY 30, 2005 initial version
00005 
00006 #include "DataFormats/METReco/interface/MET.h"
00007 #include "TVector.h"
00008 
00009 using namespace std;
00010 using namespace reco;
00011 
00012 //---------------------------------------------------------------------------
00013 // Default Constructor;
00014 //-----------------------------------
00015 MET::MET()
00016 {
00017   sumet = 0.0;
00018   elongit = 0.0;
00019   signif_dxx=signif_dyy=signif_dyx=signif_dxy=0.;
00020 }
00021 //---------------------------------------------------------------------------
00022 
00023 //---------------------------------------------------------------------------
00024 // Constructer for the case when only p4_ =  (mEx, mEy, 0, mEt) is known.
00025 // The vertex information is currently not used (but may be in the future)
00026 // and is required by the RecoCandidate constructer.
00027 //-----------------------------------
00028 MET::MET( const LorentzVector& p4_, const Point& vtx_ ) : 
00029   RecoCandidate( 0, p4_, vtx_ )
00030 {
00031   sumet = 0.0;
00032   elongit = 0.0;
00033   signif_dxx=signif_dyy=signif_dyx=signif_dxy=0.;
00034 }
00035 //---------------------------------------------------------------------------
00036 
00037 //---------------------------------------------------------------------------
00038 // Constructer for the case when the SumET is known in addition to 
00039 // p4_ = (mEx, mEy, 0, mEt).  The vertex information is currently not used
00040 // (but see above).
00041 //-----------------------------------
00042 MET::MET( double sumet_, const LorentzVector& p4_, const Point& vtx_ ) : 
00043   RecoCandidate( 0, p4_, vtx_ ) 
00044 {
00045   sumet = sumet_;
00046   elongit = 0.0;
00047   signif_dxx=signif_dyy=signif_dyx=signif_dxy=0.;
00048 }
00049 //---------------------------------------------------------------------------
00050 
00051 //---------------------------------------------------------------------------
00052 // Constructor for the case when the SumET, the corrections which
00053 // were applied to the MET, as well the MET itself p4_ = (mEx, mEy, 0, mEt)
00054 // are all known.  See above concerning the vertex information. 
00055 //-----------------------------------
00056 MET::MET( double sumet_, std::vector<CorrMETData> corr_, 
00057           const LorentzVector& p4_, const Point& vtx_ ) : 
00058   RecoCandidate( 0, p4_, vtx_ )
00059 {
00060   sumet = sumet_;
00061   elongit = 0.0;
00062   signif_dxx=signif_dyy=signif_dyx=signif_dxy=0.;
00063   //-----------------------------------
00064   // Fill the vector containing the corrections (corr) with vector of 
00065   // known corrections (corr_) passed in via the constructor.
00066   std::vector<CorrMETData>::const_iterator i;
00067   for( i = corr_.begin(); i != corr_.end();  i++ ) 
00068     {
00069       corr.push_back( *i );
00070     }
00071 }
00072 //---------------------------------------------------------------------------
00073 
00074 
00075 //----------------------------------------------------------------- 
00076 //explicit clone function
00077 MET * MET::clone() const {
00078      return new MET( * this );
00079 }
00080 //----------------------------------------------------------------- 
00081 
00082 // function that calculates the MET significance from the vector information.
00083 double MET::significance() const {
00084   if(signif_dxx==0 && signif_dyy==0 && signif_dxy==0 && signif_dyx==0)
00085     return -1;
00086   TMatrixD metmat = getSignificanceMatrix();
00087   TVectorD metvec(2);
00088   metvec(0)=this->px();
00089   metvec(1)=this->py();
00090   double signif = -1;
00091   if(std::fabs(metmat.Determinant())>0.000001){
00092     metmat.Invert();
00093     signif = metvec * (metmat * metvec);
00094   }
00095   return signif;
00096 }
00097 
00098 //---------------------------------------------------------------------------
00099 // Returns the vector of all corrections applied to the x component of the
00100 // missing transverse momentum, mEx
00101 //-----------------------------------
00102 std::vector<double> MET::dmEx() const 
00103 {
00104   std::vector<double> deltas;
00105   std::vector<CorrMETData>::const_iterator i;
00106   for( i = corr.begin(); i != corr.end(); i++ )
00107     {
00108       deltas.push_back( i->mex );
00109     }
00110   return deltas;
00111 }
00112 //---------------------------------------------------------------------------
00113 
00114 //---------------------------------------------------------------------------
00115 // Returns the vector of all corrections applied to the y component of the
00116 // missing transverse momentum, mEy
00117 //-----------------------------------
00118 std::vector<double> MET::dmEy() const 
00119 {
00120   std::vector<double> deltas;
00121   std::vector<CorrMETData>::const_iterator i;
00122   for( i = corr.begin(); i != corr.end(); i++ )
00123     {
00124       deltas.push_back( i->mey );
00125     }
00126   return deltas;
00127 }
00128 //---------------------------------------------------------------------------
00129 
00130 //---------------------------------------------------------------------------
00131 // Returns the vector of all corrections applied to the scalar sum of the 
00132 // transverse energy (over all objects)
00133 //-----------------------------------
00134 std::vector<double> MET::dsumEt() const 
00135 {
00136   std::vector<double> deltas;
00137   std::vector<CorrMETData>::const_iterator i;
00138   for( i = corr.begin(); i != corr.end(); i++ )
00139     {
00140       deltas.push_back( i->sumet );
00141     }
00142   return deltas;
00143 }
00144 
00145 //---------------------------------------------------------------------------
00146 // returns the significance matrix
00147 //---------------------------------------------------------------------------
00148 
00149 TMatrixD MET::getSignificanceMatrix(void) const
00150 {
00151   TMatrixD result(2,2);
00152   result(0,0)=signif_dxx;
00153   result(0,1)=signif_dxy;
00154   result(1,0)=signif_dyx;
00155   result(1,1)=signif_dyy;
00156   return result;
00157 }
00158 
00159 //---------------------------------------------------------------------------
00160 // Returns the vector of all corrections applied to the significance of the 
00161 // transverse energy (over all objects)
00162 //-----------------------------------
00163 std::vector<double> MET::dSignificance() const
00164 {
00165   std::vector<double> deltas;
00166   std::vector<CorrMETData>::const_iterator i;
00167   for( i = corr.begin(); i != corr.end(); i++ )
00168     {
00169       deltas.push_back( i->significance );
00170     }
00171   return deltas;
00172 }
00173 //---------------------------------------------------------------------------
00174 
00175 //---------------------------------------------------------------------------
00176 // Required RecoCandidate polymorphism
00177 //-----------------------------------
00178 bool MET::overlap( const Candidate & ) const 
00179 {
00180   return false;
00181 }
00182 //---------------------------------------------------------------------------
00183 
00184 void
00185 MET::setSignificanceMatrix(const TMatrixD &inmatrix){
00186   signif_dxx=inmatrix(0,0);
00187   signif_dxy=inmatrix(0,1);
00188   signif_dyx=inmatrix(1,0);
00189   signif_dyy=inmatrix(1,1);
00190   return;
00191 }