CMS 3D CMS Logo

PFEnergyCalibration Class Reference

#include <RecoParticleFlow/PFClusterTools/interface/PFEnergyCalibration.h>

List of all members.

Public Member Functions

double energyEm (const reco::PFCluster &clusterEcal, std::vector< double > &EclustersPS1, std::vector< double > &EclustersPS2)
double energyEm (double uncalibratedEnergyECAL, double eta=0, double phi=0) const
double energyEmHad (double uncalibratedEnergyECAL, double uncalibratedEnergyHCAL, double eta=0, double phi=0) const
double energyHad (double uncalibratedEnergyHCAL, double eta=0, double phi=0) const
double paramECAL_offset () const
double paramECAL_slope () const
double paramECALplusHCAL_offset () const
double paramECALplusHCAL_slopeECAL () const
double paramECALplusHCAL_slopeHCAL () const
double paramHCAL_damping () const
double paramHCAL_offset () const
double paramHCAL_slope () const
 PFEnergyCalibration (double e_slope, double e_offset, double eh_eslope, double eh_hslope, double eh_offset, double h_slope, double h_offset, double h_damping)
 PFEnergyCalibration ()
void setCalibrationParametersEm (double paramECAL_slope, double paramECAL_offset)
 ~PFEnergyCalibration ()

Protected Attributes

double paramECAL_offset_
double paramECAL_slope_
double paramECALplusHCAL_offset_
double paramECALplusHCAL_slopeECAL_
double paramECALplusHCAL_slopeHCAL_
double paramHCAL_damping_
double paramHCAL_offset_
double paramHCAL_slope_

Private Member Functions

double Alpha (double eta)
double Beta (double E, double eta)
double CorrBarrel (double E, double eta)
double CorrEta (double eta)
double CorrPhi (double phi, double eta)
double dCrackPhi (double phi, double eta)
double Ecorr (double eEcal, double ePS1, double ePS2, double eta, double phi)
double EcorrBarrel (double E, double eta, double phi)
double EcorrPS (double eEcal, double ePS1, double ePS2, double etaEcal)
double EcorrPS_ePSNil (double eEcal, double eta)
double EcorrZoneAfterPS (double E, double eta)
double EcorrZoneBeforePS (double E, double eta)
double Gamma (double etaEcal)
double minimum (double a, double b)

Friends

std::ostream & operator<< (std::ostream &out, const PFEnergyCalibration &calib)


Detailed Description

Definition at line 34 of file PFEnergyCalibration.h.


Constructor & Destructor Documentation

PFEnergyCalibration::PFEnergyCalibration (  ) 

Definition at line 8 of file PFEnergyCalibration.cc.

References paramECAL_offset_, paramECAL_slope_, paramECALplusHCAL_offset_, paramECALplusHCAL_slopeECAL_, paramECALplusHCAL_slopeHCAL_, paramHCAL_damping_, paramHCAL_offset_, and paramHCAL_slope_.

00008                                          {
00009 
00010 //--- initialize calibration parameters
00011 //    for energy correction applied to energy deposits of electrons 
00012 //    and photons in ECAL
00013 //  paramECAL_slope_ = 1.;
00014 //  paramECAL_offset_ = 0.;
00015 //  
00018 //  paramHCAL_slope_ = 2.17;
00019 //  paramHCAL_offset_ = 1.73;
00020 //  paramHCAL_damping_ = 2.49;
00021 //
00024 //  paramECALplusHCAL_slopeECAL_ = 1.05;
00025 //  paramECALplusHCAL_slopeHCAL_ = 1.06;
00026 //  paramECALplusHCAL_offset_ = 6.11;
00027         
00028   paramECAL_slope_ = 1.;
00029   paramECAL_offset_ = 0.;
00030   
00031 //--- initialize calibration parameters
00032 //    for energy correction applied to energy deposits of hadrons in HCAL
00033   paramHCAL_slope_ = 1.0;
00034   paramHCAL_offset_ = 0.0;
00035   paramHCAL_damping_ = 1.0;
00036 
00037 //--- initialize calibration parameters
00038 //    for energy correction applied to combined energy deposits of hadrons in HCAL and ECAL
00039   paramECALplusHCAL_slopeECAL_ = 1.0;
00040   paramECALplusHCAL_slopeHCAL_ = 1.0;
00041   paramECALplusHCAL_offset_ = 0.0;
00042 }

