CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
SCEnergyCorrectorSemiParm Class Reference

#include <SCEnergyCorrectorSemiParm.h>

Classes

class  RegParam
 

Public Member Functions

std::pair< double, double > getCorrections (const reco::SuperCluster &sc) const
 
std::vector< float > getRegData (const reco::SuperCluster &sc) const
 
void modifyObject (reco::SuperCluster &sc) const
 
 SCEnergyCorrectorSemiParm ()
 
 SCEnergyCorrectorSemiParm (const edm::ParameterSet &iConfig, edm::ConsumesCollector cc)
 
void setEvent (const edm::Event &e)
 
void setEventSetup (const edm::EventSetup &es)
 
template<edm::Transition tr = edm::Transition::BeginLuminosityBlock>
void setTokens (const edm::ParameterSet &iConfig, edm::ConsumesCollector cc)
 

Static Public Member Functions

static void fillPSetDescription (edm::ParameterSetDescription &desc)
 
static edm::ParameterSetDescription makePSetDescription ()
 

Private Member Functions

std::vector< float > getRegDataECALHLTV1 (const reco::SuperCluster &sc) const
 
std::vector< float > getRegDataECALV1 (const reco::SuperCluster &sc) const
 
std::vector< float > getRegDataHGCALHLTV1 (const reco::SuperCluster &sc) const
 
std::vector< float > getRegDataHGCALV1 (const reco::SuperCluster &sc) const
 
const RegParamgetRegParam (const DetId &detId) const
 

Private Attributes

bool applySigmaIetaIphiBug_
 
const CaloGeometrycaloGeom_
 
edm::ESGetToken< CaloGeometry, CaloGeometryRecordcaloGeomToken_
 
const CaloTopologycaloTopo_
 
edm::ESGetToken< CaloTopology, CaloTopologyRecordcaloTopoToken_
 
float hgcalCylinderR_
 
HGCalShowerShapeHelper hgcalShowerShapes_
 
float hitsEnergyThreshold_
 
bool isHLT_
 
bool isPhaseII_
 
int nHitsAboveThresholdEB_
 
int nHitsAboveThresholdEE_
 
int nHitsAboveThresholdHG_
 
edm::Handle< EcalRecHitCollectionrecHitsEB_
 
edm::Handle< EcalRecHitCollectionrecHitsEE_
 
edm::Handle< reco::PFRecHitCollectionrecHitsHgcal_
 
RegParam regParamBarrel_
 
RegParam regParamEndcap_
 
bool regTrainedWithPS_
 
edm::EDGetTokenT< EcalRecHitCollectiontokenEBRecHits_
 
edm::EDGetTokenT< EcalRecHitCollectiontokenEERecHits_
 
edm::EDGetTokenT< reco::PFRecHitCollectiontokenHgcalRecHits_
 
edm::EDGetTokenT< reco::VertexCollectiontokenVertices_
 
edm::Handle< reco::VertexCollectionvertices_
 

Detailed Description

Definition at line 37 of file SCEnergyCorrectorSemiParm.h.

Constructor & Destructor Documentation

◆ SCEnergyCorrectorSemiParm() [1/2]

SCEnergyCorrectorSemiParm::SCEnergyCorrectorSemiParm ( )

Definition at line 49 of file SCEnergyCorrectorSemiParm.cc.

50  : caloTopo_(nullptr),
51  caloGeom_(nullptr),
52  isHLT_(false),
53  isPhaseII_(false),
54  regTrainedWithPS_(true),
60  hgcalCylinderR_(0.) {}

◆ SCEnergyCorrectorSemiParm() [2/2]

SCEnergyCorrectorSemiParm::SCEnergyCorrectorSemiParm ( const edm::ParameterSet iConfig,
edm::ConsumesCollector  cc 
)

Definition at line 62 of file SCEnergyCorrectorSemiParm.cc.

References gpuPixelDoublets::cc, and setTokens().

64  setTokens(iConfig, cc);
65 }
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
void setTokens(const edm::ParameterSet &iConfig, edm::ConsumesCollector cc)

Member Function Documentation

◆ fillPSetDescription()

void SCEnergyCorrectorSemiParm::fillPSetDescription ( edm::ParameterSetDescription desc)
static

Definition at line 67 of file SCEnergyCorrectorSemiParm.cc.

References submitPVResolutionJobs::desc, ProducerED_cfi::InputTag, EgammaHGCALIDParamDefaults::kRCylinder, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by makePSetDescription().

