CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/RecoMET/METAlgorithms/src/significanceAlgo.cc

Go to the documentation of this file.
00001 #include "RecoMET/METAlgorithms/interface/significanceAlgo.h"
00002 #include <iostream>
00003 #include <string>
00004 #include "TMath.h"
00005 #include "math.h"
00006 #include "TROOT.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 
00009 // -*- C++ -*-
00010 //
00011 // Package:    METAlgorithms
00012 // Class:      SignAlgoResolutions
00013 // 
00021 //
00022 // Original Author:  Kyle Story, Freya Blekman (Cornell University)
00023 //         Created:  Fri Apr 18 11:58:33 CEST 2008
00024 // $Id: significanceAlgo.cc,v 1.15 2010/04/20 14:56:07 fblekman Exp $
00025 //
00026 //
00027 
00028 metsig::significanceAlgo::significanceAlgo():
00029   //  eventVec_(0),
00030   signifmatrix_(2,2),
00031   set_worker_(0),
00032   xmet_(0),
00033   ymet_(0)
00034 
00035 {
00036   //  std::cout << "in constructor ! " << std::endl;
00037   signifmatrix_(0,0)=signifmatrix_(1,0)=signifmatrix_(0,1)=signifmatrix_(1,1)=0;
00038 }
00039 
00040 
00041 //******* Add an existing significance matrix to the algo, so that the vector sum can be continued. Only makes sense if matrix is empty or you want to purposefully increase uncertainties (for example in systematic studies)!
00042 const void
00043 metsig::significanceAlgo::addSignifMatrix(const TMatrixD &input){
00044   // check that the matrix is size 2:
00045   if(input.GetNrows()==2 && input.GetNcols()==2) {
00046     signifmatrix_+=input;
00047   }
00048   return;
00049 }
00052 
00053 const void
00054 metsig::significanceAlgo::setSignifMatrix(const TMatrixD &input,const double &met_r, const double &met_phi, const double &met_set){
00055   // check that the matrix is size 2:
00056   if(input.GetNrows()==2 && input.GetNcols()==2) {
00057     signifmatrix_=input;
00058     set_worker_=met_set;
00059     xmet_=met_r*cos(met_phi);
00060     ymet_=met_r*sin(met_phi);
00061   }
00062   return;
00063 }
00064 
00065 // ********destructor ********
00066 
00067 metsig::significanceAlgo::~significanceAlgo(){
00068 }
00069 
00070 // //*** rotate a 2D matrix by angle theta **********************//
00071 
00072 void
00073 metsig::significanceAlgo::rotateMatrix( Double_t theta, TMatrixD &v)
00074 {
00075   // I suggest not using this to rotate trivial matrices.
00076   TMatrixD r(2,2);
00077   TMatrixD rInv(2,2);
00078 
00079   r(0,0) = cos(theta); r(0,1) = sin(theta); r(1,0) = -sin(theta); r(1,1) = cos(theta);
00080   rInv = r;
00081   rInv.Invert();
00082   //-- Rotate v --//
00083   v = rInv * v * r;
00084 }
00085 //************************************************************//
00086 
00087 
00088 const void 
00089 metsig::significanceAlgo::subtractObjects(const std::vector<metsig::SigInputObj>& eventVec)
00090 { 
00091   TMatrixD v_tot = signifmatrix_;
00092   //--- Loop over physics objects in the event ---//
00093   //  for(unsigned int objnum=1; objnum < EventVec.size(); objnum++ ) {
00094   for(std::vector<SigInputObj>::const_iterator obj = eventVec.begin(); obj!= eventVec.end(); ++obj){
00095     double et_tmp     = obj->get_energy();
00096     double phi_tmp    = obj->get_phi();
00097     double sigma_et   = obj->get_sigma_e();
00098     double sigma_tan  = obj->get_sigma_tan();
00099 
00100     double cosphi=cos(phi_tmp);
00101     double sinphi=sin(phi_tmp);
00102     double xval = et_tmp * cosphi;
00103     double yval = et_tmp * sinphi;
00104     xmet_ += xval;
00105     ymet_ += yval;
00106     set_worker_ -= et_tmp;
00107 
00108     double sigma0_2=sigma_et*sigma_et;
00109     double sigma1_2=sigma_tan*sigma_tan;
00110 
00111     v_tot(0,0)-= sigma0_2*cosphi*cosphi + sigma1_2*sinphi*sinphi;
00112     v_tot(0,1)-= cosphi*sinphi*(sigma0_2 - sigma1_2);
00113     v_tot(1,0)-= cosphi*sinphi*(sigma0_2 - sigma1_2);
00114     v_tot(1,1)-= sigma1_2*cosphi*cosphi + sigma0_2*sinphi*sinphi;
00115     
00116   }
00117   signifmatrix_=v_tot;
00118 }
00119 //************************************************************//
00120 
00121 
00122 const void 
00123 metsig::significanceAlgo::addObjects(const std::vector<metsig::SigInputObj>& eventVec)
00124 { 
00125   TMatrixD v_tot = signifmatrix_;
00126   //--- Loop over physics objects in the event ---//
00127   //  for(unsigned int objnum=1; objnum < EventVec.size(); objnum++ ) {
00128   for(std::vector<SigInputObj>::const_iterator obj = eventVec.begin(); obj!= eventVec.end(); ++obj){
00129     double et_tmp     = obj->get_energy();
00130     double phi_tmp    = obj->get_phi();
00131     double sigma_et   = obj->get_sigma_e();
00132     double sigma_tan  = obj->get_sigma_tan();
00133 
00134     double cosphi=cos(phi_tmp);
00135     double sinphi=sin(phi_tmp);
00136     double xval = et_tmp * cosphi;
00137     double yval = et_tmp * sinphi;
00138     xmet_ -= xval;
00139     ymet_ -= yval;
00140     set_worker_ += et_tmp;
00141 
00142     double sigma0_2=sigma_et*sigma_et;
00143     double sigma1_2=sigma_tan*sigma_tan;
00144 
00145     v_tot(0,0)+= sigma0_2*cosphi*cosphi + sigma1_2*sinphi*sinphi;
00146     v_tot(0,1)+= cosphi*sinphi*(sigma0_2 - sigma1_2);
00147     v_tot(1,0)+= cosphi*sinphi*(sigma0_2 - sigma1_2);
00148     v_tot(1,1)+= sigma1_2*cosphi*cosphi + sigma0_2*sinphi*sinphi;
00149     
00150   }
00151   signifmatrix_=v_tot;
00152 }
00153 
00154 //************************************************************//
00155 const double 
00156 metsig::significanceAlgo::significance(double &met_r, double &met_phi, double &met_set) 
00157 {
00158   
00159   if(signifmatrix_(0,0)==0 && signifmatrix_(1,1)==0 && signifmatrix_(1,0)==0 && signifmatrix_(0,1)==0){
00160     //edm::LogWarning("SignCaloSpecificAlgo") << "Event Vector is empty!  Return significance -1";
00161     return(-1);
00162   } 
00163 
00164   //--- Temporary variables ---//
00165  
00166   TMatrixD v_tot(2,2);
00167   TVectorD metvec(2);
00168 
00169  //--- Initialize sum of rotated covariance matrices ---//  
00170   v_tot=signifmatrix_;
00171   //  std::cout << "INPUT:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl;
00172 
00173 
00174 
00175   //--- Calculate magnitude and angle of MET, store in returned variables ---//
00176   met_r = sqrt(xmet_*xmet_ + ymet_*ymet_);
00177   met_set = set_worker_;
00178   met_phi= TMath::ATan2(ymet_, xmet_);
00179 
00180   // one other option: if particles cancel there could be small numbers.
00181   // this check fixes this, added by F.Blekman
00182   if(fabs(v_tot.Determinant())<0.000001)
00183     return -1;
00184 
00185 
00186   // save matrix into object:
00187   //  std::cout << "SAVED:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl;
00188   v_tot.Invert();
00189   //  std::cout << "INVERTED:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl;
00190   
00191 
00192 
00193   metvec(0) = xmet_; metvec(1) = ymet_;
00194   double lnSignificance = metvec * (v_tot * metvec);
00195 
00196   //  v_tot.Invert();
00197   //  std::cout << "INVERTED AGAIN:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl;
00198   return lnSignificance;
00199 }
00200 //*** End of Significance calculation ********************************//