PFEnergyCalibration::PFEnergyCalibration ( double  e_slope,
double  e_offset,
double  eh_eslope,
double  eh_hslope,
double  eh_offset,
double  h_slope,
double  h_offset,
double  h_damping 
)

Definition at line 45 of file PFEnergyCalibration.cc.

00052                                                             :
00053   
00054   paramECAL_slope_(e_slope),
00055   paramECAL_offset_(e_offset),
00056   paramECALplusHCAL_slopeECAL_(eh_eslope),
00057   paramECALplusHCAL_slopeHCAL_(eh_hslope),
00058   paramECALplusHCAL_offset_(eh_offset),
00059   paramHCAL_slope_(h_slope),
00060   paramHCAL_offset_(h_offset),
00061   paramHCAL_damping_(h_damping) {}

PFEnergyCalibration::~PFEnergyCalibration (  ) 

Definition at line 81 of file PFEnergyCalibration.cc.

00082 {
00083 //--- nothing to be done yet  
00084 }


Member Function Documentation

double PFEnergyCalibration::Alpha ( double  eta  )  [private]

Definition at line 351 of file PFEnergyCalibration.cc.

References norm, p1, p2, and HLT_VtxMuL3::result.

Referenced by EcorrPS().

00351                                      {
00352 
00353   //Energy dependency
00354   static double p0 = 5.97621e-01;
00355 
00356   //Eta dependency
00357   static double p1 =-1.86407e-01;
00358   static double p2 = 3.85197e-01; 
00359 
00360   //so that <feta()> = 1
00361   static double norm = (p1+p2*(2.6+1.656)/2);
00362 
00363   double result = p0*(p1+p2*eta)/norm;
00364 
00365   return result;
00366 }

double PFEnergyCalibration::Beta ( double  E,
double  eta 
) [private]

Definition at line 369 of file PFEnergyCalibration.cc.

References norm, p1, p2, p3, p4, p5, and HLT_VtxMuL3::result.

Referenced by EcorrPS().

00369                                               {
00370 
00371  //Energy dependency
00372   static double p0 = 0.032;
00373   static double p1 = 9.70394e-02;
00374   static double p2 = 2.23072e+01;
00375   static double p3 = 100;
00376 
00377   //Eta dependency
00378   static double p4 = 1.02496e+00 ;
00379   static double p5 = -4.40176e-03 ;
00380 
00381   //so that <feta()> = 1
00382   static double norm = (p4+p5*(2.6+1.656)/2);
00383 
00384   double result = (1.0012+p0*TMath::Exp(-E/p3)+p1*TMath::Exp(-E/p2))*(p4+p5*eta)/norm;                    
00385   return result;
00386 }

double PFEnergyCalibration::CorrBarrel ( double  E,
double  eta 
) [private]

Definition at line 315 of file PFEnergyCalibration.cc.

References p1, p2, p3, p4, p5, p6, p7, p8, and HLT_VtxMuL3::result.

Referenced by EcorrBarrel().

00315                                                     {
00316 
00317   //Energy dependency
00318   static double p0=1.00000e+00;
00319   static double p1=3.27753e+01;
00320   static double p2=2.28552e-02;
00321   static double p3=3.06139e+00;
00322   static double p4=2.25135e-01;
00323   static double p5=1.47824e+00;
00324   static double p6=1.09e-02;
00325   static double p7=4.19343e+01;
00326 
00327   //Eta dependency
00328   static double p8=2.705593e-03;
00329 
00330   double result = (p0+1/(p1+p2*TMath::Power(E,p3))+p4*TMath::Exp(-E/p5)+p6*TMath::Exp(-E*E/(p7*p7)))*(1+p8*eta*eta);
00331 
00332   return result;
00333 }

double PFEnergyCalibration::CorrEta ( double  eta  )  [private]

Definition at line 288 of file PFEnergyCalibration.cc.

References a, e, i, m, HLT_VtxMuL3::result, s, and ss.

Referenced by EcorrBarrel().