67  {
68  desc.add<bool>("isHLT", false);
69  desc.add<bool>("isPhaseII", false);
70  desc.add<bool>("regTrainedWithPS", true);
71  desc.add<bool>("applySigmaIetaIphiBug", false);
72  desc.add<edm::InputTag>("ecalRecHitsEE", edm::InputTag("ecalRecHit", "EcalRecHitsEE"));
73  desc.add<edm::InputTag>("ecalRecHitsEB", edm::InputTag("ecalRecHit", "EcalRecHitsEB"));
74  desc.add<std::string>("regressionKeyEB", "pfscecal_EBCorrection_offline_v2");
75  desc.add<std::string>("regressionKeyEE", "pfscecal_EECorrection_offline_v2");
76  desc.add<std::string>("uncertaintyKeyEB", "pfscecal_EBUncertainty_offline_v2");
77  desc.add<std::string>("uncertaintyKeyEE", "pfscecal_EEUncertainty_offline_v2");
78  desc.add<double>("regressionMinEB", 0.2);
79  desc.add<double>("regressionMaxEB", 2.0);
80  desc.add<double>("regressionMinEE", 0.2);
81  desc.add<double>("regressionMaxEE", 2.0);
82  desc.add<double>("uncertaintyMinEB", 0.0002);
83  desc.add<double>("uncertaintyMaxEB", 0.5);
84  desc.add<double>("uncertaintyMinEE", 0.0002);
85  desc.add<double>("uncertaintyMaxEE", 0.5);
86  desc.add<edm::InputTag>("vertexCollection", edm::InputTag("offlinePrimaryVertices"));
87  desc.add<double>("eRecHitThreshold", 1.);
88  desc.add<edm::InputTag>("hgcalRecHits", edm::InputTag());
89  desc.add<double>("hgcalCylinderR", EgammaHGCALIDParamDefaults::kRCylinder);
90 }

◆ getCorrections()

std::pair< double, double > SCEnergyCorrectorSemiParm::getCorrections ( const reco::SuperCluster sc) const

Definition at line 130 of file SCEnergyCorrectorSemiParm.cc.

References DetId::det(), DetId::Ecal, EcalEndcap, getRegData(), getRegParam(), isHLT_, SiStripPI::mean, reco::SuperCluster::preshowerEnergy(), reco::SuperCluster::rawEnergy(), regTrainedWithPS_, reco::SuperCluster::seed(), reco::CaloCluster::seed(), and DetId::subdetId().

Referenced by modifyObject().

130  {
131  std::pair<double, double> corrEnergyAndRes = {-1, -1};
132 
133  const auto regData = getRegData(sc);
134  if (regData.empty()) {
135  //supercluster has no valid regression, return default values
136  return corrEnergyAndRes;
137  }
138  DetId seedId = sc.seed()->seed();
139  const auto& regParam = getRegParam(seedId);
140 
141  double mean = regParam.mean(regData);
142  double sigma = regParam.sigma(regData);
143 
144  double energyCorr = mean * sc.rawEnergy();
145  if (isHLT_ && sc.seed()->seed().det() == DetId::Ecal && seedId.subdetId() == EcalEndcap && !regTrainedWithPS_) {
146  energyCorr += sc.preshowerEnergy();
147  }
148  double resolutionEst = sigma * energyCorr;
149 
150  corrEnergyAndRes.first = energyCorr;
151  corrEnergyAndRes.second = resolutionEst;
152 
153  return corrEnergyAndRes;
154 }
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:219
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:60
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
const RegParam & getRegParam(const DetId &detId) const
Definition: DetId.h:17
std::vector< float > getRegData(const reco::SuperCluster &sc) const
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:79
double preshowerEnergy() const
energy deposited in preshower
Definition: SuperCluster.h:63

◆ getRegData()

std::vector< float > SCEnergyCorrectorSemiParm::getRegData ( const reco::SuperCluster sc) const

Definition at line 167 of file SCEnergyCorrectorSemiParm.cc.

References DetId::det(), DetId::Ecal, EcalEndcap, Exception, getRegDataECALHLTV1(), getRegDataECALV1(), getRegDataHGCALHLTV1(), getRegDataHGCALV1(), DetId::HGCalEE, isHLT_, isPhaseII_, reco::SuperCluster::seed(), reco::CaloCluster::seed(), and DetId::subdetId().

Referenced by getCorrections(), and SCEnergyCorrectorProducer::produce().

167  {
168  switch (sc.seed()->seed().det()) {
169  case DetId::Ecal:
170  if (isPhaseII_ && sc.seed()->seed().subdetId() == EcalEndcap) {
171  throw cms::Exception("ConfigError") << " Error in SCEnergyCorrectorSemiParm: "
172  << " running over events with EcalEndcap clusters while enabling "
173  "isPhaseII, please set isPhaseII = False in regression config";
174  }
175  return isHLT_ ? getRegDataECALHLTV1(sc) : getRegDataECALV1(sc);
176  case DetId::HGCalEE:
177  if (!isPhaseII_) {
178  throw cms::Exception("ConfigError") << " Error in SCEnergyCorrectorSemiParm: "
179  << " running over PhaseII events without enabling isPhaseII, please set "
180  "isPhaseII = True in regression config";
181  }
182  return isHLT_ ? getRegDataHGCALHLTV1(sc) : getRegDataHGCALV1(sc);
183  default:
184  return std::vector<float>();
185  }
186 }
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:219
std::vector< float > getRegDataECALV1(const reco::SuperCluster &sc) const
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
std::vector< float > getRegDataECALHLTV1(const reco::SuperCluster &sc) const
std::vector< float > getRegDataHGCALHLTV1(const reco::SuperCluster &sc) const
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
std::vector< float > getRegDataHGCALV1(const reco::SuperCluster &sc) const
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:79

