CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/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 "CommonTools/Utils/interface/StringObjectFunction.h"
00007 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
00008 
00009 #include "DataFormats/PatCandidates/interface/PATObject.h"
00010 #include "DataFormats/PatCandidates/interface/Electron.h"
00011 #include "DataFormats/PatCandidates/interface/Muon.h"
00012 #include "DataFormats/PatCandidates/interface/Jet.h"
00013 #include "DataFormats/PatCandidates/interface/MET.h"
00014 
00015 #include "TopQuarkAnalysis/TopObjectResolutions/interface/MET.h"
00016 #include "TopQuarkAnalysis/TopObjectResolutions/interface/Jet.h"
00017 #include "TopQuarkAnalysis/TopObjectResolutions/interface/Muon.h"
00018 #include "TopQuarkAnalysis/TopObjectResolutions/interface/Electron.h"
00019 #include "TopQuarkAnalysis/TopKinFitter/interface/TopKinFitter.h"
00020 
00021 class CovarianceMatrix {
00022     
00023  private:
00024 
00026   std::vector<std::string> binsUdsc_, binsB_, binsLep_, binsMet_;
00028   std::vector<std::string> funcEtUdsc_ , funcEtB_ , funcEtLep_ , funcEtMet_;
00029   std::vector<std::string> funcEtaUdsc_, funcEtaB_, funcEtaLep_, funcEtaMet_;
00030   std::vector<std::string> funcPhiUdsc_, funcPhiB_, funcPhiLep_, funcPhiMet_;
00031 
00032  public:
00033   
00034   CovarianceMatrix(){};
00035   CovarianceMatrix(const std::vector<edm::ParameterSet> udscResolutions, const std::vector<edm::ParameterSet> bResolutions){
00036     for(std::vector<edm::ParameterSet>::const_iterator iSet = udscResolutions.begin(); iSet != udscResolutions.end(); ++iSet){
00037       if(iSet->exists("bin")) binsUdsc_.push_back(iSet->getParameter<std::string>("bin"));
00038       else if(udscResolutions.size()==1) binsUdsc_.push_back("");
00039       else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
00040       
00041       funcEtUdsc_.push_back(iSet->getParameter<std::string>("et"));
00042       funcEtaUdsc_.push_back(iSet->getParameter<std::string>("eta"));
00043       funcPhiUdsc_.push_back(iSet->getParameter<std::string>("phi"));
00044     }
00045     for(std::vector<edm::ParameterSet>::const_iterator iSet = bResolutions.begin(); iSet != bResolutions.end(); ++iSet){
00046       if(iSet->exists("bin")) binsB_.push_back(iSet->getParameter<std::string>("bin"));
00047       else if(bResolutions.size()==1) binsB_.push_back("");
00048       else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
00049     
00050       funcEtB_.push_back(iSet->getParameter<std::string>("et"));
00051       funcEtaB_.push_back(iSet->getParameter<std::string>("eta"));
00052       funcPhiB_.push_back(iSet->getParameter<std::string>("phi"));
00053     }
00054   };
00055   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){
00056     for(std::vector<edm::ParameterSet>::const_iterator iSet = udscResolutions.begin(); iSet != udscResolutions.end(); ++iSet){
00057       if(iSet->exists("bin")) binsUdsc_.push_back(iSet->getParameter<std::string>("bin"));
00058       else if(udscResolutions.size()==1) binsUdsc_.push_back("");
00059       else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
00060 
00061       funcEtUdsc_.push_back(iSet->getParameter<std::string>("et"));
00062       funcEtaUdsc_.push_back(iSet->getParameter<std::string>("eta"));
00063       funcPhiUdsc_.push_back(iSet->getParameter<std::string>("phi"));
00064     }
00065     for(std::vector<edm::ParameterSet>::const_iterator iSet = bResolutions.begin(); iSet != bResolutions.end(); ++iSet){
00066       if(iSet->exists("bin")) binsB_.push_back(iSet->getParameter<std::string>("bin"));
00067       else if(bResolutions.size()==1) binsB_.push_back("");
00068       else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
00069 
00070       funcEtB_.push_back(iSet->getParameter<std::string>("et"));
00071       funcEtaB_.push_back(iSet->getParameter<std::string>("eta"));
00072       funcPhiB_.push_back(iSet->getParameter<std::string>("phi"));
00073     }
00074     for(std::vector<edm::ParameterSet>::const_iterator iSet = lepResolutions.begin(); iSet != lepResolutions.end(); ++iSet){
00075       if(iSet->exists("bin")) binsLep_.push_back(iSet->getParameter<std::string>("bin"));
00076       else if(lepResolutions.size()==1) binsLep_.push_back("");
00077       else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
00078 
00079       funcEtLep_.push_back(iSet->getParameter<std::string>("et"));
00080       funcEtaLep_.push_back(iSet->getParameter<std::string>("eta"));
00081       funcPhiLep_.push_back(iSet->getParameter<std::string>("phi"));
00082     }
00083     for(std::vector<edm::ParameterSet>::const_iterator iSet = metResolutions.begin(); iSet != metResolutions.end(); ++iSet){
00084       if(iSet->exists("bin")) binsMet_.push_back(iSet->getParameter<std::string>("bin"));
00085       else if(metResolutions.size()==1) binsMet_.push_back("");
00086       else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
00087 
00088       funcEtMet_.push_back(iSet->getParameter<std::string>("et"));
00089       funcEtaMet_.push_back(iSet->getParameter<std::string>("eta"));
00090       funcPhiMet_.push_back(iSet->getParameter<std::string>("phi"));
00091     }
00092    };
00093    ~CovarianceMatrix(){};
00094  
00095   template <class ObjectType>
00096     TMatrixD setupMatrix(const pat::PATObject<ObjectType>& object, TopKinFitter::Param param, std::string resolutionProvider);
00097 
00098   template <class ObjectType>
00099     double getResolution(const pat::PATObject<ObjectType>& object, const std::string whichResolution, bool isBJet);
00100 };
00101 
00102 template <class ObjectType>
00103 double CovarianceMatrix::getResolution(const pat::PATObject<ObjectType>& object, const std::string whichResolution = "", bool isBJet = false)
00104 {
00105   std::vector<std::string> * bins_, * funcEt_, * funcEta_, * funcPhi_;
00106 
00107   if( dynamic_cast<const reco::Jet*>(&object) && !isBJet ) {
00108     bins_    = &binsUdsc_;
00109     funcEt_  = &funcEtUdsc_;
00110     funcEta_ = &funcEtaUdsc_;
00111     funcPhi_ = &funcPhiUdsc_;
00112   }
00113   else if( dynamic_cast<const reco::Jet*>(&object) && isBJet ) {
00114     bins_    = &binsB_;
00115     funcEt_  = &funcEtB_;
00116     funcEta_ = &funcEtaB_;
00117     funcPhi_ = &funcPhiB_;
00118   }
00119   else if( dynamic_cast<const reco::Muon*>(&object) || dynamic_cast<const reco::GsfElectron*>(&object) ) {
00120     bins_    = &binsLep_;
00121     funcEt_  = &funcEtLep_;
00122     funcEta_ = &funcEtaLep_;
00123     funcPhi_ = &funcPhiLep_;
00124   }
00125   else if( dynamic_cast<const reco::MET*>(&object) ) {
00126     bins_    = &binsMet_;
00127     funcEt_  = &funcEtMet_;
00128     funcEta_ = &funcEtaMet_;
00129     funcPhi_ = &funcPhiMet_;
00130   }
00131   else{
00132     throw cms::Exception("UnsupportedObject") << "The object given is not supported!\n";
00133   }
00134 
00135   int selectedBin=-1;
00136   reco::LeafCandidate candidate;
00137   for(unsigned int i=0; i<bins_->size(); ++i){
00138     StringCutObjectSelector<reco::LeafCandidate> select_(bins_->at(i));
00139     candidate = reco::LeafCandidate( 0, reco::LeafCandidate::LorentzVector(object.px(), object.py(), object.pz(), object.energy()) );
00140     if(select_(candidate)){
00141       selectedBin = i;
00142       break;
00143     }
00144   }
00145   if(selectedBin>=0){
00146     if(whichResolution == "et")       return StringObjectFunction<reco::LeafCandidate>(funcEt_ ->at(selectedBin)).operator()(candidate);
00147     else if(whichResolution == "eta") return StringObjectFunction<reco::LeafCandidate>(funcEta_->at(selectedBin)).operator()(candidate);
00148     else if(whichResolution == "phi") return StringObjectFunction<reco::LeafCandidate>(funcPhi_->at(selectedBin)).operator()(candidate);
00149     else throw cms::Exception("ProgrammingError") << "Only 'et', 'eta' and 'phi' resolutions supported!\n";
00150   }
00151   return 0;
00152 }
00153 
00154 
00155 template <class ObjectType>
00156 TMatrixD CovarianceMatrix::setupMatrix(const pat::PATObject<ObjectType>& object, TopKinFitter::Param param, std::string resolutionProvider = "")
00157 {
00158   TMatrixD CovM3 (3,3); CovM3.Zero();
00159   TMatrixD CovM4 (4,4); CovM4.Zero();
00160   TMatrixD* CovM = &CovM3;
00161   // This part is for pat objects with resolutions embedded
00162   if(object.hasKinResolution())
00163     {
00164       switch(param){
00165       case TopKinFitter::kEtEtaPhi :
00166         CovM3(0,0) = pow(object.resolEt(resolutionProvider) , 2);
00167         if( dynamic_cast<const reco::MET*>(&object) ) CovM3(1,1) = pow(9999., 2);
00168         else CovM3(1,1) = pow(object.resolEta(resolutionProvider), 2);
00169         CovM3(2,2) = pow(object.resolPhi(resolutionProvider), 2);
00170         CovM = &CovM3;
00171         break;
00172       case TopKinFitter::kEtThetaPhi :
00173         CovM3(0,0) = pow(object.resolEt(resolutionProvider)   , 2);
00174         CovM3(1,1) = pow(object.resolTheta(resolutionProvider), 2);
00175         CovM3(2,2) = pow(object.resolPhi(resolutionProvider)  , 2);
00176         CovM = &CovM3;
00177         break;
00178       case TopKinFitter::kEMom :
00179         CovM4(0,0) = pow(1, 2);
00180         CovM4(1,1) = pow(1, 2);
00181         CovM4(2,2) = pow(1, 2);
00182         CovM4(3,3) = pow(1, 2);
00183         CovM = &CovM4;
00184         break;
00185       }
00186     }
00187   // This part is for objects without resolutions embedded
00188   else
00189     {
00190       double pt = object.pt(), eta = object.eta();
00191       // if object is a jet
00192       if( dynamic_cast<const reco::Jet*>(&object) ) {
00193         res::HelperJet jetRes;
00194         switch(param){
00195         case TopKinFitter::kEMom :
00196           if(resolutionProvider == "bjets") {
00197             CovM4(0,0) = pow(jetRes.a (pt, eta, res::HelperJet::kB  ), 2); 
00198             CovM4(1,1) = pow(jetRes.b (pt, eta, res::HelperJet::kB  ), 2); 
00199             CovM4(2,2) = pow(jetRes.c (pt, eta, res::HelperJet::kB  ), 2);
00200             CovM4(3,3) = pow(jetRes.d (pt, eta, res::HelperJet::kB  ), 2);
00201           }
00202           else {
00203             CovM4(0,0) = pow(jetRes.a (pt, eta, res::HelperJet::kUds), 2);
00204             CovM4(1,1) = pow(jetRes.b (pt, eta, res::HelperJet::kUds), 2);
00205             CovM4(2,2) = pow(jetRes.c (pt, eta, res::HelperJet::kUds), 2);
00206             CovM4(3,3) = pow(jetRes.d (pt, eta, res::HelperJet::kUds), 2);
00207           }
00208           CovM = &CovM4;
00209           break;
00210         case TopKinFitter::kEtEtaPhi : 
00211           if(resolutionProvider == "bjets") {
00212             if(!binsB_.size()){
00213               CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kB  ), 2); 
00214               CovM3(1,1) = pow(jetRes.eta(pt, eta, res::HelperJet::kB  ), 2); 
00215               CovM3(2,2) = pow(jetRes.phi(pt, eta, res::HelperJet::kB  ), 2);
00216             }
00217             else{
00218               CovM3(0,0) = pow(getResolution(object, "et" , true), 2); 
00219               CovM3(1,1) = pow(getResolution(object, "eta", true), 2); 
00220               CovM3(2,2) = pow(getResolution(object, "phi", true), 2);
00221             }
00222           }
00223           else {
00224             if(!binsUdsc_.size()){
00225               CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kUds), 2);
00226               CovM3(1,1) = pow(jetRes.eta(pt, eta, res::HelperJet::kUds), 2);
00227               CovM3(2,2) = pow(jetRes.phi(pt, eta, res::HelperJet::kUds), 2);
00228             }
00229             else{
00230               CovM3(0,0) = pow(getResolution(object, "et") , 2); 
00231               CovM3(1,1) = pow(getResolution(object, "eta"), 2); 
00232               CovM3(2,2) = pow(getResolution(object, "phi"), 2);
00233             }
00234           }       
00235           CovM = &CovM3;
00236           break;
00237         case TopKinFitter::kEtThetaPhi :
00238           if(resolutionProvider == "bjets") {
00239             CovM3(0,0) = pow(jetRes.et   (pt, eta, res::HelperJet::kB  ), 2); 
00240             CovM3(1,1) = pow(jetRes.theta(pt, eta, res::HelperJet::kB  ), 2); 
00241             CovM3(2,2) = pow(jetRes.phi  (pt, eta, res::HelperJet::kB  ), 2);
00242           }
00243           else {
00244             CovM3(0,0) = pow(jetRes.et   (pt, eta, res::HelperJet::kUds), 2);
00245             CovM3(1,1) = pow(jetRes.theta(pt, eta, res::HelperJet::kUds), 2);
00246             CovM3(2,2) = pow(jetRes.phi  (pt, eta, res::HelperJet::kUds), 2);
00247           }
00248           CovM = &CovM3;
00249           break;
00250         } 
00251       }
00252       // if object is an electron
00253       else if( dynamic_cast<const reco::GsfElectron*>(&object) ) {
00254         res::HelperElectron elecRes;
00255         switch(param){
00256         case TopKinFitter::kEMom :
00257           CovM3(0,0) = pow(elecRes.a (pt, eta), 2);
00258           CovM3(1,1) = pow(elecRes.b (pt, eta), 2); 
00259           CovM3(2,2) = pow(elecRes.c (pt, eta), 2);
00260           CovM = &CovM3;
00261           break;
00262         case TopKinFitter::kEtEtaPhi :
00263           if(!binsLep_.size()){
00264             CovM3(0,0) = pow(elecRes.et (pt, eta), 2);
00265             CovM3(1,1) = pow(elecRes.eta(pt, eta), 2); 
00266             CovM3(2,2) = pow(elecRes.phi(pt, eta), 2);
00267           }
00268           else{
00269             CovM3(0,0) = pow(getResolution(object, "et") , 2);
00270             CovM3(1,1) = pow(getResolution(object, "eta"), 2);
00271             CovM3(2,2) = pow(getResolution(object, "phi"), 2);
00272           }
00273           CovM = &CovM3;
00274           break;
00275         case TopKinFitter::kEtThetaPhi :
00276           CovM3(0,0) = pow(elecRes.et   (pt, eta), 2);
00277           CovM3(1,1) = pow(elecRes.theta(pt, eta), 2); 
00278           CovM3(2,2) = pow(elecRes.phi  (pt, eta), 2);
00279           CovM = &CovM3;
00280           break;
00281         }
00282       }
00283       // if object is a muon
00284       else if( dynamic_cast<const reco::Muon*>(&object) ) {
00285         res::HelperMuon muonRes;
00286         switch(param){
00287         case TopKinFitter::kEMom :
00288           CovM3(0,0) = pow(muonRes.a (pt, eta), 2);
00289           CovM3(1,1) = pow(muonRes.b (pt, eta), 2); 
00290           CovM3(2,2) = pow(muonRes.c (pt, eta), 2);
00291           CovM = &CovM3;
00292           break;
00293         case TopKinFitter::kEtEtaPhi :
00294           if(!binsLep_.size()){
00295             CovM3(0,0) = pow(muonRes.et (pt, eta), 2);
00296             CovM3(1,1) = pow(muonRes.eta(pt, eta), 2); 
00297             CovM3(2,2) = pow(muonRes.phi(pt, eta), 2);
00298           }
00299           else{
00300             CovM3(0,0) = pow(getResolution(object, "et") , 2);
00301             CovM3(1,1) = pow(getResolution(object, "eta"), 2);
00302             CovM3(2,2) = pow(getResolution(object, "phi"), 2);
00303           }
00304           CovM = &CovM3;
00305           break;
00306         case TopKinFitter::kEtThetaPhi :
00307           CovM3(0,0) = pow(muonRes.et   (pt, eta), 2);
00308           CovM3(1,1) = pow(muonRes.theta(pt, eta), 2); 
00309           CovM3(2,2) = pow(muonRes.phi  (pt, eta), 2);
00310           CovM = &CovM3;
00311           break;
00312         }
00313       }
00314       // if object is met
00315       else if( dynamic_cast<const reco::MET*>(&object) ) {
00316         res::HelperMET metRes;
00317         switch(param){
00318         case TopKinFitter::kEMom :
00319           CovM3(0,0) = pow(metRes.a(pt), 2);
00320           CovM3(1,1) = pow(metRes.b(pt), 2);
00321           CovM3(2,2) = pow(metRes.c(pt), 2);
00322           CovM = &CovM3;
00323           break;
00324         case TopKinFitter::kEtEtaPhi :
00325           if(!binsMet_.size()){
00326             CovM3(0,0) = pow(metRes.et(pt) , 2);
00327             CovM3(1,1) = pow(        9999. , 2);
00328             CovM3(2,2) = pow(metRes.phi(pt), 2);
00329           }
00330           else{
00331             CovM3(0,0) = pow(getResolution(object, "et") , 2);
00332             CovM3(1,1) = pow(getResolution(object, "eta"), 2);
00333             CovM3(2,2) = pow(getResolution(object, "phi"), 2);
00334           }
00335           CovM = &CovM3;
00336           break;
00337         case TopKinFitter::kEtThetaPhi :
00338           CovM3(0,0) = pow(metRes.et(pt) , 2);
00339           CovM3(1,1) = pow(        9999. , 2);
00340           CovM3(2,2) = pow(metRes.phi(pt), 2);
00341           CovM = &CovM3;
00342           break;
00343         }
00344       }
00345     }
00346   return *CovM;
00347 }
00348 
00349 #endif