00288                                       {
00289   
00290   // we use a gaussian with a screwness for each of the 5 |eta|-cracks
00291   static std::vector<double> a;  //amplitude
00292   static std::vector<double> m;  //mean
00293   static std::vector<double> s;  //sigma
00294   static std::vector<double> sa; // screwness amplitude
00295   static std::vector<double> ss; // screwness sigma
00296 
00297   if(a.size()==0)
00298     {
00299       a.push_back(6.13349e-01) ;a.push_back(5.08146e-01)  ;a.push_back(4.44480e-01) ;a.push_back(3.3487e-01)   ;a.push_back(7.65627e-01) ;
00300       m.push_back(-1.79514e-02);m.push_back(4.44747e-01)  ;m.push_back(7.92824e-01) ;m.push_back(1.14090e+00)  ;m.push_back(1.47464e+00) ;
00301       s.push_back(7.92382e-03) ;s.push_back(3.06028e-03)  ;s.push_back(3.36139e-03) ;s.push_back(3.94521e-03)  ;s.push_back(8.63950e-04) ;
00302       sa.push_back(1.27228e+01);sa.push_back(3.81517e-02) ;sa.push_back(1.63507e-01);sa.push_back(-6.56480e-02);sa.push_back(1.87160e-01);
00303       ss.push_back(5.48753e-02);ss.push_back(-1.00223e-02);ss.push_back(2.22866e-03);ss.push_back(4.26288e-04) ;ss.push_back(2.67937e-03);
00304     }
00305  double result = 1;
00306 
00307  for(unsigned i=0;i<=4;i++) result+=a[i]*TMath::Gaus(eta,m[i],s[i])*(1+sa[i]*TMath::Sign(1.,eta-m[i])*TMath::Exp(-TMath::Abs(eta-m[i])/ss[i]));
00308 
00309   return result;
00310 }

double PFEnergyCalibration::CorrPhi ( double  phi,
double  eta 
) [private]

Definition at line 262 of file PFEnergyCalibration.cc.

References dCrackPhi(), p1, p2, p3, p4, p5, p6, p7, p8, and HLT_VtxMuL3::result.

Referenced by EcorrBarrel().

00262                                                    {
00263 
00264   // we use 3 gaussians to correct the phi-cracks effect
00265   static double p1=   5.59379e-01;
00266   static double p2=   -1.26607e-03;
00267   static double p3=  9.61133e-04;
00268 
00269   static double p4=   1.81691e-01;
00270   static double p5=   -4.97535e-03;
00271   static double p6=   1.31006e-03;
00272 
00273   static double p7=   1.38498e-01;
00274   static double p8=   1.18599e-04;
00275   static double p9= 2.01858e-03;
00276   
00277 
00278   double dminphi = dCrackPhi(phi,eta);
00279   
00280   double result = (1+p1*TMath::Gaus(dminphi,p2,p3)+p4*TMath::Gaus(dminphi,p5,p6)+p7*TMath::Gaus(dminphi,p8,p9));
00281 
00282   return result;
00283 }   

double PFEnergyCalibration::dCrackPhi ( double  phi,
double  eta 
) [private]

Definition at line 210 of file PFEnergyCalibration.cc.

References GenMuonPlsPt100GeV_cfg::cout, lat::endl(), i, m, minimum(), and pi.

Referenced by CorrPhi().

00210                                                     {
00211 
00212   static double pi= M_PI;// 3.14159265358979323846;
00213   
00214   //Location of the 18 phi-cracks
00215   static std::vector<double> cPhi;
00216   if(cPhi.size()==0)
00217     {
00218       cPhi.resize(18,0);
00219       cPhi[0]=2.97025;
00220       for(unsigned i=1;i<=17;++i) cPhi[i]=cPhi[0]-2*i*pi/18;
00221     }
00222 
00223   //Shift of this location if eta<0
00224   static double delta_cPhi=0.00638;
00225 
00226   double m; //the result
00227 
00228   //the location is shifted
00229   if(eta<0) phi +=delta_cPhi;
00230 
00231   if (phi>=-pi && phi<=pi){
00232 
00233     //the problem of the extrema
00234     if (phi<cPhi[17] || phi>=cPhi[0]){
00235       if (phi<0) phi+= 2*pi;
00236       m = minimum(phi -cPhi[0],phi-cPhi[17]-2*pi);              
00237     }
00238 
00239     //between these extrema...
00240     else{
00241       bool OK = false;
00242       unsigned i=16;
00243       while(!OK){
00244         if (phi<cPhi[i]){
00245           m=minimum(phi-cPhi[i+1],phi-cPhi[i]);
00246           OK=true;
00247         }
00248         else i-=1;
00249       }
00250     }
00251   }
00252   else{
00253     m=0.;        //if there is a problem, we assum that we are in a crack
00254     std::cout<<"Problem in dminphi"<<std::endl;
00255   }
00256   if(eta<0) m=-m;   //because of the disymetry
00257   return m;
00258 }