◆ getRegDataECALHLTV1()

std::vector< float > SCEnergyCorrectorSemiParm::getRegDataECALHLTV1 ( const reco::SuperCluster sc) const
private

Definition at line 317 of file SCEnergyCorrectorSemiParm.cc.

References caloTopo_, reco::SuperCluster::clusters(), EcalBarrel, reco::CaloCluster::eta(), reco::CaloCluster::hitsAndFractions(), SiStripPI::max, nHitsAboveThresholdEB_, nHitsAboveThresholdEE_, reco::SuperCluster::phiWidth(), reco::SuperCluster::rawEnergy(), FastTrackerRecHitMaskProducer_cfi::recHits, recHitsEB_, recHitsEE_, reco::SuperCluster::seed(), and edm::PtrVectorBase::size().

Referenced by getRegData().

317  {
318  std::vector<float> eval(7, 0.);
319  auto maxDRNonSeedClus = getMaxDRNonSeedCluster(sc);
320  const float clusterMaxDR = maxDRNonSeedClus.first ? maxDRNonSeedClus.second : 999.;
321 
322  const reco::CaloCluster& seedCluster = *(sc.seed());
323  const bool iseb = seedCluster.hitsAndFractions()[0].first.subdetId() == EcalBarrel;
324  const EcalRecHitCollection* recHits = iseb ? recHitsEB_.product() : recHitsEE_.product();
325 
327  eval[1] = sc.eta();
328  eval[2] = sc.phiWidth();
329  eval[3] = EcalClusterTools::e3x3(seedCluster, recHits, caloTopo_) / sc.rawEnergy();
330  eval[4] = std::max(0, static_cast<int>(sc.clusters().size()) - 1);
331  eval[5] = clusterMaxDR;
332  eval[6] = sc.rawEnergy();
333 
334  return eval;
335 }
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:75
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:60
double phiWidth() const
obtain phi and eta width of the Super Cluster
Definition: SuperCluster.h:68
edm::Handle< EcalRecHitCollection > recHitsEE_
edm::Handle< EcalRecHitCollection > recHitsEB_
const CaloClusterPtrVector & clusters() const
const access to the cluster list itself
Definition: SuperCluster.h:82
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:79
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181

◆ getRegDataECALV1()

std::vector< float > SCEnergyCorrectorSemiParm::getRegDataECALV1 ( const reco::SuperCluster sc) const
private

Definition at line 201 of file SCEnergyCorrectorSemiParm.cc.

References applySigmaIetaIphiBug_, caloTopo_, reco::SuperCluster::clusters(), reco::SuperCluster::clustersBegin(), reco::SuperCluster::clustersEnd(), reco::deltaPhi(), reco::deltaR(), EcalBarrel, cosmicPhotonAnalyzer_cfi::eMax, reco::CaloCluster::energy(), reco::CaloCluster::eta(), reco::SuperCluster::etaWidth(), reco::CaloCluster::hitsAndFractions(), EBDetId::ieta(), edm::isNotFinite(), EEDetId::ix(), SiStripPI::max, L1TPhase2MuonOffline_cfi::maxDR, reco::CaloCluster::phi(), reco::SuperCluster::phiWidth(), reco::SuperCluster::rawEnergy(), FastTrackerRecHitMaskProducer_cfi::recHits, recHitsEB_, recHitsEE_, reco::SuperCluster::seed(), reco::CaloCluster::seed(), edm::PtrVectorBase::size(), mathSSE::sqrt(), and vertices_.

Referenced by getRegData().

