#include <CovarianceMatrix.h>
Public Types | |
enum | ObjectType { kUdscJet, kBJet, kMuon, kElectron, kMet } |
Public Member Functions | |
CovarianceMatrix () | |
default constructor | |
CovarianceMatrix (const std::vector< edm::ParameterSet > &udscResolutions, const std::vector< edm::ParameterSet > &bResolutions, const std::vector< double > &jetEnergyResolutionScaleFactors, const std::vector< double > &jetEnergyResolutionEtaBinning) | |
constructor for the fully-hadronic channel | |
CovarianceMatrix (const std::vector< edm::ParameterSet > &udscResolutions, const std::vector< edm::ParameterSet > &bResolutions, const std::vector< edm::ParameterSet > &lepResolutions, const std::vector< edm::ParameterSet > &metResolutions, const std::vector< double > &jetEnergyResolutionScaleFactors, const std::vector< double > &jetEnergyResolutionEtaBinning) | |
constructor for the lepton+jets channel | |
template<class T > | |
double | getResolution (const pat::PATObject< T > &object, const std::string &whichResolution, const bool isBJet=false) |
get resolution for a given PAT object | |
double | getResolution (const TLorentzVector &object, const ObjectType objType, const std::string &whichResolution="") |
get resolution for a given component of an object | |
TMatrixD | setupMatrix (const TLorentzVector &object, const ObjectType objType, const TopKinFitter::Param param) |
return covariance matrix for a plain 4-vector | |
template<class T > | |
TMatrixD | setupMatrix (const pat::PATObject< T > &object, const TopKinFitter::Param param, const std::string &resolutionProvider="") |
return covariance matrix for a PAT object | |
~CovarianceMatrix () | |
Private Member Functions | |
template<class T > | |
double | getEtaDependentScaleFactor (const pat::PATObject< T > &object) |
get eta dependent smear factor for a PAT object | |
double | getEtaDependentScaleFactor (const TLorentzVector &object) |
get eta-dependent scale factor for a plain 4-vector | |
template<class T > | |
ObjectType | getObjectType (const pat::PATObject< T > &object, const bool isBJet=false) |
determine type for a given PAT object | |
Private Attributes | |
std::vector< std::string > | binsB_ |
std::vector< std::string > | binsLep_ |
std::vector< std::string > | binsMet_ |
std::vector< std::string > | binsUdsc_ |
vector of strings for the binning of the resolutions | |
std::vector< std::string > | funcEtaB_ |
std::vector< std::string > | funcEtaLep_ |
std::vector< std::string > | funcEtaMet_ |
std::vector< std::string > | funcEtaUdsc_ |
std::vector< std::string > | funcEtB_ |
std::vector< std::string > | funcEtLep_ |
std::vector< std::string > | funcEtMet_ |
std::vector< std::string > | funcEtUdsc_ |
vectors for the resolution functions | |
std::vector< std::string > | funcPhiB_ |
std::vector< std::string > | funcPhiLep_ |
std::vector< std::string > | funcPhiMet_ |
std::vector< std::string > | funcPhiUdsc_ |
const std::vector< double > | jetEnergyResolutionEtaBinning_ |
const std::vector< double > | jetEnergyResolutionScaleFactors_ |
scale factors for the jet energy resolution |
Definition at line 27 of file CovarianceMatrix.h.
CovarianceMatrix::CovarianceMatrix | ( | ) | [inline] |
CovarianceMatrix::CovarianceMatrix | ( | const std::vector< edm::ParameterSet > & | udscResolutions, |
const std::vector< edm::ParameterSet > & | bResolutions, | ||
const std::vector< double > & | jetEnergyResolutionScaleFactors, | ||
const std::vector< double > & | jetEnergyResolutionEtaBinning | ||
) |
constructor for the fully-hadronic channel
Definition at line 6 of file CovarianceMatrix.cc.
References binsB_, binsUdsc_, Exception, funcEtaB_, funcEtaUdsc_, funcEtB_, funcEtUdsc_, funcPhiB_, funcPhiUdsc_, and AlCaHLTBitMon_QueryRunRegistry::string.
: jetEnergyResolutionScaleFactors_(jetEnergyResolutionScaleFactors), jetEnergyResolutionEtaBinning_(jetEnergyResolutionEtaBinning) { for(std::vector<edm::ParameterSet>::const_iterator iSet = udscResolutions.begin(); iSet != udscResolutions.end(); ++iSet){ if(iSet->exists("bin")) binsUdsc_.push_back(iSet->getParameter<std::string>("bin")); else if(udscResolutions.size()==1) binsUdsc_.push_back(""); else throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n"; funcEtUdsc_.push_back(iSet->getParameter<std::string>("et")); funcEtaUdsc_.push_back(iSet->getParameter<std::string>("eta")); funcPhiUdsc_.push_back(iSet->getParameter<std::string>("phi")); } for(std::vector<edm::ParameterSet>::const_iterator iSet = bResolutions.begin(); iSet != bResolutions.end(); ++iSet){ if(iSet->exists("bin")) binsB_.push_back(iSet->getParameter<std::string>("bin")); else if(bResolutions.size()==1) binsB_.push_back(""); else throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n"; funcEtB_.push_back(iSet->getParameter<std::string>("et")); funcEtaB_.push_back(iSet->getParameter<std::string>("eta")); funcPhiB_.push_back(iSet->getParameter<std::string>("phi")); } }
CovarianceMatrix::CovarianceMatrix | ( | const std::vector< edm::ParameterSet > & | udscResolutions, |
const std::vector< edm::ParameterSet > & | bResolutions, | ||
const std::vector< edm::ParameterSet > & | lepResolutions, | ||
const std::vector< edm::ParameterSet > & | metResolutions, | ||
const std::vector< double > & | jetEnergyResolutionScaleFactors, | ||
const std::vector< double > & | jetEnergyResolutionEtaBinning | ||
) |
constructor for the lepton+jets channel
Definition at line 30 of file CovarianceMatrix.cc.
References binsB_, binsLep_, binsMet_, binsUdsc_, Exception, funcEtaB_, funcEtaLep_, funcEtaMet_, funcEtaUdsc_, funcEtB_, funcEtLep_, funcEtMet_, funcEtUdsc_, funcPhiB_, funcPhiLep_, funcPhiMet_, funcPhiUdsc_, i, jetEnergyResolutionEtaBinning_, jetEnergyResolutionScaleFactors_, and AlCaHLTBitMon_QueryRunRegistry::string.
: jetEnergyResolutionScaleFactors_(jetEnergyResolutionScaleFactors), jetEnergyResolutionEtaBinning_(jetEnergyResolutionEtaBinning) { if(jetEnergyResolutionScaleFactors_.size()+1!=jetEnergyResolutionEtaBinning_.size()) throw cms::Exception("Configuration") << "The number of scale factors does not fit to the number of eta bins!\n"; for(unsigned int i=0; i<jetEnergyResolutionEtaBinning_.size(); i++) if(jetEnergyResolutionEtaBinning_[i]<0. && i<jetEnergyResolutionEtaBinning_.size()-1) throw cms::Exception("Configuration") << "eta binning in absolut values required!\n"; for(std::vector<edm::ParameterSet>::const_iterator iSet = udscResolutions.begin(); iSet != udscResolutions.end(); ++iSet){ if(iSet->exists("bin")) binsUdsc_.push_back(iSet->getParameter<std::string>("bin")); else if(udscResolutions.size()==1) binsUdsc_.push_back(""); else throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n"; funcEtUdsc_.push_back(iSet->getParameter<std::string>("et")); funcEtaUdsc_.push_back(iSet->getParameter<std::string>("eta")); funcPhiUdsc_.push_back(iSet->getParameter<std::string>("phi")); } for(std::vector<edm::ParameterSet>::const_iterator iSet = bResolutions.begin(); iSet != bResolutions.end(); ++iSet){ if(iSet->exists("bin")) binsB_.push_back(iSet->getParameter<std::string>("bin")); else if(bResolutions.size()==1) binsB_.push_back(""); else throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n"; funcEtB_.push_back(iSet->getParameter<std::string>("et")); funcEtaB_.push_back(iSet->getParameter<std::string>("eta")); funcPhiB_.push_back(iSet->getParameter<std::string>("phi")); } for(std::vector<edm::ParameterSet>::const_iterator iSet = lepResolutions.begin(); iSet != lepResolutions.end(); ++iSet){ if(iSet->exists("bin")) binsLep_.push_back(iSet->getParameter<std::string>("bin")); else if(lepResolutions.size()==1) binsLep_.push_back(""); else throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n"; funcEtLep_.push_back(iSet->getParameter<std::string>("et")); funcEtaLep_.push_back(iSet->getParameter<std::string>("eta")); funcPhiLep_.push_back(iSet->getParameter<std::string>("phi")); } for(std::vector<edm::ParameterSet>::const_iterator iSet = metResolutions.begin(); iSet != metResolutions.end(); ++iSet){ if(iSet->exists("bin")) binsMet_.push_back(iSet->getParameter<std::string>("bin")); else if(metResolutions.size()==1) binsMet_.push_back(""); else throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n"; funcEtMet_.push_back(iSet->getParameter<std::string>("et")); funcEtaMet_.push_back(iSet->getParameter<std::string>("eta")); funcPhiMet_.push_back(iSet->getParameter<std::string>("phi")); } }
CovarianceMatrix::~CovarianceMatrix | ( | ) | [inline] |
Definition at line 43 of file CovarianceMatrix.h.
{};
double CovarianceMatrix::getEtaDependentScaleFactor | ( | const pat::PATObject< T > & | object | ) | [private] |
get eta dependent smear factor for a PAT object
Definition at line 153 of file CovarianceMatrix.h.
References abs, eta(), i, jetEnergyResolutionEtaBinning_, and jetEnergyResolutionScaleFactors_.
Referenced by setupMatrix().
{ double etaDependentScaleFactor = 1.; for(unsigned int i=0; i<jetEnergyResolutionEtaBinning_.size(); i++){ if(std::abs(object.eta())>=jetEnergyResolutionEtaBinning_[i] && jetEnergyResolutionEtaBinning_[i]>=0.){ if(i==jetEnergyResolutionEtaBinning_.size()-1) { edm::LogWarning("CovarianceMatrix") << "object eta ("<<std::abs(object.eta())<<") beyond last eta bin ("<<jetEnergyResolutionEtaBinning_[i]<<") using scale factor 1.0!"; etaDependentScaleFactor=1.; break; } etaDependentScaleFactor=jetEnergyResolutionScaleFactors_[i]; } else break; } return etaDependentScaleFactor; }
double CovarianceMatrix::getEtaDependentScaleFactor | ( | const TLorentzVector & | object | ) | [private] |
get eta-dependent scale factor for a plain 4-vector
Definition at line 308 of file CovarianceMatrix.cc.
References abs, reco::tau::disc::Eta(), i, jetEnergyResolutionEtaBinning_, and jetEnergyResolutionScaleFactors_.
{ double etaDependentScaleFactor = 1.; for(unsigned int i=0; i<jetEnergyResolutionEtaBinning_.size(); i++){ if(std::abs(object.Eta())>=jetEnergyResolutionEtaBinning_[i] && jetEnergyResolutionEtaBinning_[i]>=0.){ if(i==jetEnergyResolutionEtaBinning_.size()-1) { edm::LogWarning("CovarianceMatrix") << "object eta ("<<std::abs(object.Eta())<<") beyond last eta bin ("<<jetEnergyResolutionEtaBinning_[i]<<") using scale factor 1.0!"; etaDependentScaleFactor=1.; break; } etaDependentScaleFactor=jetEnergyResolutionScaleFactors_[i]; } else break; } return etaDependentScaleFactor; }
CovarianceMatrix::ObjectType CovarianceMatrix::getObjectType | ( | const pat::PATObject< T > & | object, |
const bool | isBJet = false |
||
) | [private] |
determine type for a given PAT object
Definition at line 127 of file CovarianceMatrix.h.
References Exception, kBJet, kElectron, kMet, kMuon, and kUdscJet.
Referenced by getResolution(), and setupMatrix().
{ ObjectType objType; // jets if( dynamic_cast<const reco::Jet*>(&object) ) { if(isBJet) objType = kBJet; else objType = kUdscJet; } // muons else if( dynamic_cast<const reco::Muon*>(&object) ) objType = kMuon; // electrons else if( dynamic_cast<const reco::GsfElectron*>(&object) ) objType = kElectron; // MET else if( dynamic_cast<const reco::MET*>(&object) ) objType = kMet; // catch anything else else throw cms::Exception("UnsupportedObject") << "The object given is not supported!\n"; return objType; }
double CovarianceMatrix::getResolution | ( | const pat::PATObject< T > & | object, |
const std::string & | whichResolution, | ||
const bool | isBJet = false |
||
) | [inline] |
get resolution for a given PAT object
Definition at line 54 of file CovarianceMatrix.h.
References relval_parameters_module::energy, getObjectType(), and getResolution().
{ return getResolution(TLorentzVector(object.px(), object.py(), object.pz(), object.energy()), getObjectType(object, isBJet), whichResolution); }
double CovarianceMatrix::getResolution | ( | const TLorentzVector & | object, |
const ObjectType | objType, | ||
const std::string & | whichResolution = "" |
||
) |
get resolution for a given component of an object
Definition at line 79 of file CovarianceMatrix.cc.
References binsB_, binsLep_, binsMet_, binsUdsc_, Exception, funcEtaB_, funcEtaLep_, funcEtaMet_, funcEtaUdsc_, funcEtB_, funcEtLep_, funcEtMet_, funcEtUdsc_, funcPhiB_, funcPhiLep_, funcPhiMet_, funcPhiUdsc_, i, kBJet, kElectron, kMet, kMuon, and kUdscJet.
Referenced by getResolution(), and setupMatrix().
{ std::vector<std::string> * bins_, * funcEt_, * funcEta_, * funcPhi_; switch(objType) { case kUdscJet: bins_ = &binsUdsc_; funcEt_ = &funcEtUdsc_; funcEta_ = &funcEtaUdsc_; funcPhi_ = &funcPhiUdsc_; break; case kBJet: bins_ = &binsB_; funcEt_ = &funcEtB_; funcEta_ = &funcEtaB_; funcPhi_ = &funcPhiB_; break; case kMuon: bins_ = &binsLep_; funcEt_ = &funcEtLep_; funcEta_ = &funcEtaLep_; funcPhi_ = &funcPhiLep_; break; case kElectron: bins_ = &binsLep_; funcEt_ = &funcEtLep_; funcEta_ = &funcEtaLep_; funcPhi_ = &funcPhiLep_; break; case kMet: bins_ = &binsMet_; funcEt_ = &funcEtMet_; funcEta_ = &funcEtaMet_; funcPhi_ = &funcPhiMet_; break; default: throw cms::Exception("UnsupportedObject") << "The object given is not supported!\n"; } int selectedBin=-1; reco::LeafCandidate candidate; for(unsigned int i=0; i<bins_->size(); ++i){ StringCutObjectSelector<reco::LeafCandidate> select_(bins_->at(i)); candidate = reco::LeafCandidate( 0, reco::LeafCandidate::LorentzVector(object.Px(), object.Py(), object.Pz(), object.Energy()) ); if(select_(candidate)){ selectedBin = i; break; } } double res = 0; if(selectedBin>=0){ if(whichResolution == "et") res = StringObjectFunction<reco::LeafCandidate>(funcEt_ ->at(selectedBin)).operator()(candidate); else if(whichResolution == "eta") res = StringObjectFunction<reco::LeafCandidate>(funcEta_->at(selectedBin)).operator()(candidate); else if(whichResolution == "phi") res = StringObjectFunction<reco::LeafCandidate>(funcPhi_->at(selectedBin)).operator()(candidate); else throw cms::Exception("ProgrammingError") << "Only 'et', 'eta' and 'phi' resolutions supported!\n"; } return res; }
TMatrixD CovarianceMatrix::setupMatrix | ( | const pat::PATObject< T > & | object, |
const TopKinFitter::Param | param, | ||
const std::string & | resolutionProvider = "" |
||
) |
return covariance matrix for a PAT object
Definition at line 81 of file CovarianceMatrix.h.
References relval_parameters_module::energy, getEtaDependentScaleFactor(), getObjectType(), TopKinFitter::kEMom, TopKinFitter::kEtEtaPhi, TopKinFitter::kEtThetaPhi, p4, and funct::pow().
Referenced by TtFullHadKinFitter::fit(), and TtSemiLepKinFitter::fit().
{ // This part is for pat objects with resolutions embedded if(object.hasKinResolution()) { TMatrixD CovM3 (3,3); CovM3.Zero(); TMatrixD CovM4 (4,4); CovM4.Zero(); TMatrixD* CovM = &CovM3; switch(param){ case TopKinFitter::kEtEtaPhi : CovM3(0,0) = pow(object.resolEt(resolutionProvider) , 2); if( dynamic_cast<const reco::Jet*>(&object) ) CovM3(0,0)*=getEtaDependentScaleFactor(object); if( dynamic_cast<const reco::MET*>(&object) ) CovM3(1,1) = pow(9999., 2); else CovM3(1,1) = pow(object.resolEta(resolutionProvider), 2); CovM3(2,2) = pow(object.resolPhi(resolutionProvider), 2); CovM = &CovM3; break; case TopKinFitter::kEtThetaPhi : CovM3(0,0) = pow(object.resolEt(resolutionProvider) , 2); if( dynamic_cast<const reco::Jet*>(&object) ) CovM3(0,0)*=getEtaDependentScaleFactor(object); CovM3(1,1) = pow(object.resolTheta(resolutionProvider), 2); CovM3(2,2) = pow(object.resolPhi(resolutionProvider) , 2); CovM = &CovM3; break; case TopKinFitter::kEMom : CovM4(0,0) = pow(1, 2); CovM4(1,1) = pow(1, 2); CovM4(2,2) = pow(1, 2); CovM4(3,3) = pow(1, 2); CovM = &CovM4; break; } return *CovM; } // This part is for objects without resolutions embedded else { const ObjectType objType = getObjectType(object, (resolutionProvider=="bjets")); const TLorentzVector p4(object.px(), object.py(), object.pz(), object.energy()); return setupMatrix(p4, objType, param); } }
TMatrixD CovarianceMatrix::setupMatrix | ( | const TLorentzVector & | object, |
const ObjectType | objType, | ||
const TopKinFitter::Param | param | ||
) |
return covariance matrix for a plain 4-vector
Definition at line 138 of file CovarianceMatrix.cc.
References res::HelperJet::a(), res::HelperMuon::a(), res::HelperElectron::a(), res::HelperMET::a(), res::HelperJet::b(), res::HelperMuon::b(), res::HelperElectron::b(), res::HelperMET::b(), binsLep_, binsMet_, binsUdsc_, res::HelperJet::c(), res::HelperMuon::c(), res::HelperElectron::c(), res::HelperMET::c(), res::HelperJet::d(), res::HelperElectron::et(), res::HelperJet::et(), res::HelperMuon::et(), res::HelperMET::et(), res::HelperJet::eta(), eta(), res::HelperElectron::eta(), res::HelperMuon::eta(), Exception, getEtaDependentScaleFactor(), getResolution(), res::HelperJet::kB, kBJet, kElectron, TopKinFitter::kEMom, TopKinFitter::kEtEtaPhi, TopKinFitter::kEtThetaPhi, kMet, kMuon, res::HelperJet::kUds, kUdscJet, res::HelperJet::phi(), res::HelperElectron::phi(), res::HelperMET::phi(), res::HelperMuon::phi(), funct::pow(), res::HelperJet::theta(), res::HelperMuon::theta(), and res::HelperElectron::theta().
{ TMatrixD CovM3 (3,3); CovM3.Zero(); TMatrixD CovM4 (4,4); CovM4.Zero(); const double pt = object.Pt(); const double eta = object.Eta(); switch(objType) { case kUdscJet: // if object is a non-b jet { res::HelperJet jetRes; switch(param) { case TopKinFitter::kEMom : CovM4(0,0) = pow(jetRes.a (pt, eta, res::HelperJet::kUds), 2); CovM4(1,1) = pow(jetRes.b (pt, eta, res::HelperJet::kUds), 2); CovM4(2,2) = pow(jetRes.c (pt, eta, res::HelperJet::kUds), 2); CovM4(3,3) = pow(jetRes.d (pt, eta, res::HelperJet::kUds), 2); return CovM4; case TopKinFitter::kEtEtaPhi : if(!binsUdsc_.size()){ CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kUds), 2); CovM3(0,0)*= pow(getEtaDependentScaleFactor(object) , 2); CovM3(1,1) = pow(jetRes.eta(pt, eta, res::HelperJet::kUds), 2); CovM3(2,2) = pow(jetRes.phi(pt, eta, res::HelperJet::kUds), 2); } else{ CovM3(0,0) = pow(getResolution(object, objType, "et") , 2); CovM3(0,0)*= pow(getEtaDependentScaleFactor(object) , 2); CovM3(1,1) = pow(getResolution(object, objType, "eta"), 2); CovM3(2,2) = pow(getResolution(object, objType, "phi"), 2); } return CovM3; case TopKinFitter::kEtThetaPhi : CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kUds), 2); CovM3(0,0)*= pow(getEtaDependentScaleFactor(object) , 2); CovM3(1,1) = pow(jetRes.theta(pt, eta, res::HelperJet::kUds), 2); CovM3(2,2) = pow(jetRes.phi (pt, eta, res::HelperJet::kUds), 2); return CovM3; } } break; case kBJet: // if object is a b jet { res::HelperJet jetRes; switch(param) { case TopKinFitter::kEMom : CovM4(0,0) = pow(jetRes.a (pt, eta, res::HelperJet::kB), 2); CovM4(1,1) = pow(jetRes.b (pt, eta, res::HelperJet::kB), 2); CovM4(2,2) = pow(jetRes.c (pt, eta, res::HelperJet::kB), 2); CovM4(3,3) = pow(jetRes.d (pt, eta, res::HelperJet::kB), 2); return CovM4; case TopKinFitter::kEtEtaPhi : if(!binsUdsc_.size()){ CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kB), 2); CovM3(0,0)*= pow(getEtaDependentScaleFactor(object) , 2); CovM3(1,1) = pow(jetRes.eta(pt, eta, res::HelperJet::kB), 2); CovM3(2,2) = pow(jetRes.phi(pt, eta, res::HelperJet::kB), 2); } else{ CovM3(0,0) = pow(getResolution(object, objType, "et") , 2); CovM3(0,0)*= pow(getEtaDependentScaleFactor(object) , 2); CovM3(1,1) = pow(getResolution(object, objType, "eta"), 2); CovM3(2,2) = pow(getResolution(object, objType, "phi"), 2); } return CovM3; case TopKinFitter::kEtThetaPhi : CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kB), 2); CovM3(0,0)*= pow(getEtaDependentScaleFactor(object) , 2); CovM3(1,1) = pow(jetRes.theta(pt, eta, res::HelperJet::kB), 2); CovM3(2,2) = pow(jetRes.phi (pt, eta, res::HelperJet::kB), 2); return CovM3; } } break; case kMuon: // if object is a muon { res::HelperMuon muonRes; switch(param){ case TopKinFitter::kEMom : CovM3(0,0) = pow(muonRes.a (pt, eta), 2); CovM3(1,1) = pow(muonRes.b (pt, eta), 2); CovM3(2,2) = pow(muonRes.c (pt, eta), 2); return CovM3; case TopKinFitter::kEtEtaPhi : if(!binsLep_.size()){ CovM3(0,0) = pow(muonRes.et (pt, eta), 2); CovM3(1,1) = pow(muonRes.eta(pt, eta), 2); CovM3(2,2) = pow(muonRes.phi(pt, eta), 2); } else{ CovM3(0,0) = pow(getResolution(object, objType, "et") , 2); CovM3(1,1) = pow(getResolution(object, objType, "eta"), 2); CovM3(2,2) = pow(getResolution(object, objType, "phi"), 2); } return CovM3; case TopKinFitter::kEtThetaPhi : CovM3(0,0) = pow(muonRes.et (pt, eta), 2); CovM3(1,1) = pow(muonRes.theta(pt, eta), 2); CovM3(2,2) = pow(muonRes.phi (pt, eta), 2); return CovM3; } } break; case kElectron: { // if object is an electron res::HelperElectron elecRes; switch(param){ case TopKinFitter::kEMom : CovM3(0,0) = pow(elecRes.a (pt, eta), 2); CovM3(1,1) = pow(elecRes.b (pt, eta), 2); CovM3(2,2) = pow(elecRes.c (pt, eta), 2); return CovM3; case TopKinFitter::kEtEtaPhi : if(!binsLep_.size()){ CovM3(0,0) = pow(elecRes.et (pt, eta), 2); CovM3(1,1) = pow(elecRes.eta(pt, eta), 2); CovM3(2,2) = pow(elecRes.phi(pt, eta), 2); } else{ CovM3(0,0) = pow(getResolution(object, objType, "et") , 2); CovM3(1,1) = pow(getResolution(object, objType, "eta"), 2); CovM3(2,2) = pow(getResolution(object, objType, "phi"), 2); } return CovM3; case TopKinFitter::kEtThetaPhi : CovM3(0,0) = pow(elecRes.et (pt, eta), 2); CovM3(1,1) = pow(elecRes.theta(pt, eta), 2); CovM3(2,2) = pow(elecRes.phi (pt, eta), 2); return CovM3; } } break; case kMet: // if object is met { res::HelperMET metRes; switch(param){ case TopKinFitter::kEMom : CovM3(0,0) = pow(metRes.a(pt), 2); CovM3(1,1) = pow(metRes.b(pt), 2); CovM3(2,2) = pow(metRes.c(pt), 2); return CovM3; case TopKinFitter::kEtEtaPhi : if(!binsMet_.size()){ CovM3(0,0) = pow(metRes.et(pt) , 2); CovM3(1,1) = pow( 9999. , 2); CovM3(2,2) = pow(metRes.phi(pt), 2); } else{ CovM3(0,0) = pow(getResolution(object, objType, "et") , 2); CovM3(1,1) = pow(getResolution(object, objType, "eta"), 2); CovM3(2,2) = pow(getResolution(object, objType, "phi"), 2); } return CovM3; case TopKinFitter::kEtThetaPhi : CovM3(0,0) = pow(metRes.et(pt) , 2); CovM3(1,1) = pow( 9999. , 2); CovM3(2,2) = pow(metRes.phi(pt), 2); return CovM3; } } break; } cms::Exception("Logic") << "Something went wrong while trying to setup a covariance matrix!\n"; return CovM4; //should never get here }
std::vector<std::string> CovarianceMatrix::binsB_ [private] |
Definition at line 60 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::binsLep_ [private] |
Definition at line 60 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), getResolution(), and setupMatrix().
std::vector<std::string> CovarianceMatrix::binsMet_ [private] |
Definition at line 60 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), getResolution(), and setupMatrix().
std::vector<std::string> CovarianceMatrix::binsUdsc_ [private] |
vector of strings for the binning of the resolutions
Definition at line 60 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), getResolution(), and setupMatrix().
std::vector<std::string> CovarianceMatrix::funcEtaB_ [private] |
Definition at line 63 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcEtaLep_ [private] |
Definition at line 63 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcEtaMet_ [private] |
Definition at line 63 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcEtaUdsc_ [private] |
Definition at line 63 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcEtB_ [private] |
Definition at line 62 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcEtLep_ [private] |
Definition at line 62 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcEtMet_ [private] |
Definition at line 62 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcEtUdsc_ [private] |
vectors for the resolution functions
Definition at line 62 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcPhiB_ [private] |
Definition at line 64 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcPhiLep_ [private] |
Definition at line 64 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcPhiMet_ [private] |
Definition at line 64 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
std::vector<std::string> CovarianceMatrix::funcPhiUdsc_ [private] |
Definition at line 64 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getResolution().
const std::vector<double> CovarianceMatrix::jetEnergyResolutionEtaBinning_ [private] |
Definition at line 67 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getEtaDependentScaleFactor().
const std::vector<double> CovarianceMatrix::jetEnergyResolutionScaleFactors_ [private] |
scale factors for the jet energy resolution
Definition at line 66 of file CovarianceMatrix.h.
Referenced by CovarianceMatrix(), and getEtaDependentScaleFactor().