double PFEnergyCalibration::Ecorr ( double  eEcal,
double  ePS1,
double  ePS2,
double  eta,
double  phi 
) [private]

Definition at line 532 of file PFEnergyCalibration.cc.

References EcorrBarrel(), EcorrPS(), EcorrPS_ePSNil(), EcorrZoneAfterPS(), EcorrZoneBeforePS(), and HLT_VtxMuL3::result.

Referenced by energyEm().

00532                                                                                      {
00533 
00534   static double endBarrel=1.48;
00535   static double beginingPS=1.65;
00536   static double endPS=2.6;
00537   static double endEndCap=2.98;
00538  
00539   double result=0;
00540 
00541   eta=TMath::Abs(eta);
00542 
00543   if(eEcal>0){
00544     if(eta <= endBarrel)                         result = EcorrBarrel(eEcal,eta,phi);
00545     else if(eta <= beginingPS)                   result = EcorrZoneBeforePS(eEcal,eta);
00546     else if((eta < endPS) && ePS1==0 && ePS2==0) result = EcorrPS_ePSNil(eEcal,eta);
00547     else if(eta < endPS)                         result = EcorrPS(eEcal,ePS1,ePS2,eta);
00548     else if(eta < endEndCap)                     result = EcorrZoneAfterPS(eEcal,eta); 
00549     else result =eEcal;
00550   }
00551   else result = eEcal;// useful if eEcal=0 or eta>2.98
00552   return result;
00553 }

double PFEnergyCalibration::EcorrBarrel ( double  E,
double  eta,
double  phi 
) [private]

Definition at line 419 of file PFEnergyCalibration.cc.

References CorrBarrel(), CorrEta(), CorrPhi(), and HLT_VtxMuL3::result.

Referenced by Ecorr().

00419                                                                 {
00420 
00421   double result = E*CorrBarrel(E,eta)*CorrEta(eta)*CorrPhi(phi,eta);
00422 
00423   return result;
00424 }

double PFEnergyCalibration::EcorrPS ( double  eEcal,
double  ePS1,
double  ePS2,
double  etaEcal 
) [private]

Definition at line 455 of file PFEnergyCalibration.cc.

References Alpha(), Beta(), e, Gamma(), p1, p2, p3, p4, and HLT_VtxMuL3::result.

Referenced by Ecorr().

00455                                                                                 {
00456 
00457   // gives the good weights to each subdetector
00458   double E = Beta(1.0155*eEcal+0.025*(ePS1+0.5976*ePS2)/9e-5,etaEcal)*eEcal+Gamma(etaEcal)*(ePS1+Alpha(etaEcal)*ePS2)/9e-5 ;
00459 
00460   //Correction of the residual energy dependency
00461   static double p0 = 1.00;
00462   static double p1 = 2.18;
00463   static double p2 =1.94;
00464   static double p3 =4.13;
00465   static double p4 =1.127;
00466 
00467   double result = E*(p0+p1*TMath::Exp(-E/p2)-p3*TMath::Exp(-E/p4));
00468 
00469   return result;
00470 } 

double PFEnergyCalibration::EcorrPS_ePSNil ( double  eEcal,
double  eta 
) [private]

Definition at line 476 of file PFEnergyCalibration.cc.

References norm, p1, p2, p3, p4, p5, and HLT_VtxMuL3::result.

Referenced by Ecorr().

00476                                                           {
00477 
00478   //Energy dependency
00479   static double p0= 1.02;
00480   static double p1= 0.165;
00481   static double p2= 6.5 ;
00482   static double p3=  2.1 ;
00483 
00484   //Eta dependency
00485   static double p4 = 1.02496e+00 ;
00486   static double p5 = -4.40176e-03 ;
00487 
00488   //so that <feta()> = 1
00489   static double norm = (p4+p5*(2.6+1.656)/2);
00490 
00491   double result = eEcal*(p0+p1*TMath::Exp(-TMath::Abs(eEcal-p3)/p2))*(p4+p5*eta)/norm;
00492                   
00493   return result;
00494 }