201  {
202  std::vector<float> eval(30, 0.);
203 
204  const reco::CaloCluster& seedCluster = *(sc.seed());
205  const bool iseb = seedCluster.hitsAndFractions()[0].first.subdetId() == EcalBarrel;
206  const EcalRecHitCollection* recHits = iseb ? recHitsEB_.product() : recHitsEE_.product();
207 
208  const double raw_energy = sc.rawEnergy();
209  const int numberOfClusters = sc.clusters().size();
210 
211  const auto& localCovariances = EcalClusterTools::localCovariances(seedCluster, recHits, caloTopo_);
212 
213  const float eLeft = EcalClusterTools::eLeft(seedCluster, recHits, caloTopo_);
214  const float eRight = EcalClusterTools::eRight(seedCluster, recHits, caloTopo_);
215  const float eTop = EcalClusterTools::eTop(seedCluster, recHits, caloTopo_);
216  const float eBottom = EcalClusterTools::eBottom(seedCluster, recHits, caloTopo_);
217 
218  float sigmaIetaIeta = sqrt(localCovariances[0]);
219  float sigmaIetaIphi = std::numeric_limits<float>::max();
220  float sigmaIphiIphi = std::numeric_limits<float>::max();
221 
222  if (!edm::isNotFinite(localCovariances[2]))
223  sigmaIphiIphi = sqrt(localCovariances[2]);
224 
225  // extra shower shapes
226  const float see_by_spp = sigmaIetaIeta * (applySigmaIetaIphiBug_ ? std::numeric_limits<float>::max() : sigmaIphiIphi);
227  if (see_by_spp > 0) {
228  sigmaIetaIphi = localCovariances[1] / see_by_spp;
229  } else if (localCovariances[1] > 0) {
230  sigmaIetaIphi = 1.f;
231  } else {
232  sigmaIetaIphi = -1.f;
233  }
234 
235  // calculate sub-cluster variables
236  std::vector<float> clusterRawEnergy;
237  clusterRawEnergy.resize(std::max(3, numberOfClusters), 0);
238  std::vector<float> clusterDEtaToSeed;
239  clusterDEtaToSeed.resize(std::max(3, numberOfClusters), 0);
240  std::vector<float> clusterDPhiToSeed;
241  clusterDPhiToSeed.resize(std::max(3, numberOfClusters), 0);
242  float clusterMaxDR = 999.;
243  float clusterMaxDRDPhi = 999.;
244  float clusterMaxDRDEta = 999.;
245  float clusterMaxDRRawEnergy = 0.;
246 
247  size_t iclus = 0;
248  float maxDR = 0;
250  const edm::Ptr<reco::CaloCluster>& theseed = sc.seed();
251  // loop over all clusters that aren't the seed
252  auto clusend = sc.clustersEnd();
253  for (auto clus = sc.clustersBegin(); clus != clusend; ++clus) {
254  pclus = *clus;
255 
256  if (theseed == pclus)
257  continue;
258  clusterRawEnergy[iclus] = pclus->energy();
259  clusterDPhiToSeed[iclus] = reco::deltaPhi(pclus->phi(), theseed->phi());
260  clusterDEtaToSeed[iclus] = pclus->eta() - theseed->eta();
261 
262  // find cluster with max dR
263  const auto the_dr = reco::deltaR(*pclus, *theseed);
264  if (the_dr > maxDR) {
265  maxDR = the_dr;
266  clusterMaxDR = maxDR;
267  clusterMaxDRDPhi = clusterDPhiToSeed[iclus];
268  clusterMaxDRDEta = clusterDEtaToSeed[iclus];
269  clusterMaxDRRawEnergy = clusterRawEnergy[iclus];
270  }
271  ++iclus;
272  }
273 
274  eval[0] = vertices_->size();
275  eval[1] = raw_energy;
276  eval[2] = sc.etaWidth();
277  eval[3] = sc.phiWidth();
278  eval[4] = EcalClusterTools::e3x3(seedCluster, recHits, caloTopo_) / raw_energy;
279  eval[5] = seedCluster.energy() / raw_energy;
280  eval[6] = EcalClusterTools::eMax(seedCluster, recHits) / raw_energy;
281  eval[7] = EcalClusterTools::e2nd(seedCluster, recHits) / raw_energy;
282  eval[8] = (eLeft + eRight != 0.f ? (eLeft - eRight) / (eLeft + eRight) : 0.f);
283  eval[9] = (eTop + eBottom != 0.f ? (eTop - eBottom) / (eTop + eBottom) : 0.f);
284  eval[10] = sigmaIetaIeta;
285  eval[11] = sigmaIetaIphi;
286  eval[12] = sigmaIphiIphi;
287  eval[13] = std::max(0, numberOfClusters - 1);
288  eval[14] = clusterMaxDR;
289  eval[15] = clusterMaxDRDPhi;
290  eval[16] = clusterMaxDRDEta;
291  eval[17] = clusterMaxDRRawEnergy / raw_energy;
292  eval[18] = clusterRawEnergy[0] / raw_energy;
293  eval[19] = clusterRawEnergy[1] / raw_energy;
294  eval[20] = clusterRawEnergy[2] / raw_energy;
295  eval[21] = clusterDPhiToSeed[0];
296  eval[22] = clusterDPhiToSeed[1];
297  eval[23] = clusterDPhiToSeed[2];
298  eval[24] = clusterDEtaToSeed[0];
299  eval[25] = clusterDEtaToSeed[1];
300  eval[26] = clusterDEtaToSeed[2];
301  if (iseb) {
302  EBDetId ebseedid(seedCluster.seed());
303  eval[27] = ebseedid.ieta();
304  eval[28] = ebseedid.iphi();
305  } else {
306  EEDetId eeseedid(seedCluster.seed());
307  eval[27] = eeseedid.ix();
308  eval[28] = eeseedid.iy();
309  //seed cluster eta is only needed for the 106X Ultra Legacy regressions
310  //and was not used in the 74X regression however as its just an extra varaible
311  //at the end, its harmless to add for the 74X regression
312  eval[29] = seedCluster.eta();
313  }
314  return eval;
315 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:75
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:219
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:60
CaloCluster_iterator clustersBegin() const
fist iterator over BasicCluster constituents
Definition: SuperCluster.h:88
int ix() const
Definition: EEDetId.h:77
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:184
double phiWidth() const
obtain phi and eta width of the Super Cluster
Definition: SuperCluster.h:68
CaloCluster_iterator clustersEnd() const
last iterator over BasicCluster constituents
Definition: SuperCluster.h:91
int ieta() const
get the crystal ieta
Definition: EBDetId.h:49
T sqrt(T t)
Definition: SSEVec.h:19
edm::Handle< EcalRecHitCollection > recHitsEE_
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
double energy() const
cluster energy
Definition: CaloCluster.h:149
double etaWidth() const
Definition: SuperCluster.h:69
edm::Handle< EcalRecHitCollection > recHitsEB_
const CaloClusterPtrVector & clusters() const
const access to the cluster list itself
Definition: SuperCluster.h:82
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:79
edm::Handle< reco::VertexCollection > vertices_
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181

◆ getRegDataHGCALHLTV1()

std::vector< float > SCEnergyCorrectorSemiParm::getRegDataHGCALHLTV1 ( const reco::SuperCluster sc) const
private

Definition at line 368 of file SCEnergyCorrectorSemiParm.cc.

References reco::SuperCluster::clusters(), HGCalShowerShapeHelper::createCalc(), reco::CaloCluster::eta(), hgcalCylinderR_, hgcalShowerShapes_, SiStripPI::max, nHitsAboveThresholdEB_, nHitsAboveThresholdHG_, reco::SuperCluster::phiWidth(), reco::SuperCluster::rawEnergy(), and edm::PtrVectorBase::size().

Referenced by getRegData().

368  {
369  std::vector<float> eval(7, 0.);
370  const float clusterMaxDR = getMaxDRNonSeedCluster(sc).second;
371 
372  auto ssCalc = hgcalShowerShapes_.createCalc(sc);
373 
374  eval[0] = sc.rawEnergy();
375  eval[1] = sc.eta();
376  eval[2] = sc.phiWidth();
377  eval[3] = std::max(0, static_cast<int>(sc.clusters().size()) - 1);
378  eval[4] = ssCalc.getRvar(hgcalCylinderR_);
379  eval[5] = clusterMaxDR;
381 
382  return eval;
383 }
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:75
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:60
double phiWidth() const
obtain phi and eta width of the Super Cluster
Definition: SuperCluster.h:68
const CaloClusterPtrVector & clusters() const
const access to the cluster list itself
Definition: SuperCluster.h:82
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181
HGCalShowerShapeHelper hgcalShowerShapes_
HGCalShowerShapeHelper::ShowerShapeCalc createCalc(const std::vector< std::pair< DetId, float > > &hitsAndFracs, double rawEnergy, double minHitE=0, double minHitET=0, int minLayer=1, int maxLayer=-1, DetId::Detector subDet=DetId::HGCalEE) const

◆ getRegDataHGCALV1()

std::vector< float > SCEnergyCorrectorSemiParm::getRegDataHGCALV1 ( const reco::SuperCluster sc) const
private

Definition at line 337 of file SCEnergyCorrectorSemiParm.cc.

References reco::SuperCluster::clusters(), HGCalShowerShapeHelper::createCalc(), reco::deltaPhi(), reco::CaloCluster::energy(), reco::CaloCluster::eta(), reco::SuperCluster::etaWidth(), HGCalShowerShapeHelper::ShowerShapeCalc::getPCAWidths(), hgcalCylinderR_, hgcalShowerShapes_, reco::CaloCluster::hitsAndFractions(), nHitsAboveThresholdEB_, nHitsAboveThresholdHG_, reco::CaloCluster::phi(), reco::SuperCluster::phiWidth(), reco::SuperCluster::rawEnergy(), reco::SuperCluster::seed(), edm::PtrVectorBase::size(), and mathSSE::sqrt().

Referenced by getRegData().

337  {
338  std::vector<float> eval(17, 0.);
339 
340  auto ssCalc = hgcalShowerShapes_.createCalc(sc);
341  auto pcaWidths = ssCalc.getPCAWidths(hgcalCylinderR_);
342  auto energyHighestHits = ssCalc.getEnergyHighestHits(2);
343 
344  auto maxDRNonSeedClus = getMaxDRNonSeedCluster(sc);
345  const float clusterMaxDR = maxDRNonSeedClus.first ? maxDRNonSeedClus.second : 999.;
346 
347  eval[0] = sc.rawEnergy();
348  eval[1] = sc.eta();
349  eval[2] = sc.etaWidth();
350  eval[3] = sc.phiWidth();
351  eval[4] = sc.clusters().size();
352  eval[5] = sc.hitsAndFractions().size();
353  eval[6] = clusterMaxDR;
354  eval[7] = sc.eta() - sc.seed()->eta();
355  eval[8] = reco::deltaPhi(sc.phi(), sc.seed()->phi());
356  eval[9] = energyHighestHits[0] / sc.rawEnergy();
357  eval[10] = energyHighestHits[1] / sc.rawEnergy();
358  eval[11] = std::sqrt(pcaWidths.sigma2uu);
359  eval[12] = std::sqrt(pcaWidths.sigma2vv);
360  eval[13] = std::sqrt(pcaWidths.sigma2ww);
361  eval[14] = ssCalc.getRvar(hgcalCylinderR_, sc.rawEnergy());
362  eval[15] = sc.seed()->energy() / sc.rawEnergy();
364 
365  return eval;
366 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:75
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:60
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:184
double phiWidth() const
obtain phi and eta width of the Super Cluster
Definition: SuperCluster.h:68
T sqrt(T t)
Definition: SSEVec.h:19
double energy() const
cluster energy
Definition: CaloCluster.h:149
double etaWidth() const
Definition: SuperCluster.h:69
const CaloClusterPtrVector & clusters() const
const access to the cluster list itself
Definition: SuperCluster.h:82
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:79
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181
HGCalShowerShapeHelper hgcalShowerShapes_
ShowerWidths getPCAWidths(double cylinderR, bool useFractions=false) const
HGCalShowerShapeHelper::ShowerShapeCalc createCalc(const std::vector< std::pair< DetId, float > > &hitsAndFracs, double rawEnergy, double minHitE=0, double minHitET=0, int minLayer=1, int maxLayer=-1, DetId::Detector subDet=DetId::HGCalEE) const

◆ getRegParam()

const RegParam& SCEnergyCorrectorSemiParm::getRegParam ( const DetId detId) const
inlineprivate

Definition at line 100 of file SCEnergyCorrectorSemiParm.h.

References DetId::det(), DetId::Ecal, EcalBarrel, regParamBarrel_, regParamEndcap_, and DetId::subdetId().

Referenced by getCorrections().

100  {
101  return detId.det() == DetId::Ecal && detId.subdetId() == EcalBarrel ? regParamBarrel_ : regParamEndcap_;
102  }
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48

◆ makePSetDescription()

edm::ParameterSetDescription SCEnergyCorrectorSemiParm::makePSetDescription ( )
static

◆ modifyObject()

void SCEnergyCorrectorSemiParm::modifyObject ( reco::SuperCluster sc) const

Definition at line 156 of file SCEnergyCorrectorSemiParm.cc.

References getCorrections(), reco::CaloCluster::setCorrectedEnergy(), reco::CaloCluster::setCorrectedEnergyUncertainty(), and reco::CaloCluster::setEnergy().

Referenced by SCEnergyCorrectorProducer::produce().

156  {
157  std::pair<double, double> cor = getCorrections(sc);
158  if (cor.first < 0)
159  return;
160  sc.setEnergy(cor.first);
161  sc.setCorrectedEnergy(cor.first);
162  if (cor.second >= 0) {
163  sc.setCorrectedEnergyUncertainty(cor.second);
164  }
165 }
void setEnergy(double energy)
Definition: CaloCluster.h:136
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:137
std::pair< double, double > getCorrections(const reco::SuperCluster &sc) const
void setCorrectedEnergyUncertainty(float energyerr)
Definition: CaloCluster.h:138

◆ setEvent()

void SCEnergyCorrectorSemiParm::setEvent ( const edm::Event e)

Definition at line 110 of file SCEnergyCorrectorSemiParm.cc.

References hgcalShowerShapes_, hitsEnergyThreshold_, HGCalShowerShapeHelper::initPerEvent(), isHLT_, isPhaseII_, nHitsAboveThresholdEB_, nHitsAboveThresholdEE_, nHitsAboveThresholdHG_, recHitsEB_, recHitsEE_, recHitsHgcal_, tokenEBRecHits_, tokenEERecHits_, tokenHgcalRecHits_, tokenVertices_, and vertices_.

Referenced by SCEnergyCorrectorProducer::produce().

110  {
111  event.getByToken(tokenEBRecHits_, recHitsEB_);
112  if (!isPhaseII_) {
113  event.getByToken(tokenEERecHits_, recHitsEE_);
114  } else {
115  event.getByToken(tokenHgcalRecHits_, recHitsHgcal_);
117  }
118  if (isHLT_ || isPhaseII_) {
119  //note countRecHits checks the validity of the handle and returns 0
120  //if invalid so its okay to call on all rec-hit collections here
124  }
125  if (!isHLT_) {
126  event.getByToken(tokenVertices_, vertices_);
127  }
128 }
edm::Handle< reco::PFRecHitCollection > recHitsHgcal_
edm::EDGetTokenT< EcalRecHitCollection > tokenEERecHits_
void initPerEvent(const std::vector< reco::PFRecHit > &recHits)
edm::Handle< EcalRecHitCollection > recHitsEE_
edm::EDGetTokenT< EcalRecHitCollection > tokenEBRecHits_
edm::Handle< EcalRecHitCollection > recHitsEB_
edm::Handle< reco::VertexCollection > vertices_
edm::EDGetTokenT< reco::PFRecHitCollection > tokenHgcalRecHits_
edm::EDGetTokenT< reco::VertexCollection > tokenVertices_
HGCalShowerShapeHelper hgcalShowerShapes_

◆ setEventSetup()

void SCEnergyCorrectorSemiParm::setEventSetup ( const edm::EventSetup es)

Definition at line 98 of file SCEnergyCorrectorSemiParm.cc.

References caloGeom_, caloGeomToken_, caloTopo_, caloTopoToken_, edm::EventSetup::getData(), hgcalShowerShapes_, HGCalShowerShapeHelper::initPerSetup(), isPhaseII_, regParamBarrel_, regParamEndcap_, and SCEnergyCorrectorSemiParm::RegParam::setForests().

Referenced by SCEnergyCorrectorProducer::beginLuminosityBlock().

98  {
101 
104 
105  if (isPhaseII_) {
107  }
108 }
edm::ESGetToken< CaloTopology, CaloTopologyRecord > caloTopoToken_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
void setForests(const edm::EventSetup &setup)
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
void initPerSetup(const edm::EventSetup &iSetup)
HGCalShowerShapeHelper hgcalShowerShapes_

◆ setTokens()

template<edm::Transition esTransition>
void SCEnergyCorrectorSemiParm::setTokens ( const edm::ParameterSet iConfig,
edm::ConsumesCollector  cc 
)

Definition at line 147 of file SCEnergyCorrectorSemiParm.h.

References applySigmaIetaIphiBug_, caloGeomToken_, caloTopoToken_, gpuPixelDoublets::cc, edm::ParameterSet::getParameter(), hgcalCylinderR_, hgcalShowerShapes_, hitsEnergyThreshold_, isHLT_, isPhaseII_, regParamBarrel_, regParamEndcap_, regTrainedWithPS_, SCEnergyCorrectorSemiParm::RegParam::setTokens(), HGCalShowerShapeHelper::setTokens(), AlCaHLTBitMon_QueryRunRegistry::string, tokenEBRecHits_, tokenEERecHits_, tokenHgcalRecHits_, and tokenVertices_.

Referenced by SCEnergyCorrectorSemiParm().

147  {
148  isHLT_ = iConfig.getParameter<bool>("isHLT");
149  isPhaseII_ = iConfig.getParameter<bool>("isPhaseII");
150  regTrainedWithPS_ = iConfig.getParameter<bool>("regTrainedWithPS");
151  applySigmaIetaIphiBug_ = iConfig.getParameter<bool>("applySigmaIetaIphiBug");
152  tokenEBRecHits_ = cc.consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("ecalRecHitsEB"));
153  if (not isPhaseII_) {
154  tokenEERecHits_ = cc.consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("ecalRecHitsEE"));
155  } else {
156  tokenHgcalRecHits_ = cc.consumes<reco::PFRecHitCollection>(iConfig.getParameter<edm::InputTag>("hgcalRecHits"));
157  hgcalCylinderR_ = iConfig.getParameter<double>("hgcalCylinderR");
158  hgcalShowerShapes_.setTokens<esTransition>(cc);
159  }
160  caloGeomToken_ = cc.esConsumes<CaloGeometry, CaloGeometryRecord, esTransition>();
161  caloTopoToken_ = cc.esConsumes<CaloTopology, CaloTopologyRecord, esTransition>();
162 
163  regParamBarrel_ = RegParam(iConfig.getParameter<std::string>("regressionKeyEB"),
164  iConfig.getParameter<double>("regressionMinEB"),
165  iConfig.getParameter<double>("regressionMaxEB"),
166  iConfig.getParameter<std::string>("uncertaintyKeyEB"),
167  iConfig.getParameter<double>("uncertaintyMinEB"),
168  iConfig.getParameter<double>("uncertaintyMaxEB"));
169  regParamBarrel_.setTokens<esTransition>(cc);
170  regParamEndcap_ = RegParam(iConfig.getParameter<std::string>("regressionKeyEE"),
171  iConfig.getParameter<double>("regressionMinEE"),
172  iConfig.getParameter<double>("regressionMaxEE"),
173  iConfig.getParameter<std::string>("uncertaintyKeyEE"),
174  iConfig.getParameter<double>("uncertaintyMinEE"),
175  iConfig.getParameter<double>("uncertaintyMaxEE"));
176  regParamEndcap_.setTokens<esTransition>(cc);
177  hitsEnergyThreshold_ = iConfig.getParameter<double>("eRecHitThreshold");
178  if (not isHLT_) {
179  tokenVertices_ = cc.consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexCollection"));
180  }
181 }
edm::ESGetToken< CaloTopology, CaloTopologyRecord > caloTopoToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< PFRecHit > PFRecHitCollection
collection of PFRecHit objects
Definition: PFRecHitFwd.h:9
edm::EDGetTokenT< EcalRecHitCollection > tokenEERecHits_
void setTokens(edm::ConsumesCollector cc)
edm::EDGetTokenT< EcalRecHitCollection > tokenEBRecHits_
void setTokens(edm::ConsumesCollector consumesCollector)
edm::EDGetTokenT< reco::PFRecHitCollection > tokenHgcalRecHits_
edm::EDGetTokenT< reco::VertexCollection > tokenVertices_
HGCalShowerShapeHelper hgcalShowerShapes_

