#include <significanceAlgo.h>
Public Member Functions | |
const void | addObjects (const std::vector< metsig::SigInputObj > &EventVec) |
const void | addSignifMatrix (const TMatrixD &input) |
TMatrixD | getSignifMatrix () const |
const void | setSignifMatrix (const TMatrixD &input, const double &met_r, const double &met_phi, const double &met_set) |
const double | significance (double &met_r, double &met_phi, double &met_set) |
significanceAlgo () | |
const void | subtractObjects (const std::vector< metsig::SigInputObj > &EventVec) |
~significanceAlgo () | |
Private Member Functions | |
void | rotateMatrix (Double_t theta, TMatrixD &v) |
Private Attributes | |
double | set_worker_ |
TMatrixD | signifmatrix_ |
double | xmet_ |
double | ymet_ |
Definition at line 88 of file significanceAlgo.h.
metsig::significanceAlgo::significanceAlgo | ( | ) |
Definition at line 28 of file significanceAlgo.cc.
References signifmatrix_.
: // eventVec_(0), signifmatrix_(2,2), set_worker_(0), xmet_(0), ymet_(0) { // std::cout << "in constructor ! " << std::endl; signifmatrix_(0,0)=signifmatrix_(1,0)=signifmatrix_(0,1)=signifmatrix_(1,1)=0; }
metsig::significanceAlgo::~significanceAlgo | ( | ) |
Definition at line 67 of file significanceAlgo.cc.
{ }
const void metsig::significanceAlgo::addObjects | ( | const std::vector< metsig::SigInputObj > & | EventVec | ) |
Definition at line 123 of file significanceAlgo.cc.
References funct::cos(), VarParsing::obj, and funct::sin().
Referenced by SignCaloSpecificAlgo::calculateBaseCaloMET(), and pat::PATMHTProducer::produce().
{ TMatrixD v_tot = signifmatrix_; //--- Loop over physics objects in the event ---// // for(unsigned int objnum=1; objnum < EventVec.size(); objnum++ ) { for(std::vector<SigInputObj>::const_iterator obj = eventVec.begin(); obj!= eventVec.end(); ++obj){ double et_tmp = obj->get_energy(); double phi_tmp = obj->get_phi(); double sigma_et = obj->get_sigma_e(); double sigma_tan = obj->get_sigma_tan(); double cosphi=cos(phi_tmp); double sinphi=sin(phi_tmp); double xval = et_tmp * cosphi; double yval = et_tmp * sinphi; xmet_ -= xval; ymet_ -= yval; set_worker_ += et_tmp; double sigma0_2=sigma_et*sigma_et; double sigma1_2=sigma_tan*sigma_tan; v_tot(0,0)+= sigma0_2*cosphi*cosphi + sigma1_2*sinphi*sinphi; v_tot(0,1)+= cosphi*sinphi*(sigma0_2 - sigma1_2); v_tot(1,0)+= cosphi*sinphi*(sigma0_2 - sigma1_2); v_tot(1,1)+= sigma1_2*cosphi*cosphi + sigma0_2*sinphi*sinphi; } signifmatrix_=v_tot; }
const void metsig::significanceAlgo::addSignifMatrix | ( | const TMatrixD & | input | ) |
Definition at line 43 of file significanceAlgo.cc.
References collect_tpl::input.
Referenced by SignCaloSpecificAlgo::calculateBaseCaloMET().
{ // check that the matrix is size 2: if(input.GetNrows()==2 && input.GetNcols()==2) { signifmatrix_+=input; } return; }
TMatrixD metsig::significanceAlgo::getSignifMatrix | ( | ) | const [inline] |
Definition at line 98 of file significanceAlgo.h.
References signifmatrix_.
Referenced by SignCaloSpecificAlgo::calculateBaseCaloMET(), and metsig::SignPFSpecificAlgo::getSignifMatrix().
{return signifmatrix_;}
void metsig::significanceAlgo::rotateMatrix | ( | Double_t | theta, |
TMatrixD & | v | ||
) | [private] |
Definition at line 73 of file significanceAlgo.cc.
References funct::cos(), csvReporter::r, and funct::sin().
const void metsig::significanceAlgo::setSignifMatrix | ( | const TMatrixD & | 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
Definition at line 54 of file significanceAlgo.cc.
References funct::cos(), collect_tpl::input, and funct::sin().
{ // check that the matrix is size 2: if(input.GetNrows()==2 && input.GetNcols()==2) { signifmatrix_=input; set_worker_=met_set; xmet_=met_r*cos(met_phi); ymet_=met_r*sin(met_phi); } return; }
const double metsig::significanceAlgo::significance | ( | double & | met_r, |
double & | met_phi, | ||
double & | met_set | ||
) |
Definition at line 156 of file significanceAlgo.cc.
References mathSSE::sqrt().
Referenced by SignCaloSpecificAlgo::calculateBaseCaloMET(), and pat::PATMHTProducer::produce().
{ if(signifmatrix_(0,0)==0 && signifmatrix_(1,1)==0 && signifmatrix_(1,0)==0 && signifmatrix_(0,1)==0){ //edm::LogWarning("SignCaloSpecificAlgo") << "Event Vector is empty! Return significance -1"; return(-1); } //--- Temporary variables ---// TMatrixD v_tot(2,2); TVectorD metvec(2); //--- Initialize sum of rotated covariance matrices ---// v_tot=signifmatrix_; // std::cout << "INPUT:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl; //--- Calculate magnitude and angle of MET, store in returned variables ---// met_r = sqrt(xmet_*xmet_ + ymet_*ymet_); met_set = set_worker_; met_phi= TMath::ATan2(ymet_, xmet_); // one other option: if particles cancel there could be small numbers. // this check fixes this, added by F.Blekman if(fabs(v_tot.Determinant())<0.000001) return -1; // save matrix into object: // std::cout << "SAVED:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl; v_tot.Invert(); // std::cout << "INVERTED:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl; metvec(0) = xmet_; metvec(1) = ymet_; double lnSignificance = metvec * (v_tot * metvec); // v_tot.Invert(); // std::cout << "INVERTED AGAIN:\n"<< v_tot(0,0) << "," << v_tot(0,1) << "\n" << v_tot(1,0) << "," << v_tot(1,1) << std::endl; return lnSignificance; }
const void metsig::significanceAlgo::subtractObjects | ( | const std::vector< metsig::SigInputObj > & | EventVec | ) |
Definition at line 89 of file significanceAlgo.cc.
References funct::cos(), VarParsing::obj, and funct::sin().
{ TMatrixD v_tot = signifmatrix_; //--- Loop over physics objects in the event ---// // for(unsigned int objnum=1; objnum < EventVec.size(); objnum++ ) { for(std::vector<SigInputObj>::const_iterator obj = eventVec.begin(); obj!= eventVec.end(); ++obj){ double et_tmp = obj->get_energy(); double phi_tmp = obj->get_phi(); double sigma_et = obj->get_sigma_e(); double sigma_tan = obj->get_sigma_tan(); double cosphi=cos(phi_tmp); double sinphi=sin(phi_tmp); double xval = et_tmp * cosphi; double yval = et_tmp * sinphi; xmet_ += xval; ymet_ += yval; set_worker_ -= et_tmp; double sigma0_2=sigma_et*sigma_et; double sigma1_2=sigma_tan*sigma_tan; v_tot(0,0)-= sigma0_2*cosphi*cosphi + sigma1_2*sinphi*sinphi; v_tot(0,1)-= cosphi*sinphi*(sigma0_2 - sigma1_2); v_tot(1,0)-= cosphi*sinphi*(sigma0_2 - sigma1_2); v_tot(1,1)-= sigma1_2*cosphi*cosphi + sigma0_2*sinphi*sinphi; } signifmatrix_=v_tot; }
double metsig::significanceAlgo::set_worker_ [private] |
Definition at line 106 of file significanceAlgo.h.
TMatrixD metsig::significanceAlgo::signifmatrix_ [private] |
Definition at line 104 of file significanceAlgo.h.
Referenced by getSignifMatrix(), and significanceAlgo().
double metsig::significanceAlgo::xmet_ [private] |
Definition at line 107 of file significanceAlgo.h.
double metsig::significanceAlgo::ymet_ [private] |
Definition at line 108 of file significanceAlgo.h.