double PFEnergyCalibration::EcorrZoneAfterPS ( double  E,
double  eta 
) [private]

Definition at line 499 of file PFEnergyCalibration.cc.

References norm, p1, p2, p3, p4, p5, p6, p7, p8, and HLT_VtxMuL3::result.

Referenced by Ecorr().

00499                                                          {
00500 
00501   //Energy dependency
00502   static double p0 =1; 
00503   static double p1 = 0.058;
00504   static double p2 =12.5;
00505   static double p3 =-1.05444e+00;
00506   static double p4 =-5.39557e+00;
00507   static double p5 =8.38444e+00;
00508   static double p6 = 6.10998e-01  ;
00509 
00510   //Eta dependency
00511   static double p7 =1.06161e+00;
00512   static double p8 = 0.41;
00513   static double p9 =2.918;
00514   static double p10 =0.0181;
00515   static double p11= 2.05;
00516   static double p12 =2.99;
00517   static double p13=0.0287;
00518 
00519   //so that <feta()> = 1
00520   static double norm=1.045;
00521 
00522   double result = E*(p0+p1*TMath::Exp(-(E-p3)/p2)+1/(p4+p5*TMath::Power(E,p6)))*(p7+p8*TMath::Gaus(eta,p9,p10)+p11*TMath::Gaus(eta,p12,p13))/norm;
00523   return result;
00524 }

double PFEnergyCalibration::EcorrZoneBeforePS ( double  E,
double  eta 
) [private]

Definition at line 429 of file PFEnergyCalibration.cc.

References norm, p1, p2, p3, p4, p5, p6, p7, and HLT_VtxMuL3::result.

Referenced by Ecorr().

00429                                                           {
00430 
00431  //Energy dependency
00432   static double p0 =1; 
00433   static double p1 =0.18;
00434   static double p2 =8.;
00435 
00436   //Eta dependency
00437   static double p3 =0.3;
00438   static double p4 =1.11;
00439   static double p5 =0.025;
00440   static double p6 =1.49;
00441   static double p7 =0.6;
00442 
00443   //so that <feta()> = 1
00444   static double norm = 1.21;
00445 
00446   double result = E*(p0+p1*TMath::Exp(-E/p2))*(p3+p4*TMath::Gaus(eta,p6,p5)+p7*eta)/norm;
00447 
00448   return result;
00449 }

double PFEnergyCalibration::energyEm ( const reco::PFCluster clusterEcal,
std::vector< double > &  EclustersPS1,
std::vector< double > &  EclustersPS2 
)

Definition at line 103 of file PFEnergyCalibration.cc.

References reco::PFCluster::calculatePositionREP(), GenMuonPlsPt100GeV_cfg::cout, Ecorr(), lat::endl(), reco::PFCluster::energy(), eta, i, phi, and reco::PFCluster::positionREP().

00103                                                                                                                              {
00104   double eEcal = clusterEcal.energy();
00105   //temporaty ugly fix
00106   reco::PFCluster myPFCluster=clusterEcal;
00107   myPFCluster.calculatePositionREP();
00108   double eta = myPFCluster.positionREP().eta();
00109   double phi = myPFCluster.positionREP().phi();
00110 
00111   double ePS1 = 0;
00112   double ePS2 = 0;
00113 
00114   for(unsigned i=0;i<EclustersPS1.size();i++) ePS1 += EclustersPS1[i];
00115   for(unsigned i=0;i<EclustersPS2.size();i++) ePS2 += EclustersPS2[i];
00116 
00117   double calibrated = Ecorr(eEcal,ePS1,ePS2,eta,phi);
00118   if(eEcal!=0 && calibrated==0) std::cout<<"Eecal = "<<eEcal<<"  eta = "<<eta<<"  phi = "<<phi<<std::endl; 
00119   return calibrated; 
00120 }

double PFEnergyCalibration::energyEm ( double  uncalibratedEnergyECAL,
double  eta = 0,
double  phi = 0 
) const

Definition at line 88 of file PFEnergyCalibration.cc.

References paramECAL_offset_, and paramECAL_slope_.

Referenced by PFElectronAlgo::SetCandidates(), and PFElectronAlgo::SetIDOutputs().