Member Data Documentation

◆ applySigmaIetaIphiBug_

bool SCEnergyCorrectorSemiParm::applySigmaIetaIphiBug_
private

Definition at line 131 of file SCEnergyCorrectorSemiParm.h.

Referenced by getRegDataECALV1(), and setTokens().

◆ caloGeom_

const CaloGeometry* SCEnergyCorrectorSemiParm::caloGeom_
private

Definition at line 114 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup().

◆ caloGeomToken_

edm::ESGetToken<CaloGeometry, CaloGeometryRecord> SCEnergyCorrectorSemiParm::caloGeomToken_
private

Definition at line 116 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup(), and setTokens().

◆ caloTopo_

const CaloTopology* SCEnergyCorrectorSemiParm::caloTopo_
private

◆ caloTopoToken_

edm::ESGetToken<CaloTopology, CaloTopologyRecord> SCEnergyCorrectorSemiParm::caloTopoToken_
private

Definition at line 115 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup(), and setTokens().

◆ hgcalCylinderR_

float SCEnergyCorrectorSemiParm::hgcalCylinderR_
private

Definition at line 136 of file SCEnergyCorrectorSemiParm.h.

Referenced by getRegDataHGCALHLTV1(), getRegDataHGCALV1(), and setTokens().

◆ hgcalShowerShapes_

