CMS 3D CMS Logo

PhotonEnergyCalibrator.cc
Go to the documentation of this file.
2 
6 #include <CLHEP/Random/RandGaussQ.h>
7 
10 
12  correctionRetriever_(correctionFile),
13  rng_(nullptr),
14  minEt_(1.0)
15 {
16 
17 }
18 
20 {
21  rng_ = rnd;
22 }
23 
24 std::array<float,EGEnergySysIndex::kNrSysErrs> PhotonEnergyCalibrator::
25 calibrate(reco::Photon &photon,const unsigned int runNumber,
26  const EcalRecHitCollection *recHits,
27  edm::StreamID const & id,
28  const PhotonEnergyCalibrator::EventType eventType) const
29 {
30  return calibrate(photon,runNumber,recHits,gauss(id),eventType);
31 }
32 
33 std::array<float,EGEnergySysIndex::kNrSysErrs> PhotonEnergyCalibrator::
34 calibrate(reco::Photon &photon,const unsigned int runNumber,
35  const EcalRecHitCollection *recHits,
36  const float smearNrSigma,
37  const PhotonEnergyCalibrator::EventType eventType) const
38 {
39  const float scEtaAbs = std::abs(photon.superCluster()->eta());
40  const float et = photon.getCorrectedEnergy(reco::Photon::P4type::regression2) / cosh(scEtaAbs);
41 
42  if (et < minEt_) {
43  std::array<float,EGEnergySysIndex::kNrSysErrs> retVal;
44  retVal.fill(photon.getCorrectedEnergy(reco::Photon::P4type::regression2));
45  retVal[EGEnergySysIndex::kScaleValue] = 1.0;
46  retVal[EGEnergySysIndex::kSmearValue] = 0.0;
47  retVal[EGEnergySysIndex::kSmearNrSigma] = smearNrSigma;
48  retVal[EGEnergySysIndex::kEcalErrPreCorr] = photon.getCorrectedEnergyError(reco::Photon::P4type::regression2);
49  retVal[EGEnergySysIndex::kEcalErrPostCorr] = photon.getCorrectedEnergyError(reco::Photon::P4type::regression2);
54 
55  return retVal;
56  }
57 
58  const DetId seedDetId = photon.superCluster()->seed()->seed();
59  EcalRecHitCollection::const_iterator seedRecHit = recHits->find(seedDetId);
60  unsigned int gainSeedSC = 12;
61  if (seedRecHit != recHits->end()) {
62  if(seedRecHit->checkFlag(EcalRecHit::kHasSwitchToGain6)) gainSeedSC = 6;
63  if(seedRecHit->checkFlag(EcalRecHit::kHasSwitchToGain1)) gainSeedSC = 1;
64  }
65 
66  const EnergyScaleCorrection::ScaleCorrection* scaleCorr = correctionRetriever_.getScaleCorr(runNumber, et, scEtaAbs, photon.full5x5_r9(), gainSeedSC);
67  const EnergyScaleCorrection::SmearCorrection* smearCorr = correctionRetriever_.getSmearCorr(runNumber, et, scEtaAbs, photon.full5x5_r9(), gainSeedSC);
68  if(scaleCorr==nullptr) scaleCorr=&defaultScaleCorr_;
69  if(smearCorr==nullptr) smearCorr=&defaultSmearCorr_;
70 
71 
72  std::array<float,EGEnergySysIndex::kNrSysErrs> uncertainties{};
73 
74  uncertainties[EGEnergySysIndex::kScaleValue] = scaleCorr->scale();
75  uncertainties[EGEnergySysIndex::kSmearValue] = smearCorr->sigma(et); //even though we use scale = 1.0, we still store the value returned for MC
76  uncertainties[EGEnergySysIndex::kSmearNrSigma] = smearNrSigma;
77  //MC central values are not scaled (scale = 1.0), data is not smeared (smearNrSigma = 0)
78  //smearing still has a second order effect on data as it enters the E/p combination as an
79  //extra uncertainty on the calo energy
80  //MC gets all the scale systematics
81  if(eventType == EventType::DATA){
82  setEnergyAndSystVarations(scaleCorr->scale(),0.,et,*scaleCorr,*smearCorr,photon,uncertainties);
83  }else if(eventType == EventType::MC){
84  setEnergyAndSystVarations(1.0,smearNrSigma,et,*scaleCorr,*smearCorr,photon,uncertainties);
85  }
86 
87  return uncertainties;
88 
89 }
90 
92 setEnergyAndSystVarations(const float scale,const float smearNrSigma,const float et,
96  std::array<float,EGEnergySysIndex::kNrSysErrs>& energyData)const
97 {
98 
99  const float smear = smearCorr.sigma(et);
100  const float smearRhoUp = smearCorr.sigma(et,1,0);
101  const float smearRhoDn = smearCorr.sigma(et,-1,0);
102  const float smearPhiUp = smearCorr.sigma(et,0,1);
103  const float smearPhiDn = smearCorr.sigma(et,0,-1);
104 
105  const float corr = scale + smear * smearNrSigma;
106  const float corrRhoUp = scale + smearRhoUp * smearNrSigma;
107  const float corrRhoDn = scale + smearRhoDn * smearNrSigma;
108  const float corrPhiUp = scale + smearPhiUp * smearNrSigma;
109  const float corrPhiDn = scale + smearPhiDn * smearNrSigma;
110  const float corrUp = corrRhoUp;
111  const float corrDn = corrRhoDn;
112 
113 
114  const double oldEcalEnergy = photon.getCorrectedEnergy(reco::Photon::P4type::regression2);
115  const double oldEcalEnergyError = photon.getCorrectedEnergyError(reco::Photon::P4type::regression2);
116 
117  energyData[EGEnergySysIndex::kEcalPreCorr] = oldEcalEnergy;
118  energyData[EGEnergySysIndex::kEcalErrPreCorr] = oldEcalEnergyError;
119 
120  const double newEcalEnergy = oldEcalEnergy * corr;
121  const double newEcalEnergyError = std::hypot(oldEcalEnergyError * corr, smear * newEcalEnergy);
122  photon.setCorrectedEnergy(reco::Photon::P4type::regression2, newEcalEnergy, newEcalEnergyError, true);
123 
124  energyData[EGEnergySysIndex::kScaleStatUp] = oldEcalEnergy * (corr + scaleCorr.scaleErrStat());
125  energyData[EGEnergySysIndex::kScaleStatDown] = oldEcalEnergy * (corr - scaleCorr.scaleErrStat());
126  energyData[EGEnergySysIndex::kScaleSystUp] = oldEcalEnergy * (corr + scaleCorr.scaleErrSyst());
127  energyData[EGEnergySysIndex::kScaleSystDown] = oldEcalEnergy * (corr - scaleCorr.scaleErrSyst());
128  energyData[EGEnergySysIndex::kScaleGainUp] = oldEcalEnergy * (corr + scaleCorr.scaleErrGain());
129  energyData[EGEnergySysIndex::kScaleGainDown] = oldEcalEnergy * (corr - scaleCorr.scaleErrGain());
130  energyData[EGEnergySysIndex::kSmearRhoUp] = oldEcalEnergy * corrRhoUp;
131  energyData[EGEnergySysIndex::kSmearRhoDown] = oldEcalEnergy * corrRhoDn;
132  energyData[EGEnergySysIndex::kSmearPhiUp] = oldEcalEnergy * corrPhiUp;
133  energyData[EGEnergySysIndex::kSmearPhiDown] = oldEcalEnergy * corrPhiDn;
134 
135  // The total variation
136  energyData[EGEnergySysIndex::kScaleUp] = oldEcalEnergy * (corr + scaleCorr.scaleErr(EnergyScaleCorrection::kErrStatSystGain));
137  energyData[EGEnergySysIndex::kScaleDown] = oldEcalEnergy * (corr - scaleCorr.scaleErr(EnergyScaleCorrection::kErrStatSystGain));
138  energyData[EGEnergySysIndex::kSmearUp] = oldEcalEnergy * corrUp;
139  energyData[EGEnergySysIndex::kSmearDown] = oldEcalEnergy * corrDn;
140 
141 
142  energyData[EGEnergySysIndex::kEcalPostCorr] = photon.getCorrectedEnergy(reco::Photon::P4type::regression2);
143  energyData[EGEnergySysIndex::kEcalErrPostCorr] = photon.getCorrectedEnergyError(reco::Photon::P4type::regression2);
144 }
145 
146 
148 {
149  if (rng_) {
150  return rng_->Gaus();
151  } else {
153  if ( !rng.isAvailable() ) {
154  throw cms::Exception("Configuration")
155  << "XXXXXXX requires the RandomNumberGeneratorService\n"
156  "which is not present in the configuration file. You must add the service\n"
157  "in the configuration file or remove the modules that require it.";
158  }
159  CLHEP::RandGaussQ gaussDistribution(rng->getEngine(id), 0.0, 1.0);
160  return gaussDistribution.fire();
161  }
162 }
163 
EnergyScaleCorrection correctionRetriever_
void initPrivateRng(TRandom *rnd)
void setCorrectedEnergy(P4type type, float E, float dE, bool toCand=true)
std::vector< EcalRecHit >::const_iterator const_iterator
#define nullptr
reco::SuperClusterRef superCluster() const override
Ref to SuperCluster.
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
float full5x5_r9() const
Definition: Photon.h:240
const ScaleCorrection * getScaleCorr(unsigned int runnr, double et, double eta, double r9, unsigned int gainSeed) const
bool isAvailable() const
Definition: Service.h:46
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const SmearCorrection * getSmearCorr(unsigned int runnr, double et, double eta, double r9, unsigned int gainSeed) const
static const EnergyScaleCorrection::ScaleCorrection defaultScaleCorr_
static const EnergyScaleCorrection::SmearCorrection defaultSmearCorr_
float getCorrectedEnergyError(P4type type) const
JetCorrectorParameters corr
Definition: classes.h:5
float scaleErr(const std::bitset< kErrNrBits > &uncBitMask) const
float sigma(const float et, const float nrSigmaRho=0., const float nrSigmaPhi=0.) const
const_iterator end() const
Definition: DetId.h:18
et
define resolution functions of each parameter
float getCorrectedEnergy(P4type type) const
double gauss(edm::StreamID const &id) const
iterator find(key_type k)
void setEnergyAndSystVarations(const float scale, const float smearNrSigma, const float et, const EnergyScaleCorrection::ScaleCorrection &scaleCorr, const EnergyScaleCorrection::SmearCorrection &smearCorr, reco::Photon &photon, std::array< float, EGEnergySysIndex::kNrSysErrs > &energyData) const
std::array< float, EGEnergySysIndex::kNrSysErrs > calibrate(reco::Photon &photon, const unsigned int runNumber, const EcalRecHitCollection *recHits, edm::StreamID const &id, const EventType eventType) const