00089                                                             {
00090 
00091   //--- apply calibration correction
00092   //    for energy deposits of electrons and photons in ECAL
00093   //    (eta and phi dependence not implemented yet)
00094   
00095   double calibrated = paramECAL_slope_*uncalibratedEnergyECAL;
00096   calibrated += paramECAL_offset_;
00097 
00098   return calibrated;
00099 }

double PFEnergyCalibration::energyEmHad ( double  uncalibratedEnergyECAL,
double  uncalibratedEnergyHCAL,
double  eta = 0,
double  phi = 0 
) const

Definition at line 142 of file PFEnergyCalibration.cc.

References paramECALplusHCAL_offset_, paramECALplusHCAL_slopeECAL_, and paramECALplusHCAL_slopeHCAL_.

00144                                                                {
00145 //--- apply calibration correction
00146 //    for energy deposits of hadrons in ECAL and HCAL
00147 //    (eta and phi dependence not implemented yet)
00148 
00149   double calibrated = paramECALplusHCAL_slopeECAL_*uncalibratedEnergyECAL;
00150   calibrated += paramECALplusHCAL_slopeHCAL_*uncalibratedEnergyHCAL;
00151   calibrated += paramECALplusHCAL_offset_;
00152 
00153   return calibrated;
00154 }

double PFEnergyCalibration::energyHad ( double  uncalibratedEnergyHCAL,
double  eta = 0,
double  phi = 0 
) const

Definition at line 124 of file PFEnergyCalibration.cc.

References funct::exp(), paramHCAL_damping_, paramHCAL_offset_, and paramHCAL_slope_.

00125                                                              {
00126   
00127   //--- apply calibration correction
00128   //    for energy deposits of hadrons in HCAL
00129   //    (eta and phi dependence not implemented yet)
00130   
00131   double numerator = paramHCAL_slope_*uncalibratedEnergyHCAL;
00132   numerator += paramHCAL_offset_;
00133 
00134   double denominator = 1 + exp(paramHCAL_damping_/uncalibratedEnergyHCAL);
00135 
00136 
00137   return numerator/denominator;
00138 }

double PFEnergyCalibration::Gamma ( double  etaEcal  )  [private]

Definition at line 390 of file PFEnergyCalibration.cc.

References norm, p1, p2, and HLT_VtxMuL3::result.

Referenced by EcorrPS().

00390                                          {
00391 
00392  //Energy dependency
00393   static double p0 = 2.49752e-02;
00394 
00395   //Eta dependency
00396   static double p1 = 6.48816e-02;
00397   static double p2 = -1.59517e-02; 
00398  
00399   //so that <feta()> = 1
00400   static double norm = (p1+p2*(2.6+1.656)/2);
00401 
00402   double result = p0*(p1+p2*etaEcal)/norm;                                        
00403 
00404   return result;
00405 }

double PFEnergyCalibration::minimum ( double  a,
double  b 
) [private]

Definition at line 202 of file PFEnergyCalibration.cc.

Referenced by dCrackPhi().

00202                                              {
00203   if(TMath::Abs(b)<TMath::Abs(a)) a=b;
00204   return a;
00205 }

double PFEnergyCalibration::paramECAL_offset (  )  const [inline]

Definition at line 73 of file PFEnergyCalibration.h.

References paramECAL_offset_.

00073 {return paramECAL_offset_;} 

double PFEnergyCalibration::paramECAL_slope (  )  const [inline]

Definition at line 71 of file PFEnergyCalibration.h.

References paramECAL_slope_.

00071 {return  paramECAL_slope_;} 

double PFEnergyCalibration::paramECALplusHCAL_offset (  )  const [inline]

Definition at line 83 of file PFEnergyCalibration.h.

References paramECALplusHCAL_offset_.

00083 {return paramECALplusHCAL_offset_;} 

double PFEnergyCalibration::paramECALplusHCAL_slopeECAL (  )  const [inline]

Definition at line 75 of file PFEnergyCalibration.h.

References paramECALplusHCAL_slopeECAL_.

00075                                              {
00076     return paramECALplusHCAL_slopeECAL_;
00077   } 

double PFEnergyCalibration::paramECALplusHCAL_slopeHCAL (  )  const [inline]

Definition at line 79 of file PFEnergyCalibration.h.

References paramECALplusHCAL_slopeHCAL_.

00079                                              {
00080     return paramECALplusHCAL_slopeHCAL_;
00081   } 

