CMS 3D CMS Logo

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