CMS 3D CMS Logo

Public Member Functions | Private Attributes

CovarianceMatrix Class Reference

#include <CovarianceMatrix.h>

List of all members.

Public Member Functions

 CovarianceMatrix ()
 CovarianceMatrix (const std::vector< edm::ParameterSet > udscResolutions, const std::vector< edm::ParameterSet > bResolutions)
 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)
template<class ObjectType >
double getResolution (const pat::PATObject< ObjectType > &object, const std::string whichResolution, bool isBJet)
template<class ObjectType >
TMatrixD setupMatrix (const pat::PATObject< ObjectType > &object, TopKinFitter::Param param, std::string resolutionProvider)
 ~CovarianceMatrix ()

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_

Detailed Description

Definition at line 21 of file CovarianceMatrix.h.


Constructor & Destructor Documentation

CovarianceMatrix::CovarianceMatrix ( ) [inline]

Definition at line 34 of file CovarianceMatrix.h.

{};
CovarianceMatrix::CovarianceMatrix ( const std::vector< edm::ParameterSet udscResolutions,
const std::vector< edm::ParameterSet bResolutions 
) [inline]

Definition at line 35 of file CovarianceMatrix.h.

References binsB_, binsUdsc_, Exception, funcEtaB_, funcEtaUdsc_, funcEtB_, funcEtUdsc_, funcPhiB_, and funcPhiUdsc_.

                                                                                                                   {
    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("WrongConfig") << "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("WrongConfig") << "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 
) [inline]

Definition at line 55 of file CovarianceMatrix.h.

References binsB_, binsLep_, binsMet_, binsUdsc_, Exception, funcEtaB_, funcEtaLep_, funcEtaMet_, funcEtaUdsc_, funcEtB_, funcEtLep_, funcEtMet_, funcEtUdsc_, funcPhiB_, funcPhiLep_, funcPhiMet_, and funcPhiUdsc_.

                                                                                                                                                                                                                     {
    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("WrongConfig") << "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("WrongConfig") << "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("WrongConfig") << "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("WrongConfig") << "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 93 of file CovarianceMatrix.h.

{};

Member Function Documentation

template<class ObjectType >
double CovarianceMatrix::getResolution ( const pat::PATObject< ObjectType > &  object,
const std::string  whichResolution = "",
bool  isBJet = false 
)

Definition at line 103 of file CovarianceMatrix.h.

References binsB_, binsLep_, binsMet_, binsUdsc_, relval_parameters_module::energy, Exception, funcEtaB_, funcEtaLep_, funcEtaMet_, funcEtaUdsc_, funcEtB_, funcEtLep_, funcEtMet_, funcEtUdsc_, funcPhiB_, funcPhiLep_, funcPhiMet_, funcPhiUdsc_, and i.

Referenced by setupMatrix().

{
  std::vector<std::string> * bins_, * funcEt_, * funcEta_, * funcPhi_;

  if( dynamic_cast<const reco::Jet*>(&object) && !isBJet ) {
    bins_    = &binsUdsc_;
    funcEt_  = &funcEtUdsc_;
    funcEta_ = &funcEtaUdsc_;
    funcPhi_ = &funcPhiUdsc_;
  }
  else if( dynamic_cast<const reco::Jet*>(&object) && isBJet ) {
    bins_    = &binsB_;
    funcEt_  = &funcEtB_;
    funcEta_ = &funcEtaB_;
    funcPhi_ = &funcPhiB_;
  }
  else if( dynamic_cast<const reco::Muon*>(&object) || dynamic_cast<const reco::GsfElectron*>(&object) ) {
    bins_    = &binsLep_;
    funcEt_  = &funcEtLep_;
    funcEta_ = &funcEtaLep_;
    funcPhi_ = &funcPhiLep_;
  }
  else if( dynamic_cast<const reco::MET*>(&object) ) {
    bins_    = &binsMet_;
    funcEt_  = &funcEtMet_;
    funcEta_ = &funcEtaMet_;
    funcPhi_ = &funcPhiMet_;
  }
  else{
    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;
    }
  }
  if(selectedBin>=0){
    if(whichResolution == "et")       return StringObjectFunction<reco::LeafCandidate>(funcEt_ ->at(selectedBin)).operator()(candidate);
    else if(whichResolution == "eta") return StringObjectFunction<reco::LeafCandidate>(funcEta_->at(selectedBin)).operator()(candidate);
    else if(whichResolution == "phi") return StringObjectFunction<reco::LeafCandidate>(funcPhi_->at(selectedBin)).operator()(candidate);
    else throw cms::Exception("ProgrammingError") << "Only 'et', 'eta' and 'phi' resolutions supported!\n";
  }
  return 0;
}
template<class ObjectType >
TMatrixD CovarianceMatrix::setupMatrix ( const pat::PATObject< ObjectType > &  object,
TopKinFitter::Param  param,
std::string  resolutionProvider = "" 
)