HGCalShowerShapeHelper SCEnergyCorrectorSemiParm::hgcalShowerShapes_
private

◆ hitsEnergyThreshold_

float SCEnergyCorrectorSemiParm::hitsEnergyThreshold_
private

Definition at line 135 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

◆ isHLT_

bool SCEnergyCorrectorSemiParm::isHLT_
private

Definition at line 128 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), getRegData(), setEvent(), and setTokens().

◆ isPhaseII_

bool SCEnergyCorrectorSemiParm::isPhaseII_
private

Definition at line 129 of file SCEnergyCorrectorSemiParm.h.

Referenced by getRegData(), setEvent(), setEventSetup(), and setTokens().

◆ nHitsAboveThresholdEB_

int SCEnergyCorrectorSemiParm::nHitsAboveThresholdEB_
private

◆ nHitsAboveThresholdEE_

int SCEnergyCorrectorSemiParm::nHitsAboveThresholdEE_
private

Definition at line 133 of file SCEnergyCorrectorSemiParm.h.

Referenced by getRegDataECALHLTV1(), and setEvent().

◆ nHitsAboveThresholdHG_

int SCEnergyCorrectorSemiParm::nHitsAboveThresholdHG_
private

Definition at line 134 of file SCEnergyCorrectorSemiParm.h.