double PFEnergyCalibration::paramHCAL_damping (  )  const [inline]

Definition at line 87 of file PFEnergyCalibration.h.

References paramHCAL_damping_.

00087 {return paramHCAL_damping_;} 

double PFEnergyCalibration::paramHCAL_offset (  )  const [inline]

Definition at line 86 of file PFEnergyCalibration.h.

References paramHCAL_offset_.

00086 {return paramHCAL_offset_;} 

double PFEnergyCalibration::paramHCAL_slope (  )  const [inline]

Definition at line 85 of file PFEnergyCalibration.h.

References paramHCAL_slope_.

00085 {return paramHCAL_slope_;} 

void PFEnergyCalibration::setCalibrationParametersEm ( double  paramECAL_slope,
double  paramECAL_offset 
)

Definition at line 68 of file PFEnergyCalibration.cc.

References paramECAL_offset_, and paramECAL_slope_.

00069                                                                               {
00070 
00071 //--- set calibration parameters for energy deposits of electrons 
00072 //    and photons in ECAL;
00073 //    this member function is needed by PFRootEvent
00074 
00075   paramECAL_slope_ = paramECAL_slope;
00076   paramECAL_offset_ = paramECAL_offset;
00077 }


Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  out,
const PFEnergyCalibration calib 
) [friend]

Definition at line 156 of file PFEnergyCalibration.cc.

00157                                                            {
00158 
00159   if(!out ) return out;
00160 
00161   out<<"PFEnergyCalibration -- "<<endl;
00162   out<<"ecal      = "<<calib.paramECAL_slope_
00163      <<" x E + "<< calib.paramECAL_offset_<<endl;
00164   out<<"hcal only = <add formula>"
00165      <<calib.paramHCAL_slope_<<","
00166      <<calib.paramHCAL_offset_<<","
00167      <<calib.paramHCAL_damping_<<endl;
00168   out<<"ecal+hcal = "<<calib.paramECALplusHCAL_slopeECAL_<<" x E_e + "
00169      <<calib.paramECALplusHCAL_slopeHCAL_<<" x E_h + "
00170      <<calib.paramECALplusHCAL_offset_<<endl;
00171 
00172   return out;
00173 }


Member Data Documentation

double PFEnergyCalibration::paramECAL_offset_ [protected]

Definition at line 95 of file PFEnergyCalibration.h.

Referenced by energyEm(), operator<<(), paramECAL_offset(), PFEnergyCalibration(), and setCalibrationParametersEm().

double PFEnergyCalibration::paramECAL_slope_ [protected]

Definition at line 94 of file PFEnergyCalibration.h.

Referenced by energyEm(), operator<<(), paramECAL_slope(), PFEnergyCalibration(), and setCalibrationParametersEm().

double PFEnergyCalibration::paramECALplusHCAL_offset_ [protected]

Definition at line 99 of file PFEnergyCalibration.h.

Referenced by energyEmHad(), operator<<(), paramECALplusHCAL_offset(), and PFEnergyCalibration().

double PFEnergyCalibration::paramECALplusHCAL_slopeECAL_ [protected]

Definition at line 97 of file PFEnergyCalibration.h.

Referenced by energyEmHad(), operator<<(), paramECALplusHCAL_slopeECAL(), and PFEnergyCalibration().

double PFEnergyCalibration::paramECALplusHCAL_slopeHCAL_ [protected]

Definition at line 98 of file PFEnergyCalibration.h.

Referenced by energyEmHad(), operator<<(), paramECALplusHCAL_slopeHCAL(), and PFEnergyCalibration().

double PFEnergyCalibration::paramHCAL_damping_ [protected]

Definition at line 103 of file PFEnergyCalibration.h.

Referenced by energyHad(), operator<<(), paramHCAL_damping(), and PFEnergyCalibration().

double PFEnergyCalibration::paramHCAL_offset_ [protected]

Definition at line 102 of file PFEnergyCalibration.h.

Referenced by energyHad(), operator<<(), paramHCAL_offset(), and PFEnergyCalibration().

double PFEnergyCalibration::paramHCAL_slope_ [protected]

Definition at line 101 of file PFEnergyCalibration.h.

Referenced by energyHad(), operator<<(), paramHCAL_slope(), and PFEnergyCalibration().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:29:43 2009 for CMSSW by  doxygen 1.5.4