CMS 3D CMS Logo

significanceAlgo.cc
Go to the documentation of this file.
2 #include <iostream>
3 #include <string>
4 #include "TMath.h"
5 #include <cmath>
6 #include "TROOT.h"
8 
9 // -*- C++ -*-
10 //
11 // Package: METAlgorithms
12 // Class: SignAlgoResolutions
13 //
21 //
22 // Original Author: Kyle Story, Freya Blekman (Cornell University)
23 // Created: Fri Apr 18 11:58:33 CEST 2008
24 //
25 //
26 
28  // eventVec_(0),
29  set_worker_(0),
30  xmet_(0),
31  ymet_(0)
32 
33 {
34  // std::cout << "in constructor ! " << std::endl;
36 }
37 
38 
39 //******* 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)!
40 const void
43  return;
44 }
47 
48 const void
49 metsig::significanceAlgo::setSignifMatrix(const reco::METCovMatrix &input,const double &met_r, const double &met_phi, const double &met_set){
51  set_worker_=met_set;
52  xmet_=met_r*cos(met_phi);
53  ymet_=met_r*sin(met_phi);
54 
55  return;
56 }
57 
58 // ********destructor ********
59 
61 }
62 
63 // //*** rotate a 2D matrix by angle theta **********************//
64 
65 // void
66 // metsig::significanceAlgo::rotateMatrix( Double_t theta, reco::METCovMatrix &v)
67 // {
68 // // I suggest not using this to rotate trivial matrices.
69 // reco::METCovMatrix r;
70 // reco::METCovMatrix rInv;
71 
72 // r(0,0) = cos(theta); r(0,1) = sin(theta); r(1,0) = -sin(theta); r(1,1) = cos(theta);
73 // rInv = r;
74 // rInv.Invert();
75 // //-- Rotate v --//
76 // v = rInv * v * r;
77 // }
78 //************************************************************//
79 
80 
81 const void
82 metsig::significanceAlgo::subtractObjects(const std::vector<metsig::SigInputObj>& eventVec)
83 {
85  //--- Loop over physics objects in the event ---//
86  // for(unsigned int objnum=1; objnum < EventVec.size(); objnum++ ) {
87  for(std::vector<SigInputObj>::const_iterator obj = eventVec.begin(); obj!= eventVec.end(); ++obj){
88  double et_tmp = obj->get_energy();
89  double phi_tmp = obj->get_phi();
90  double sigma_et = obj->get_sigma_e();
91  double sigma_tan = obj->get_sigma_tan();
92 
93  double cosphi=cos(phi_tmp);
94  double sinphi=sin(phi_tmp);
95  double xval = et_tmp * cosphi;
96  double yval = et_tmp * sinphi;
97  xmet_ += xval;
98  ymet_ += yval;
99  set_worker_ -= et_tmp;
100 
101  double sigma0_2=sigma_et*sigma_et;
102  double sigma1_2=sigma_tan*sigma_tan;
103 
104  v_tot(0,0)-= sigma0_2*cosphi*cosphi + sigma1_2*sinphi*sinphi;
105  v_tot(0,1)-= cosphi*sinphi*(sigma0_2 - sigma1_2);
106  v_tot(1,0)-= cosphi*sinphi*(sigma0_2 - sigma1_2);
107  v_tot(1,1)-= sigma1_2*cosphi*cosphi + sigma0_2*sinphi*sinphi;
108 
109  }
110  signifmatrix_=v_tot;
111 }
112 //************************************************************//
113 
114 
115 const void
116 metsig::significanceAlgo::addObjects(const std::vector<metsig::SigInputObj>& eventVec)
117 {
119  //--- Loop over physics objects in the event ---//
120  // for(unsigned int objnum=1; objnum < EventVec.size(); objnum++ ) {
121  for(std::vector<SigInputObj>::const_iterator obj = eventVec.begin(); obj!= eventVec.end(); ++obj){
122  double et_tmp = obj->get_energy();
123  double phi_tmp = obj->get_phi();
124  double sigma_et = obj->get_sigma_e();
125  double sigma_tan = obj->get_sigma_tan();
126 
127  double cosphi=cos(phi_tmp);
128  double sinphi=sin(phi_tmp);
129  double xval = et_tmp * cosphi;
130  double yval = et_tmp * sinphi;
131  xmet_ -= xval;
132  ymet_ -= yval;
133  set_worker_ += et_tmp;
134 
135  double sigma0_2=sigma_et*sigma_et;
136  double sigma1_2=sigma_tan*sigma_tan;
137 
138  v_tot(0,0)+= sigma0_2*cosphi*cosphi + sigma1_2*sinphi*sinphi;
139  v_tot(0,1)+= cosphi*sinphi*(sigma0_2 - sigma1_2);
140  v_tot(1,0)+= cosphi*sinphi*(sigma0_2 - sigma1_2);
141  v_tot(1,1)+= sigma1_2*cosphi*cosphi + sigma0_2*sinphi*sinphi;
142 
143  }
144  signifmatrix_=v_tot;
145 }
146 
147 //************************************************************//
148 const double
149 metsig::significanceAlgo::significance(double &met_r, double &met_phi, double &met_set)
150 {
151 
152  if(signifmatrix_(0,0)==0 && signifmatrix_(1,1)==0 && signifmatrix_(1,0)==0 && signifmatrix_(0,1)==0){
153  //edm::LogWarning("SignCaloSpecificAlgo") << "Event Vector is empty! Return significance -1";
154  return(-1);
155  }
156 
157  //--- Temporary variables ---//
158 
159  reco::METCovMatrix v_tot;
160  ROOT::Math::SVector<double, 2> metvec;
161 
162  //--- Initialize sum of rotated covariance matrices ---//
163  v_tot=signifmatrix_;
164  // std::cout << "INPUT:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl;
165 
166 
167 
168  //--- Calculate magnitude and angle of MET, store in returned variables ---//
169  met_r = sqrt(xmet_*xmet_ + ymet_*ymet_);
170  met_set = set_worker_;
171  met_phi= TMath::ATan2(ymet_, xmet_);
172 
173  // one other option: if particles cancel there could be small numbers.
174  // this check fixes this, added by F.Blekman
175  double det=0;
176  v_tot.Det(det);
177  if(fabs(det)<0.000001)
178  return -1;
179 
180 
181  // save matrix into object:
182  // std::cout << "SAVED:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl;
183  v_tot.Invert();
184  // std::cout << "INVERTED:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl;
185 
186 
187 
188  metvec(0) = xmet_; metvec(1) = ymet_;
189  double lnSignificance = ROOT::Math::Dot(metvec, (v_tot * metvec) );
190 
191  // v_tot.Invert();
192  // std::cout << "INVERTED AGAIN:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl;
193  return lnSignificance;
194 }
195 //*** End of Significance calculation ********************************//
const void addObjects(const std::vector< metsig::SigInputObj > &EventVec)
const double significance(double &met_r, double &met_phi, double &met_set)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const void setSignifMatrix(const reco::METCovMatrix &input, const double &met_r, const double &met_phi, const double &met_set)
reset the signficance matrix (this is the most likely case), so that the vector sum can be continued ...
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:40
reco::METCovMatrix signifmatrix_
const void subtractObjects(const std::vector< metsig::SigInputObj > &EventVec)
static std::string const input
Definition: EdmProvDump.cc:45
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const void addSignifMatrix(const reco::METCovMatrix &input)