Definition at line 156 of file CovarianceMatrix.h.

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(), binsB_, 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::HelperMET::et(), res::HelperMuon::et(), res::HelperJet::eta(), eta(), res::HelperElectron::eta(), res::HelperMuon::eta(), getResolution(), res::HelperJet::kB, TopKinFitter::kEMom, TopKinFitter::kEtEtaPhi, TopKinFitter::kEtThetaPhi, res::HelperJet::kUds, res::HelperJet::phi(), res::HelperElectron::phi(), res::HelperMET::phi(), res::HelperMuon::phi(), funct::pow(), res::HelperJet::theta(), res::HelperMuon::theta(), and res::HelperElectron::theta().

Referenced by TtFullHadKinFitter::fit(), and TtSemiLepKinFitter::fit().

{
  TMatrixD CovM3 (3,3); CovM3.Zero();
  TMatrixD CovM4 (4,4); CovM4.Zero();
  TMatrixD* CovM = &CovM3;
  // This part is for pat objects with resolutions embedded
  if(object.hasKinResolution())
    {
      switch(param){
      case TopKinFitter::kEtEtaPhi :
        CovM3(0,0) = pow(object.resolEt(resolutionProvider) , 2);
        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);
        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;
      }
    }
  // This part is for objects without resolutions embedded
  else
    {
      double pt = object.pt(), eta = object.eta();
      // if object is a jet
      if( dynamic_cast<const reco::Jet*>(&object) ) {
        res::HelperJet jetRes;
        switch(param){
        case TopKinFitter::kEMom :
          if(resolutionProvider == "bjets") {
            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);
          }
          else {
            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);
          }
          CovM = &CovM4;
          break;
        case TopKinFitter::kEtEtaPhi : 
          if(resolutionProvider == "bjets") {
            if(!binsB_.size()){
              CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kB  ), 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, "et" , true), 2); 
              CovM3(1,1) = pow(getResolution(object, "eta", true), 2); 
              CovM3(2,2) = pow(getResolution(object, "phi", true), 2);
            }
          }
          else {
            if(!binsUdsc_.size()){
              CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kUds), 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, "et") , 2); 
              CovM3(1,1) = pow(getResolution(object, "eta"), 2); 
              CovM3(2,2) = pow(getResolution(object, "phi"), 2);
            }
          }       
          CovM = &CovM3;
          break;
        case TopKinFitter::kEtThetaPhi :
          if(resolutionProvider == "bjets") {
            CovM3(0,0) = pow(jetRes.et   (pt, eta, res::HelperJet::kB  ), 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);
          }
          else {
            CovM3(0,0) = pow(jetRes.et   (pt, eta, res::HelperJet::kUds), 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);
          }
          CovM = &CovM3;
          break;
        } 
      }
      // if object is an electron
      else if( dynamic_cast<const reco::GsfElectron*>(&object) ) {
        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);
          CovM = &CovM3;
          break;
        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, "et") , 2);
            CovM3(1,1) = pow(getResolution(object, "eta"), 2);
            CovM3(2,2) = pow(getResolution(object, "phi"), 2);
          }
          CovM = &CovM3;
          break;
        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);
          CovM = &CovM3;
          break;
        }
      }
      // if object is a muon
      else if( dynamic_cast<const reco::Muon*>(&object) ) {
        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);
          CovM = &CovM3;
          break;
        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, "et") , 2);
            CovM3(1,1) = pow(getResolution(object, "eta"), 2);
            CovM3(2,2) = pow(getResolution(object, "phi"), 2);
          }
          CovM = &CovM3;
          break;
        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);
          CovM = &CovM3;
          break;
        }
      }
      // if object is met
      else if( dynamic_cast<const reco::MET*>(&object) ) {
        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);
          CovM = &CovM3;
          break;
        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, "et") , 2);
            CovM3(1,1) = pow(getResolution(object, "eta"), 2);
            CovM3(2,2) = pow(getResolution(object, "phi"), 2);
          }
          CovM = &CovM3;
          break;
        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);
          CovM = &CovM3;
          break;
        }
      }
    }
  return *CovM;
}

Member Data Documentation

std::vector<std::string> CovarianceMatrix::binsB_ [private]

Definition at line 26 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), getResolution(), and setupMatrix().

std::vector<std::string> CovarianceMatrix::binsLep_ [private]

Definition at line 26 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), getResolution(), and setupMatrix().

std::vector<std::string> CovarianceMatrix::binsMet_ [private]

Definition at line 26 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 26 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), getResolution(), and setupMatrix().

std::vector<std::string> CovarianceMatrix::funcEtaB_ [private]

Definition at line 29 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcEtaLep_ [private]

Definition at line 29 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcEtaMet_ [private]

Definition at line 29 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcEtaUdsc_ [private]

Definition at line 29 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcEtB_ [private]

Definition at line 28 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcEtLep_ [private]

Definition at line 28 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcEtMet_ [private]

Definition at line 28 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcEtUdsc_ [private]

vectors for the resolution functions

Definition at line 28 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcPhiB_ [private]

Definition at line 30 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcPhiLep_ [private]

Definition at line 30 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcPhiMet_ [private]

Definition at line 30 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().

std::vector<std::string> CovarianceMatrix::funcPhiUdsc_ [private]

Definition at line 30 of file CovarianceMatrix.h.

Referenced by CovarianceMatrix(), and getResolution().