Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:26

0001 #ifndef RecoParticleFlow_PFClusterTools_PFEnergyCalibration_h
0002 #define RecoParticleFlow_PFClusterTools_PFEnergyCalibration_h
0003 
0004 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0005 #include "CondFormats/PhysicsToolsObjects/interface/PerformancePayloadFromTFormula.h"
0006 
0007 #include "CondFormats/ESObjects/interface/ESEEIntercalibConstants.h"
0008 #include "CondFormats/ESObjects/interface/ESChannelStatus.h"
0009 
0010 #include <TF1.h>
0011 
0012 // -*- C++ -*-
0013 //
0014 // Package:    PFClusterTools
0015 // Class:      PFEnergyCalibration
0016 //
0017 /**\class
0018 
0019  Description: An auxiliary class of the Particle-Flow algorithm,
0020               for the Calibration of Energy Deposits in ECAL and HCAL
0021 
0022  Implementation:
0023      Original Implementation of Calibration functions in PFAlgo/PFBlock 
0024      by Colin Bernet;
0025      Code moved into separate Class by Christian Veelken 
0026 
0027  Modification 
0028      To include energy-dependent and angular-dependent calibration
0029      By Patrick Janot
0030  
0031 */
0032 //
0033 // Original Author:  Christian Veelken
0034 //         Created:  Tue Aug  8 16:26:18 CDT 2006
0035 //
0036 //
0037 
0038 #include <iostream>
0039 
0040 //#include "FWCore/ParameterSet/interface/ParameterSet.h"
0041 
0042 class PFEnergyCalibration {
0043 public:
0044   PFEnergyCalibration();
0045 
0046   // ecal calibration for photons
0047   double energyEm(const reco::PFCluster& clusterEcal, double ePS1, double ePS2, bool crackCorrection = true) const;
0048 
0049   struct CalibratedEndcapPFClusterEnergies {
0050     double clusterEnergy = 0.;
0051     double ps1Energy = 0.;
0052     double ps2Energy = 0.;
0053   };
0054 
0055   CalibratedEndcapPFClusterEnergies calibrateEndcapClusterEnergies(
0056       reco::PFCluster const& eeCluster,
0057       std::vector<reco::PFCluster const*> const& psClusterPointers,
0058       ESChannelStatus const& channelStatus,
0059       bool applyCrackCorrections) const;
0060 
0061   // ECAL+HCAL (abc) calibration, with E and eta dependent coefficients, for hadrons
0062   void energyEmHad(double t, double& e, double& h, double eta, double phi) const;
0063 
0064   // Set the run-dependent calibration functions from the global tag
0065   void setCalibrationFunctions(const PerformancePayloadFromTFormula* thePFCal) { pfCalibrations = thePFCal; }
0066 
0067   void initAlphaGamma_ESplanes_fromDB(const ESEEIntercalibConstants* esEEInterCalib) {
0068     esEEInterCalib_ = esEEInterCalib;
0069   }
0070 
0071   friend std::ostream& operator<<(std::ostream& out, const PFEnergyCalibration& calib);
0072 
0073 private:
0074   // ecal calibration for photons
0075   double energyEm(const reco::PFCluster& clusterEcal,
0076                   double ePS1,
0077                   double ePS2,
0078                   double& ps1,
0079                   double& ps2,
0080                   bool crackCorrection = true) const;
0081 
0082   // Calibration functions from global tag
0083   const PerformancePayloadFromTFormula* pfCalibrations = nullptr;
0084   const ESEEIntercalibConstants* esEEInterCalib_ = nullptr;
0085 
0086   // Barrel calibration (eta 0.00 -> 1.48)
0087   std::unique_ptr<TF1> faBarrel;
0088   std::unique_ptr<TF1> fbBarrel;
0089   std::unique_ptr<TF1> fcBarrel;
0090   std::unique_ptr<TF1> faEtaBarrelEH;
0091   std::unique_ptr<TF1> fbEtaBarrelEH;
0092   std::unique_ptr<TF1> faEtaBarrelH;
0093   std::unique_ptr<TF1> fbEtaBarrelH;
0094 
0095   // Endcap calibration (eta 1.48 -> 3.xx)
0096   std::unique_ptr<TF1> faEndcap;
0097   std::unique_ptr<TF1> fbEndcap;
0098   std::unique_ptr<TF1> fcEndcap;
0099   std::unique_ptr<TF1> faEtaEndcapEH;
0100   std::unique_ptr<TF1> fbEtaEndcapEH;
0101   std::unique_ptr<TF1> faEtaEndcapH;
0102   std::unique_ptr<TF1> fbEtaEndcapH;
0103 
0104   //added by Bhumika on 2 august 2018
0105   std::unique_ptr<TF1> fcEtaBarrelEH;
0106   std::unique_ptr<TF1> fcEtaEndcapEH;
0107   std::unique_ptr<TF1> fdEtaEndcapEH;
0108   std::unique_ptr<TF1> fcEtaBarrelH;
0109   std::unique_ptr<TF1> fcEtaEndcapH;
0110   std::unique_ptr<TF1> fdEtaEndcapH;
0111 
0112   double minimum(double a, double b) const;
0113   double dCrackPhi(double phi, double eta) const;
0114   double CorrPhi(double phi, double eta) const;
0115   double CorrEta(double eta) const;
0116   double CorrBarrel(double E, double eta) const;
0117   double Alpha(double eta) const;
0118   double Beta(double E, double eta) const;
0119   double Gamma(double etaEcal) const;
0120   double EcorrBarrel(double E, double eta, double phi, bool crackCorrection = true) const;
0121   double EcorrZoneBeforePS(double E, double eta) const;
0122   double EcorrPS(double eEcal, double ePS1, double ePS2, double etaEcal) const;
0123   double EcorrPS(double eEcal, double ePS1, double ePS2, double etaEcal, double&, double&) const;
0124   double EcorrPS_ePSNil(double eEcal, double eta) const;
0125   double EcorrZoneAfterPS(double E, double eta) const;
0126   double Ecorr(double eEcal, double ePS1, double ePS2, double eta, double phi, bool crackCorrection = true) const;
0127   double Ecorr(double eEcal,
0128                double ePS1,
0129                double ePS2,
0130                double eta,
0131                double phi,
0132                double&,
0133                double&,
0134                bool crackCorrection = true) const;
0135 
0136   // The calibration functions
0137   double aBarrel(double x) const;
0138   double bBarrel(double x) const;
0139   double cBarrel(double x) const;
0140   double aEtaBarrelEH(double x) const;
0141   double bEtaBarrelEH(double x) const;
0142   double aEtaBarrelH(double x) const;
0143   double bEtaBarrelH(double x) const;
0144   double aEndcap(double x) const;
0145   double bEndcap(double x) const;
0146   double cEndcap(double x) const;
0147   double aEtaEndcapEH(double x) const;
0148   double bEtaEndcapEH(double x) const;
0149   double aEtaEndcapH(double x) const;
0150   double bEtaEndcapH(double x) const;
0151   //added by Bhumika on 3 august 2018
0152   double cEtaBarrelEH(double x) const;
0153   double cEtaEndcapEH(double x) const;
0154   double dEtaEndcapEH(double x) const;
0155   double cEtaBarrelH(double x) const;
0156   double cEtaEndcapH(double x) const;
0157   double dEtaEndcapH(double x) const;
0158 
0159   // Threshold correction (offset)
0160   const double threshE = 3.5;
0161   const double threshH = 2.5;
0162 };
0163 
0164 #endif