CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/TopQuarkAnalysis/TopKinFitter/interface/CovarianceMatrix.h

Go to the documentation of this file.
00001 #ifndef CovarianceMatrix_h
00002 #define CovarianceMatrix_h
00003 
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 
00006 #include "DataFormats/PatCandidates/interface/PATObject.h"
00007 #include "DataFormats/PatCandidates/interface/Electron.h"
00008 #include "DataFormats/PatCandidates/interface/Muon.h"
00009 #include "DataFormats/PatCandidates/interface/Jet.h"
00010 #include "DataFormats/PatCandidates/interface/MET.h"
00011 
00012 #include "TopQuarkAnalysis/TopKinFitter/interface/TopKinFitter.h"
00013 #include "TopQuarkAnalysis/TopObjectResolutions/interface/MET.h"
00014 #include "TopQuarkAnalysis/TopObjectResolutions/interface/Jet.h"
00015 #include "TopQuarkAnalysis/TopObjectResolutions/interface/Muon.h"
00016 #include "TopQuarkAnalysis/TopObjectResolutions/interface/Electron.h"
00017 
00018 /*
00019   \class   CovarianceMatrix CovarianceMatrix.h "TopQuarkAnalysis/TopKinFitter/interface/CovarianceMatrix.h"
00020   
00021   \brief   Helper class used to setup covariance matrices for given objects and known resolutions
00022 
00023   More details to be added here...
00024   
00025 **/
00026 
00027 class CovarianceMatrix {
00028 
00029  public:
00030 
00031   enum ObjectType{kUdscJet, kBJet, kMuon, kElectron, kMet};
00032   
00034   CovarianceMatrix(){};
00036   CovarianceMatrix(const std::vector<edm::ParameterSet> udscResolutions, const std::vector<edm::ParameterSet> bResolutions,
00037                    const std::vector<double> jetEnergyResolutionScaleFactors, const std::vector<double> jetEnergyResolutionEtaBinning);
00039   CovarianceMatrix(const std::vector<edm::ParameterSet> udscResolutions, const std::vector<edm::ParameterSet> bResolutions,
00040                    const std::vector<edm::ParameterSet> lepResolutions, const std::vector<edm::ParameterSet> metResolutions,
00041                    const std::vector<double> jetEnergyResolutionScaleFactors, const std::vector<double> jetEnergyResolutionEtaBinning);
00042   // destructor
00043   ~CovarianceMatrix(){};
00044 
00046   template <class T>
00047     TMatrixD setupMatrix(const pat::PATObject<T>& object, const TopKinFitter::Param param, const std::string& resolutionProvider = "");
00049   TMatrixD setupMatrix(const TLorentzVector& object, const ObjectType objType, const TopKinFitter::Param param);
00051   double getResolution(const TLorentzVector& object, const ObjectType objType, const std::string& whichResolution = "");
00053   template <class T>
00054     double getResolution(const pat::PATObject<T>& object, const std::string& whichResolution, const bool isBJet=false) {
00055     return getResolution(TLorentzVector(object.px(), object.py(), object.pz(), object.energy()), getObjectType(object, isBJet), whichResolution); }
00056 
00057  private:
00058 
00060   std::vector<std::string> binsUdsc_, binsB_, binsLep_, binsMet_;
00062   std::vector<std::string> funcEtUdsc_ , funcEtB_ , funcEtLep_ , funcEtMet_;
00063   std::vector<std::string> funcEtaUdsc_, funcEtaB_, funcEtaLep_, funcEtaMet_;
00064   std::vector<std::string> funcPhiUdsc_, funcPhiB_, funcPhiLep_, funcPhiMet_;
00066   const std::vector<double> jetEnergyResolutionScaleFactors_;
00067   const std::vector<double> jetEnergyResolutionEtaBinning_;
00068 
00070   template <class T>
00071     ObjectType getObjectType(const pat::PATObject<T>& object, const bool isBJet=false);
00073   template <class T>
00074     double getEtaDependentScaleFactor(const pat::PATObject<T>& object);
00076   double getEtaDependentScaleFactor(const TLorentzVector& object);
00077 
00078 };
00079 
00080 template <class T>
00081 TMatrixD CovarianceMatrix::setupMatrix(const pat::PATObject<T>& object, const TopKinFitter::Param param, const std::string& resolutionProvider)
00082 {
00083   // This part is for pat objects with resolutions embedded
00084   if(object.hasKinResolution()) {
00085     TMatrixD CovM3 (3,3); CovM3.Zero();
00086     TMatrixD CovM4 (4,4); CovM4.Zero();
00087     TMatrixD* CovM = &CovM3;
00088     switch(param){
00089     case TopKinFitter::kEtEtaPhi :
00090       CovM3(0,0) = pow(object.resolEt(resolutionProvider) , 2);
00091       if( dynamic_cast<const reco::Jet*>(&object) )
00092         CovM3(0,0)*=getEtaDependentScaleFactor(object); 
00093       if( dynamic_cast<const reco::MET*>(&object) )
00094         CovM3(1,1) = pow(9999., 2);
00095       else
00096         CovM3(1,1) = pow(object.resolEta(resolutionProvider), 2);
00097       CovM3(2,2) = pow(object.resolPhi(resolutionProvider), 2);
00098       CovM = &CovM3;
00099       break;
00100     case TopKinFitter::kEtThetaPhi :
00101       CovM3(0,0) = pow(object.resolEt(resolutionProvider)   , 2);
00102       if( dynamic_cast<const reco::Jet*>(&object) )
00103         CovM3(0,0)*=getEtaDependentScaleFactor(object); 
00104       CovM3(1,1) = pow(object.resolTheta(resolutionProvider), 2);
00105       CovM3(2,2) = pow(object.resolPhi(resolutionProvider)  , 2);
00106       CovM = &CovM3;
00107       break;
00108     case TopKinFitter::kEMom :
00109       CovM4(0,0) = pow(1, 2);
00110       CovM4(1,1) = pow(1, 2);
00111       CovM4(2,2) = pow(1, 2);
00112       CovM4(3,3) = pow(1, 2);
00113       CovM = &CovM4;
00114       break;
00115     }
00116     return *CovM;
00117   }
00118   // This part is for objects without resolutions embedded
00119   else {
00120     const ObjectType objType = getObjectType(object, (resolutionProvider=="bjets"));
00121     const TLorentzVector p4(object.px(), object.py(), object.pz(), object.energy());
00122     return setupMatrix(p4, objType, param);
00123   }
00124 }
00125 
00126 template <class T>
00127 CovarianceMatrix::ObjectType CovarianceMatrix::getObjectType(const pat::PATObject<T>& object, const bool isBJet)
00128 {
00129   ObjectType objType;
00130   // jets
00131   if( dynamic_cast<const reco::Jet*>(&object) ) {
00132     if(isBJet)
00133       objType = kBJet;
00134     else
00135       objType = kUdscJet;
00136   }
00137   // muons
00138   else if( dynamic_cast<const reco::Muon*>(&object) )
00139     objType = kMuon;
00140   // electrons
00141   else if( dynamic_cast<const reco::GsfElectron*>(&object) ) 
00142     objType = kElectron;
00143   // MET
00144   else if( dynamic_cast<const reco::MET*>(&object) )
00145     objType = kMet;
00146   // catch anything else
00147   else
00148     throw cms::Exception("UnsupportedObject") << "The object given is not supported!\n";
00149   return objType;
00150 }
00151 
00152 template <class T>
00153 double CovarianceMatrix::getEtaDependentScaleFactor(const pat::PATObject<T>& object)
00154 {
00155   double etaDependentScaleFactor = 1.;
00156   for(unsigned int i=0; i<jetEnergyResolutionEtaBinning_.size(); i++){
00157     if(std::abs(object.eta())>=jetEnergyResolutionEtaBinning_[i] && jetEnergyResolutionEtaBinning_[i]>=0.){
00158       if(i==jetEnergyResolutionEtaBinning_.size()-1) {
00159         edm::LogWarning("CovarianceMatrix") << "object eta ("<<std::abs(object.eta())<<") beyond last eta bin ("<<jetEnergyResolutionEtaBinning_[i]<<") using scale factor 1.0!";
00160         etaDependentScaleFactor=1.;
00161         break;
00162       }
00163       etaDependentScaleFactor=jetEnergyResolutionScaleFactors_[i];
00164     }
00165     else
00166       break;
00167   }
00168   return etaDependentScaleFactor;
00169 }
00170 
00171 #endif