CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/RecoMET/METAlgorithms/src/SignAlgoResolutions.cc

Go to the documentation of this file.
00001 #include "RecoMET/METAlgorithms/interface/SignAlgoResolutions.h"
00002 // -*- C++ -*-
00003 //
00004 // Package:    METAlgorithms
00005 // Class:      SignAlgoResolutions
00006 // 
00014 //
00015 // Original Author:  Kyle Story, Freya Blekman (Cornell University)
00016 //         Created:  Fri Apr 18 11:58:33 CEST 2008
00017 // $Id: SignAlgoResolutions.cc,v 1.9 2011/08/15 12:32:03 akhukhun Exp $
00018 //
00019 //
00020 #include "FWCore/Framework/interface/EventSetup.h"
00021 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00022 #include "DataFormats/TrackReco/interface/Track.h"
00023 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
00024 #include "DataFormats/ParticleFlowReco/interface/PFBlockElement.h"
00025 
00026 
00027 #include <math.h>
00028 
00029 #include <cstdlib>
00030 #include <iostream>
00031 #include <sstream>
00032 #include <string>
00033 #include <sys/stat.h>
00034 
00035 metsig::SignAlgoResolutions::SignAlgoResolutions(const edm::ParameterSet &iConfig):
00036     functionmap_(),
00037     ptResol_(0),
00038     phiResol_(0)
00039 {
00040   addResolutions(iConfig);
00041 }
00042 
00043 
00044 
00045 double metsig::SignAlgoResolutions::eval(const resolutionType & type, const resolutionFunc & func, const double & et, const double & phi, const double & eta) const {
00046   // derive p from et and eta;
00047   double theta = 2*atan(exp(-eta));
00048   double p = et / sin(theta); // rough assumption: take e and p equivalent...
00049   return eval(type,func,et,phi,eta,p);
00050 
00051 }
00052 
00053 double metsig::SignAlgoResolutions::eval(const resolutionType & type, const resolutionFunc & func, const double & et, const double & phi, const double & eta, const double &p) const {
00054 
00055   functionPars x(4);
00056   x[0]=et;
00057   x[1]=phi;
00058   x[2]=eta;
00059   x[3]=p;
00060   //  std::cout << "getting function of type " << type << " " << func << " " << x[0] << " " << x[1] << " " << x[2] << " " << x[3] << std::endl;
00061   return getfunc(type,func,x);
00062 
00063 }
00064 
00065 metsig::SigInputObj  metsig::SignAlgoResolutions::evalPF(const reco::PFCandidate *candidate) const {
00066   double eta = candidate->eta();
00067   double phi = candidate->phi();
00068   double et = candidate->energy()*sin(candidate->theta());
00069   resolutionType thetype;
00070   std::string name;
00071   int type = candidate->particleId();
00072   switch (type) 
00073     {
00074     case 1: 
00075       thetype=PFtype1;name="PFChargedHadron"; break;
00076     case 2: 
00077       thetype=PFtype2;name="PFChargedEM"; break;
00078     case 3: 
00079       thetype=PFtype3;name="PFMuon"; break;
00080     case 4: 
00081       thetype=PFtype4;name="PFNeutralEM"; break;
00082     case 5: 
00083       thetype=PFtype5;name="PFNeutralHadron"; break;
00084     case 6: 
00085       thetype=PFtype6;name="PFtype6"; break;
00086     case 7:
00087       thetype=PFtype7;name="PFtype7"; break;
00088     default:
00089       thetype=PFtype7;name="PFunknown"; break;
00090   }
00091 
00092   double d_et=0, d_phi=0; //d_phi here is the error on phi component of the et
00093   reco::TrackRef trackRef = candidate->trackRef();
00094   if(!trackRef.isNull()){
00095     d_phi = et*trackRef->phiError();
00096     d_et  = (type==2) ?  ElectronPtResolution(candidate) : trackRef->ptError();
00097     //if(type==2) std::cout << eval(thetype,ET,et,phi,eta) << "    " << trackRef->ptError() << "    "<< ElectronPtResolution(candidate) << std::endl;
00098   }
00099   else{
00100     d_et = eval(thetype,ET,et,phi,eta);
00101     d_phi = eval(thetype,PHI,et,phi,eta);
00102   }
00103 
00104   metsig::SigInputObj resultingobj(name,et,phi,d_et,d_phi);
00105   return resultingobj;
00106 }
00107 
00108 
00109 metsig::SigInputObj
00110 metsig::SignAlgoResolutions::evalPFJet(const reco::PFJet *jet) const{
00111 
00112     double jpt  = jet->pt();
00113     double jphi = jet->phi();
00114     double jeta = jet->eta();
00115     double jdeltapt = 999.;
00116     double jdeltapphi = 999.;
00117 
00118     if(jpt<ptResolThreshold_ && jpt<20.){ //use temporary fix for low pT jets
00119         double feta = TMath::Abs(jeta);
00120         int ieta = feta<5.? int(feta/0.5) : 9; //bin size = 0.5 
00121         int ipt  = jpt>3. ? int(jpt-3./2) : 0; //bin size =2, starting from ptmin=3GeV
00122         jdeltapt   = jdpt[ieta][ipt];
00123         jdeltapphi = jpt*jdphi[ieta][ipt];
00124     }
00125     else{
00126         //use the resolution functions at |eta|=5 to avoid crash for jets with large eta.
00127         if(jeta>5) jeta=5;
00128         if(jeta<-5) jeta=-5;
00129         TF1* fPtEta  = ptResol_->parameterEta("sigma",jeta);
00130         TF1* fPhiEta = phiResol_->parameterEta("sigma",jeta);
00131         jdeltapt   = jpt>ptResolThreshold_ ? jpt*fPtEta->Eval(jpt)  : jpt*fPtEta->Eval(ptResolThreshold_);
00132         jdeltapphi = jpt>ptResolThreshold_ ? jpt*fPhiEta->Eval(jpt) : jpt*fPhiEta->Eval(ptResolThreshold_);
00133         delete fPtEta;
00134         delete fPhiEta;
00135     }
00136 
00137     std::string inputtype = "jet";
00138     metsig::SigInputObj obj_jet(inputtype,jpt,jphi,jdeltapt,jdeltapphi);
00139     //std::cout << "RESOLUTIONS JET: " << jpt << "   " << jphi<< "   " <<jdeltapt << "   " << jdeltapphi << std::endl;
00140     return obj_jet;
00141 }
00142 
00143 
00144 void metsig::SignAlgoResolutions::addResolutions(const edm::ParameterSet &iConfig){
00145     using namespace std;
00146 
00147   // Jet Resolutions - for now load from the files. Migrate to EventSetup asap.
00148   metsig::SignAlgoResolutions::initializeJetResolutions( iConfig );
00149   
00150   ptResolThreshold_ = iConfig.getParameter<double>("ptresolthreshold");
00151 
00152 
00153     //get temporary low pT pfjet resolutions
00154     for (int ieta=0; ieta<10; ieta++){
00155       jdpt[ieta] = iConfig.getParameter<std::vector<double> >(Form("jdpt%d", ieta));
00156       jdphi[ieta] = iConfig.getParameter<std::vector<double> >(Form("jdphi%d", ieta));
00157     }
00158 
00159 
00160   // for now: do this by hand - this can obviously also be done via ESSource etc.
00161   functionPars etparameters(3,0);
00162   functionPars phiparameters(1,0);
00163   // set the parameters per function:
00164   // ECAL, BARREL:
00165   std::vector<double> ebet = iConfig.getParameter<std::vector<double> >("EB_EtResPar");
00166   std::vector<double> ebphi = iConfig.getParameter<std::vector<double> >("EB_PhiResPar");
00167 
00168   etparameters[0]=ebet[0];
00169   etparameters[1]=ebet[1];
00170   etparameters[2]=ebet[2];
00171   phiparameters[0]=ebphi[0];
00172   addfunction(caloEB,ET,etparameters);
00173   addfunction(caloEB,PHI,phiparameters);
00174  // ECAL, ENDCAP:
00175   std::vector<double> eeet = iConfig.getParameter<std::vector<double> >("EE_EtResPar");
00176   std::vector<double> eephi = iConfig.getParameter<std::vector<double> >("EE_PhiResPar");
00177 
00178   etparameters[0]=eeet[0];
00179   etparameters[1]=eeet[1];
00180   etparameters[2]=eeet[2];
00181   phiparameters[0]=eephi[0];
00182   addfunction(caloEE,ET,etparameters);
00183   addfunction(caloEE,PHI,phiparameters);
00184  // HCAL, BARREL:
00185   std::vector<double> hbet = iConfig.getParameter<std::vector<double> >("HB_EtResPar");
00186   std::vector<double> hbphi = iConfig.getParameter<std::vector<double> >("HB_PhiResPar");
00187 
00188   etparameters[0]=hbet[0];
00189   etparameters[1]=hbet[1];
00190   etparameters[2]=hbet[2];
00191   phiparameters[0]=hbphi[0];
00192   addfunction(caloHB,ET,etparameters);
00193   addfunction(caloHB,PHI,phiparameters);
00194  // HCAL, ENDCAP:
00195   std::vector<double> heet = iConfig.getParameter<std::vector<double> >("HE_EtResPar");
00196   std::vector<double> hephi = iConfig.getParameter<std::vector<double> >("HE_PhiResPar");
00197 
00198   etparameters[0]=heet[0];
00199   etparameters[1]=heet[1];
00200   etparameters[2]=heet[2];
00201   phiparameters[0]=hephi[0];
00202   addfunction(caloHE,ET,etparameters);
00203   addfunction(caloHE,PHI,phiparameters);
00204  // HCAL, Outer
00205   std::vector<double> hoet = iConfig.getParameter<std::vector<double> >("HO_EtResPar");
00206   std::vector<double> hophi = iConfig.getParameter<std::vector<double> >("HO_PhiResPar");
00207 
00208 
00209   etparameters[0]=hoet[0];
00210   etparameters[1]=hoet[1];
00211   etparameters[2]=hoet[2];
00212   phiparameters[0]=hophi[0];
00213   addfunction(caloHO,ET,etparameters);
00214   addfunction(caloHO,PHI,phiparameters);
00215  // HCAL, Forward
00216   std::vector<double> hfet = iConfig.getParameter<std::vector<double> >("HF_EtResPar");
00217   std::vector<double> hfphi = iConfig.getParameter<std::vector<double> >("HF_PhiResPar");
00218 
00219   etparameters[0]=hfet[0];
00220   etparameters[1]=hfet[1];
00221   etparameters[2]=hfet[2];
00222   phiparameters[0]=hfphi[0];
00223   addfunction(caloHF,ET,etparameters);
00224   addfunction(caloHF,PHI,phiparameters);
00225 
00226 
00227   // PF objects:
00228   // type 1:
00229   std::vector<double> pf1et = iConfig.getParameter<std::vector<double> >("PF_EtResType1");
00230   std::vector<double> pf1phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType1");
00231   etparameters[0]=pf1et[0];
00232   etparameters[1]=pf1et[1];
00233   etparameters[2]=pf1et[2];
00234   phiparameters[0]=pf1phi[0];
00235   addfunction(PFtype1,ET,etparameters);
00236   addfunction(PFtype1,PHI,phiparameters);
00237 
00238   // PF objects:
00239   // type 2:
00240   std::vector<double> pf2et = iConfig.getParameter<std::vector<double> >("PF_EtResType2");
00241   std::vector<double> pf2phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType2");
00242   etparameters[0]=pf2et[0];
00243   etparameters[1]=pf2et[1];
00244   etparameters[2]=pf2et[2];
00245   phiparameters[0]=pf2phi[0];
00246   addfunction(PFtype2,ET,etparameters);
00247   addfunction(PFtype2,PHI,phiparameters);
00248 
00249   // PF objects:
00250   // type 3:
00251   std::vector<double> pf3et = iConfig.getParameter<std::vector<double> >("PF_EtResType3");
00252   std::vector<double> pf3phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType3");
00253   etparameters[0]=pf3et[0];
00254   etparameters[1]=pf3et[1];
00255   etparameters[2]=pf3et[2];
00256   phiparameters[0]=pf3phi[0];
00257   addfunction(PFtype3,ET,etparameters);
00258   addfunction(PFtype3,PHI,phiparameters);
00259 
00260   // PF objects:
00261   // type 4:
00262   std::vector<double> pf4et = iConfig.getParameter<std::vector<double> >("PF_EtResType4");
00263   std::vector<double> pf4phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType4");
00264   etparameters[0]=pf4et[0];
00265   etparameters[1]=pf4et[1];
00266   etparameters[2]=pf4et[2];
00267   //phiparameters[0]=pf4phi[0];
00268   addfunction(PFtype4,ET,etparameters);
00269   addfunction(PFtype4,PHI,pf4phi); //use the same functional form for photon phi error as for pT, pass whole vector
00270 
00271   // PF objects:
00272   // type 5:
00273   std::vector<double> pf5et = iConfig.getParameter<std::vector<double> >("PF_EtResType5");
00274   std::vector<double> pf5phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType5");
00275   etparameters[0]=pf5et[0];
00276   etparameters[1]=pf5et[1];
00277   etparameters[2]=pf5et[2];
00278   phiparameters[0]=pf5phi[0];
00279   addfunction(PFtype5,ET,etparameters);
00280   addfunction(PFtype5,PHI,pf5phi);
00281 
00282   // PF objects:
00283   // type 6:
00284   std::vector<double> pf6et = iConfig.getParameter<std::vector<double> >("PF_EtResType6");
00285   std::vector<double> pf6phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType6");
00286   etparameters[0]=pf6et[0];
00287   etparameters[1]=pf6et[1];
00288   etparameters[2]=pf6et[2];
00289   phiparameters[0]=pf6phi[0];
00290   addfunction(PFtype6,ET,etparameters);
00291   addfunction(PFtype6,PHI,phiparameters);
00292 
00293   
00294   // PF objects:
00295   // type 7:
00296   std::vector<double> pf7et = iConfig.getParameter<std::vector<double> >("PF_EtResType7");
00297   std::vector<double> pf7phi = iConfig.getParameter<std::vector<double> >("PF_PhiResType7");
00298   etparameters[0]=pf7et[0];
00299   etparameters[1]=pf7et[1];
00300   etparameters[2]=pf7et[2];
00301   phiparameters[0]=pf7phi[0];
00302   addfunction(PFtype7,ET,etparameters);
00303   addfunction(PFtype7,PHI,phiparameters);
00304 
00305   return;
00306 }
00307 
00308 void metsig::SignAlgoResolutions::addfunction(resolutionType type, resolutionFunc func, functionPars parameters){
00309 
00310   functionCombo mypair(type,func);
00311   functionmap_[mypair]=parameters;
00312 
00313 }
00314 
00315 double metsig::SignAlgoResolutions::getfunc(const metsig::resolutionType & type,const metsig::resolutionFunc & func, functionPars & x) const{
00316   
00317   double result=0;
00318   functionCombo mypair(type,func);
00319 
00320   if(functionmap_.count(mypair)==0){
00321         return result;
00322   }
00323   
00324   functionPars values = (functionmap_.find(mypair))->second;
00325   switch ( func ){
00326   case metsig::ET :
00327     return EtFunction(x,values);
00328   case metsig::PHI :
00329     return PhiFunction(x,values);
00330   case metsig::TRACKP :
00331     return PFunction(x,values);
00332   case metsig::CONSTPHI :
00333     return PhiConstFunction(x,values);
00334   }
00335   
00336   //  std::cout << "returning function " << type << " " << func << " " << result << " " << x[0] << std::endl; 
00337 
00338   return result;
00339 }
00340 
00341 double metsig::SignAlgoResolutions::EtFunction( const functionPars &x, const functionPars & par) const
00342 {
00343   if(par.size()<3)
00344     return 0.;
00345   if(x.size()<1)
00346     return 0.;
00347   double et=x[0];
00348   if(et<=0.)
00349     return 0.;
00350   double result = et*sqrt((par[2]*par[2])+(par[1]*par[1]/et)+(par[0]*par[0]/(et*et)));
00351   return result;
00352 }
00353 
00354 
00355 double metsig::SignAlgoResolutions::PhiFunction(const functionPars &x,const  functionPars & par) const
00356 {
00357   double et=x[0];
00358   if(et<=0.){
00359     return 0.;
00360   }
00361 
00362   //if 1 parameter is C provided, returns C*pT, if three parameters N, S, C are provided, it returns the usual resolution value, as for sigmaPt
00363   if(par.size()!=1 && par.size()!=3){//only 1 or 3 parameters supported for phi function
00364       return 0.;
00365   }
00366   else if(par.size()==1){
00367     return par[0]*et;
00368   }
00369   else{
00370     return et*sqrt((par[2]*par[2])+(par[1]*par[1]/et)+(par[0]*par[0]/(et*et)));
00371   }
00372 
00373 }
00374 double metsig::SignAlgoResolutions::PFunction(const functionPars &x, const functionPars & par) const
00375 {
00376   // not currently implemented
00377   return 0;
00378 }
00379 
00380 double metsig::SignAlgoResolutions::PhiConstFunction(const functionPars& x, const functionPars &par) const
00381 {
00382   return par[0];
00383 }
00384 
00385 void
00386 metsig::SignAlgoResolutions::initializeJetResolutions( const edm::ParameterSet &iConfig ) {
00387   
00388   using namespace std;
00389   
00390   // only reinitialize the resolutsion if the pointers are zero
00391   if ( ptResol_ == 0 ) {
00392     string resolutionsAlgo  = iConfig.getParameter<std::string>("resolutionsAlgo");     
00393     string resolutionsEra   = iConfig.getParameter<std::string>("resolutionsEra");     
00394 
00395     string cmssw_base(getenv("CMSSW_BASE"));
00396     string cmssw_release_base(getenv("CMSSW_RELEASE_BASE"));
00397     string path = cmssw_base + "/src/CondFormats/JetMETObjects/data";
00398     struct stat st;
00399     if (stat(path.c_str(),&st)!=0) {
00400       path = cmssw_release_base + "/src/CondFormats/JetMETObjects/data";
00401     }
00402     if (stat(path.c_str(),&st)!=0) {
00403       cerr<<"ERROR: tried to set path but failed, abort."<<endl;
00404     }    
00405     string era(resolutionsEra);
00406     string alg(resolutionsAlgo);
00407     string ptFileName  = path + "/" + era + "_PtResolution_" +alg+".txt";
00408     string phiFileName = path + "/" + era + "_PhiResolution_"+alg+".txt";
00409     
00410     ptResol_ = new JetResolution(ptFileName,false);
00411     phiResol_ = new JetResolution(phiFileName,false);
00412   }
00413 }
00414 
00415 double 
00416 metsig::SignAlgoResolutions::ElectronPtResolution(const reco::PFCandidate *c) const{
00417 
00418     double eta = c->eta();
00419     double energy = c->energy();
00420     double dEnergy = pfresol_->getEnergyResolutionEm(energy, eta);
00421 
00422     return dEnergy/cosh(eta);
00423 }
00424