Referenced by getRegDataHGCALHLTV1(), getRegDataHGCALV1(), and setEvent().

◆ recHitsEB_

edm::Handle<EcalRecHitCollection> SCEnergyCorrectorSemiParm::recHitsEB_
private

Definition at line 123 of file SCEnergyCorrectorSemiParm.h.

Referenced by getRegDataECALHLTV1(), getRegDataECALV1(), and setEvent().

◆ recHitsEE_

edm::Handle<EcalRecHitCollection> SCEnergyCorrectorSemiParm::recHitsEE_
private

Definition at line 124 of file SCEnergyCorrectorSemiParm.h.

Referenced by getRegDataECALHLTV1(), getRegDataECALV1(), and setEvent().

◆ recHitsHgcal_

edm::Handle<reco::PFRecHitCollection> SCEnergyCorrectorSemiParm::recHitsHgcal_
private

Definition at line 125 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent().

◆ regParamBarrel_

RegParam SCEnergyCorrectorSemiParm::regParamBarrel_
private

Definition at line 110 of file SCEnergyCorrectorSemiParm.h.

Referenced by getRegParam(), setEventSetup(), and setTokens().

◆ regParamEndcap_

RegParam SCEnergyCorrectorSemiParm::regParamEndcap_
private

Definition at line 111 of file SCEnergyCorrectorSemiParm.h.

Referenced by getRegParam(), setEventSetup(), and setTokens().

◆ regTrainedWithPS_

bool SCEnergyCorrectorSemiParm::regTrainedWithPS_
private

Definition at line 130 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setTokens().

◆ tokenEBRecHits_

edm::EDGetTokenT<EcalRecHitCollection> SCEnergyCorrectorSemiParm::tokenEBRecHits_
private

Definition at line 118 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

◆ tokenEERecHits_

edm::EDGetTokenT<EcalRecHitCollection> SCEnergyCorrectorSemiParm::tokenEERecHits_
private

Definition at line 119 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

◆ tokenHgcalRecHits_

edm::EDGetTokenT<reco::PFRecHitCollection> SCEnergyCorrectorSemiParm::tokenHgcalRecHits_
private

Definition at line 120 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

◆ tokenVertices_

edm::EDGetTokenT<reco::VertexCollection> SCEnergyCorrectorSemiParm::tokenVertices_
private

Definition at line 121 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

◆ vertices_

edm::Handle<reco::VertexCollection> SCEnergyCorrectorSemiParm::vertices_
private

Definition at line 126 of file SCEnergyCorrectorSemiParm.h.

Referenced by getRegDataECALV1(), and setEvent().