CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Attributes | Private Attributes
SCEnergyCorrectorSemiParm Class Reference

#include <SCEnergyCorrectorSemiParm.h>

Public Member Functions

std::pair< double, double > getCorrections (const reco::SuperCluster &sc) const
 
void modifyObject (reco::SuperCluster &sc)
 
 SCEnergyCorrectorSemiParm ()
 
void setEvent (const edm::Event &e)
 
void setEventSetup (const edm::EventSetup &es)
 
void setTokens (const edm::ParameterSet &iConfig, edm::ConsumesCollector &cc)
 
 ~SCEnergyCorrectorSemiParm ()
 

Protected Attributes

edm::ESHandle< CaloGeometrycalogeom_
 
edm::ESHandle< CaloTopologycalotopo_
 
edm::InputTag ecalHitsEBInputTag_
 
edm::InputTag ecalHitsEEInputTag_
 
const GBRForestDforesteb_
 
const GBRForestDforestee_
 
const GBRForestDforestsigmaeb_
 
const GBRForestDforestsigmaee_
 
edm::Handle< EcalRecHitCollectionrechitsEB_
 
edm::Handle< EcalRecHitCollectionrechitsEE_
 
std::string regressionKeyEB_
 
std::string regressionKeyEE_
 
edm::EDGetTokenT
< EcalRecHitCollection
tokenEBRecHits_
 
edm::EDGetTokenT
< EcalRecHitCollection
tokenEERecHits_
 
edm::EDGetTokenT
< reco::VertexCollection
tokenVertices_
 
std::string uncertaintyKeyEB_
 
std::string uncertaintyKeyEE_
 
edm::InputTag vertexInputTag_
 
edm::Handle
< reco::VertexCollection
vertices_
 

Private Attributes

float eThreshold_
 
bool isHLT_
 
int nHitsAboveThreshold_
 

Detailed Description

Definition at line 30 of file SCEnergyCorrectorSemiParm.h.

Constructor & Destructor Documentation

SCEnergyCorrectorSemiParm::SCEnergyCorrectorSemiParm ( )

Definition at line 16 of file SCEnergyCorrectorSemiParm.cc.

16  :
17 foresteb_(0),
18 forestee_(0),
21 calotopo_(0),
22 calogeom_(0)
23 {}
edm::ESHandle< CaloTopology > calotopo_
edm::ESHandle< CaloGeometry > calogeom_
SCEnergyCorrectorSemiParm::~SCEnergyCorrectorSemiParm ( )

Definition at line 26 of file SCEnergyCorrectorSemiParm.cc.

27 {}

Member Function Documentation

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

Definition at line 96 of file SCEnergyCorrectorSemiParm.cc.

References calotopo_, reco::SuperCluster::clusters(), reco::SuperCluster::clustersBegin(), reco::SuperCluster::clustersEnd(), constexpr, reco::deltaPhi(), reco::deltaR(), EcalBarrel, reco::CaloCluster::energy(), reco::CaloCluster::eta(), reco::SuperCluster::etaWidth(), foresteb_, forestee_, forestsigmaeb_, forestsigmaee_, DetId::Forward, GBRForestD::GetResponse(), reco::CaloCluster::hitsAndFractions(), EBDetId::ieta(), isHLT_, edm::isNotFinite(), EEDetId::ix(), bookConverter::max, timingPdfMaker::mean, nHitsAboveThreshold_, AlCaHLTBitMon_ParallelJobs::p, reco::SuperCluster::phiWidth(), reco::SuperCluster::preshowerEnergy(), edm::ESHandle< class >::product(), reco::SuperCluster::rawEnergy(), HLT_25ns10e33_v2_cff::recHits, rechitsEB_, rechitsEE_, reco::SuperCluster::seed(), reco::CaloCluster::seed(), edm::PtrVectorBase::size(), mathSSE::sqrt(), and vertices_.

Referenced by modifyObject().

96  {
97  std::pair<double, double> p;
98  p.first=-1;
99  p.second=-1;
100 
101  // protect against HGCal, don't mod the object
102  if( sc.seed()->seed().det() == DetId::Forward ) return p;
103 
104  const reco::CaloCluster &seedCluster = *(sc.seed());
105  const bool iseb = seedCluster.hitsAndFractions()[0].first.subdetId() == EcalBarrel;
106  const EcalRecHitCollection *recHits = iseb ? rechitsEB_.product() : rechitsEE_.product();
107 
108  const CaloTopology *topo = calotopo_.product();
109 
110  const double raw_energy = sc.rawEnergy();
111  const int numberOfClusters = sc.clusters().size();
112 
113  std::vector<float> localCovariances = EcalClusterTools::localCovariances(seedCluster,recHits,topo) ;
114 
115  if (not isHLT_) {
116  std::array<float, 29> eval;
117 
118  const float eLeft = EcalClusterTools::eLeft(seedCluster,recHits,topo);
119  const float eRight = EcalClusterTools::eRight(seedCluster,recHits,topo);
120  const float eTop = EcalClusterTools::eTop(seedCluster,recHits,topo);
121  const float eBottom = EcalClusterTools::eBottom(seedCluster,recHits,topo);
122 
123  float sigmaIetaIeta = sqrt(localCovariances[0]);
124  float sigmaIetaIphi = std::numeric_limits<float>::max();
125  float sigmaIphiIphi = std::numeric_limits<float>::max();
126 
127  // extra shower shapes
128  const float see_by_spp = sigmaIetaIeta*sigmaIphiIphi;
129  if( see_by_spp > 0 ) {
130  sigmaIetaIphi = localCovariances[1] / see_by_spp;
131  } else if ( localCovariances[1] > 0 ) {
132  sigmaIetaIphi = 1.f;
133  } else {
134  sigmaIetaIphi = -1.f;
135  }
136 
137  if (!edm::isNotFinite(localCovariances[2])) sigmaIphiIphi = sqrt(localCovariances[2]) ;
138 
139  // calculate sub-cluster variables
140  std::vector<float> clusterRawEnergy;
141  clusterRawEnergy.resize(std::max(3, numberOfClusters), 0);
142  std::vector<float> clusterDEtaToSeed;
143  clusterDEtaToSeed.resize(std::max(3, numberOfClusters), 0);
144  std::vector<float> clusterDPhiToSeed;
145  clusterDPhiToSeed.resize(std::max(3, numberOfClusters), 0);
146  float clusterMaxDR = 999.;
147  float clusterMaxDRDPhi = 999.;
148  float clusterMaxDRDEta = 999.;
149  float clusterMaxDRRawEnergy = 0.;
150 
151  size_t iclus = 0;
152  float maxDR = 0;
154  const edm::Ptr<reco::CaloCluster>& theseed = sc.seed();
155  // loop over all clusters that aren't the seed
156  auto clusend = sc.clustersEnd();
157  for( auto clus = sc.clustersBegin(); clus != clusend; ++clus ) {
158  pclus = *clus;
159 
160  if(theseed == pclus )
161  continue;
162  clusterRawEnergy[iclus] = pclus->energy();
163  clusterDPhiToSeed[iclus] = reco::deltaPhi(pclus->phi(),theseed->phi());
164  clusterDEtaToSeed[iclus] = pclus->eta() - theseed->eta();
165 
166  // find cluster with max dR
167  const auto the_dr = reco::deltaR(*pclus, *theseed);
168  if(the_dr > maxDR) {
169  maxDR = the_dr;
170  clusterMaxDR = maxDR;
171  clusterMaxDRDPhi = clusterDPhiToSeed[iclus];
172  clusterMaxDRDEta = clusterDEtaToSeed[iclus];
173  clusterMaxDRRawEnergy = clusterRawEnergy[iclus];
174  }
175  ++iclus;
176  }
177 
178  // SET INPUTS
179  eval[0] = vertices_->size();
180  eval[1] = raw_energy;
181  eval[2] = sc.etaWidth();
182  eval[3] = sc.phiWidth();
183  eval[4] = EcalClusterTools::e3x3(seedCluster,recHits,topo)/raw_energy;
184  eval[5] = seedCluster.energy()/raw_energy;
185  eval[6] = EcalClusterTools::eMax(seedCluster,recHits)/raw_energy;
186  eval[7] = EcalClusterTools::e2nd(seedCluster,recHits)/raw_energy;
187  eval[8] = (eLeft + eRight != 0.f ? (eLeft-eRight)/(eLeft+eRight) : 0.f);
188  eval[9] = (eTop + eBottom != 0.f ? (eTop-eBottom)/(eTop+eBottom) : 0.f);
189  eval[10] = sigmaIetaIeta;
190  eval[11] = sigmaIetaIphi;
191  eval[12] = sigmaIphiIphi;
192  eval[13] = std::max(0,numberOfClusters-1);
193  eval[14] = clusterMaxDR;
194  eval[15] = clusterMaxDRDPhi;
195  eval[16] = clusterMaxDRDEta;
196  eval[17] = clusterMaxDRRawEnergy/raw_energy;
197  eval[18] = clusterRawEnergy[0]/raw_energy;
198  eval[19] = clusterRawEnergy[1]/raw_energy;
199  eval[20] = clusterRawEnergy[2]/raw_energy;
200  eval[21] = clusterDPhiToSeed[0];
201  eval[22] = clusterDPhiToSeed[1];
202  eval[23] = clusterDPhiToSeed[2];
203  eval[24] = clusterDEtaToSeed[0];
204  eval[25] = clusterDEtaToSeed[1];
205  eval[26] = clusterDEtaToSeed[2];
206  if (iseb) {
207  EBDetId ebseedid(seedCluster.seed());
208  eval[27] = ebseedid.ieta();
209  eval[28] = ebseedid.iphi();
210  } else {
211  EEDetId eeseedid(seedCluster.seed());
212  eval[27] = eeseedid.ix();
213  eval[28] = eeseedid.iy();
214  }
215 
216  //magic numbers for MINUIT-like transformation of BDT output onto limited range
217  //(These should be stored inside the conditions object in the future as well)
218  constexpr double meanlimlow = 0.2;
219  constexpr double meanlimhigh = 2.0;
220  constexpr double meanoffset = meanlimlow + 0.5*(meanlimhigh-meanlimlow);
221  constexpr double meanscale = 0.5*(meanlimhigh-meanlimlow);
222 
223  constexpr double sigmalimlow = 0.0002;
224  constexpr double sigmalimhigh = 0.5;
225  constexpr double sigmaoffset = sigmalimlow + 0.5*(sigmalimhigh-sigmalimlow);
226  constexpr double sigmascale = 0.5*(sigmalimhigh-sigmalimlow);
227 
228  const GBRForestD *forestmean = iseb ? foresteb_ : forestee_;
229  const GBRForestD *forestsigma = iseb ? forestsigmaeb_ : forestsigmaee_;
230 
231  //these are the actual BDT responses
232  double rawmean = forestmean->GetResponse(eval.data());
233  double rawsigma = forestsigma->GetResponse(eval.data());
234 
235  //apply transformation to limited output range (matching the training)
236  double mean = meanoffset + meanscale*vdt::fast_sin(rawmean);
237  double sigma = sigmaoffset + sigmascale*vdt::fast_sin(rawsigma);
238 
239  double ecor = mean*(eval[1]);
240  const double sigmacor = sigma*ecor;
241 
242  p.first = ecor;
243  p.second = sigmacor;
244 
245  } else {
246 
247  std::array<float, 7> eval;
248  float clusterMaxDR = 999.;
249 
250  size_t iclus = 0;
251  float maxDR = 0;
253  const edm::Ptr<reco::CaloCluster>& theseed = sc.seed();
254 
255  // loop over all clusters that aren't the seed
256  auto clusend = sc.clustersEnd();
257  for( auto clus = sc.clustersBegin(); clus != clusend; ++clus ) {
258  pclus = *clus;
259 
260  if(theseed == pclus )
261  continue;
262 
263  // find cluster with max dR
264  const auto the_dr = reco::deltaR(*pclus, *theseed);
265  if(the_dr > maxDR) {
266  maxDR = the_dr;
267  clusterMaxDR = maxDR;
268  }
269  ++iclus;
270  }
271 
272  // SET INPUTS
273  eval[0] = nHitsAboveThreshold_;
274  eval[1] = sc.eta();
275  eval[2] = sc.phiWidth();
276  eval[3] = EcalClusterTools::e3x3(seedCluster,recHits,topo)/raw_energy;
277  eval[4] = std::max(0,numberOfClusters-1);
278  eval[5] = clusterMaxDR;
279  eval[6] = raw_energy;
280 
281  //magic numbers for MINUIT-like transformation of BDT output onto limited range
282  //(These should be stored inside the conditions object in the future as well)
283  constexpr double meanlimlow = 0.2;
284  constexpr double meanlimhigh = 2.0;
285  constexpr double meanoffset = meanlimlow + 0.5*(meanlimhigh-meanlimlow);
286  constexpr double meanscale = 0.5*(meanlimhigh-meanlimlow);
287 
288  const GBRForestD *forestmean = iseb ? foresteb_ : forestee_;
289 
290  double rawmean = forestmean->GetResponse(eval.data());
291  double mean = meanoffset + meanscale*vdt::fast_sin(rawmean);
292 
293  double ecor = mean*eval[6];
294  if (!iseb)
295  ecor = mean*eval[6]+sc.preshowerEnergy();
296 
297  p.first = ecor;
298  //p.second unchanged
299  }
300 
301  return p;
302 }
double GetResponse(const float *vector) const
Definition: GBRForestD.h:51
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:74
int ix() const
Definition: EEDetId.h:76
double phiWidth() const
obtain phi and eta width of the Super Cluster
Definition: SuperCluster.h:55
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:192
const CaloClusterPtrVector & clusters() const
const access to the cluster list itself
Definition: SuperCluster.h:69
#define constexpr
edm::Handle< EcalRecHitCollection > rechitsEE_
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:163
bool isNotFinite(T x)
Definition: isFinite.h:10
double etaWidth() const
Definition: SuperCluster.h:56
edm::Handle< EcalRecHitCollection > rechitsEB_
edm::ESHandle< CaloTopology > calotopo_
T sqrt(T t)
Definition: SSEVec.h:18
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
double energy() const
cluster energy
Definition: CaloCluster.h:121
int ieta() const
get the crystal ieta
Definition: EBDetId.h:51
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:47
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:202
T const * product() const
Definition: ESHandle.h:86
edm::Handle< reco::VertexCollection > vertices_
CaloCluster_iterator clustersBegin() const
fist iterator over BasicCluster constituents
Definition: SuperCluster.h:75
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:66
double preshowerEnergy() const
energy deposited in preshower
Definition: SuperCluster.h:50
CaloCluster_iterator clustersEnd() const
last iterator over BasicCluster constituents
Definition: SuperCluster.h:78
void SCEnergyCorrectorSemiParm::modifyObject ( reco::SuperCluster sc)

Definition at line 305 of file SCEnergyCorrectorSemiParm.cc.

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

305  {
306 
307  std::pair<double, double> cor = getCorrections(sc);
308  if(cor.first<0) return;
309  sc.setEnergy(cor.first);
310  sc.setCorrectedEnergy(cor.first);
311  if(! isHLT_ && cor.second>=0.) sc.setCorrectedEnergyUncertainty(cor.second);
312 }
void setEnergy(double energy)
Definition: CaloCluster.h:108
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:109
void setCorrectedEnergyUncertainty(float energyerr)
Definition: CaloCluster.h:110
std::pair< double, double > getCorrections(const reco::SuperCluster &sc) const
void SCEnergyCorrectorSemiParm::setEvent ( const edm::Event e)

Definition at line 71 of file SCEnergyCorrectorSemiParm.cc.

References edm::SortedCollection< T, SORT >::begin(), edm::SortedCollection< T, SORT >::end(), eThreshold_, edm::Event::getByToken(), isHLT_, nHitsAboveThreshold_, rechitsEB_, rechitsEE_, tokenEBRecHits_, tokenEERecHits_, tokenVertices_, and vertices_.

71  {
72 
75 
76  if (not isHLT_)
78  else {
80  const EcalRecHitCollection *recHitsEB = rechitsEB_.product();
81  const EcalRecHitCollection *recHitsEE = rechitsEE_.product();
82 
83  for (EcalRecHitCollection::const_iterator it=recHitsEB->begin(); it!=recHitsEB->end(); ++it) {
84  if (it->energy() > eThreshold_)
86  }
87 
88  for (EcalRecHitCollection::const_iterator it=recHitsEE->begin(); it!=recHitsEE->end(); ++it) {
89  if (it->energy() > eThreshold_)
91  }
92  }
93 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
std::vector< EcalRecHit >::const_iterator const_iterator
edm::Handle< EcalRecHitCollection > rechitsEE_
edm::EDGetTokenT< EcalRecHitCollection > tokenEERecHits_
edm::Handle< EcalRecHitCollection > rechitsEB_
const_iterator end() const
edm::EDGetTokenT< EcalRecHitCollection > tokenEBRecHits_
edm::Handle< reco::VertexCollection > vertices_
edm::EDGetTokenT< reco::VertexCollection > tokenVertices_
const_iterator begin() const
void SCEnergyCorrectorSemiParm::setEventSetup ( const edm::EventSetup es)

Definition at line 49 of file SCEnergyCorrectorSemiParm.cc.

References calogeom_, calotopo_, foresteb_, forestee_, forestsigmaeb_, forestsigmaee_, edm::EventSetup::get(), edm::ESHandle< class >::product(), regressionKeyEB_, regressionKeyEE_, uncertaintyKeyEB_, and uncertaintyKeyEE_.

49  {
50 
51  es.get<CaloTopologyRecord>().get(calotopo_);
52  es.get<CaloGeometryRecord>().get(calogeom_);
53 
55  edm::ESHandle<GBRForestD> readerebvar;
57  edm::ESHandle<GBRForestD> readereevar;
58 
59  es.get<GBRDWrapperRcd>().get(regressionKeyEB_, readereb);
60  es.get<GBRDWrapperRcd>().get(uncertaintyKeyEB_, readerebvar);
61  es.get<GBRDWrapperRcd>().get(regressionKeyEE_, readeree);
62  es.get<GBRDWrapperRcd>().get(uncertaintyKeyEE_, readereevar);
63 
64  foresteb_ = readereb.product();
65  forestsigmaeb_ = readerebvar.product();
66  forestee_ = readeree.product();
67  forestsigmaee_ = readereevar.product();
68 }
edm::ESHandle< CaloTopology > calotopo_
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
edm::ESHandle< CaloGeometry > calogeom_
void SCEnergyCorrectorSemiParm::setTokens ( const edm::ParameterSet iConfig,
edm::ConsumesCollector cc 
)

Definition at line 30 of file SCEnergyCorrectorSemiParm.cc.

References edm::ConsumesCollector::consumes(), eThreshold_, edm::ParameterSet::getParameter(), isHLT_, regressionKeyEB_, regressionKeyEE_, AlCaHLTBitMon_QueryRunRegistry::string, tokenEBRecHits_, tokenEERecHits_, tokenVertices_, uncertaintyKeyEB_, and uncertaintyKeyEE_.

30  {
31 
32  isHLT_ = iConfig.getParameter<bool>("isHLT");
33 
36 
37  regressionKeyEB_ = iConfig.getParameter<std::string>("regressionKeyEB");
38  uncertaintyKeyEB_ = iConfig.getParameter<std::string>("uncertaintyKeyEB");
39  regressionKeyEE_ = iConfig.getParameter<std::string>("regressionKeyEE");
40  uncertaintyKeyEE_ = iConfig.getParameter<std::string>("uncertaintyKeyEE");
41 
42  if (not isHLT_)
43  tokenVertices_ = cc.consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexCollection"));
44  else
45  eThreshold_ = iConfig.getParameter<double>("eRecHitThreshold");
46 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< EcalRecHitCollection > tokenEERecHits_
edm::EDGetTokenT< EcalRecHitCollection > tokenEBRecHits_
edm::EDGetTokenT< reco::VertexCollection > tokenVertices_

Member Data Documentation

edm::ESHandle<CaloGeometry> SCEnergyCorrectorSemiParm::calogeom_
protected

Definition at line 50 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup().

edm::ESHandle<CaloTopology> SCEnergyCorrectorSemiParm::calotopo_
protected

Definition at line 49 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setEventSetup().

edm::InputTag SCEnergyCorrectorSemiParm::ecalHitsEBInputTag_
protected

Definition at line 60 of file SCEnergyCorrectorSemiParm.h.

edm::InputTag SCEnergyCorrectorSemiParm::ecalHitsEEInputTag_
protected

Definition at line 61 of file SCEnergyCorrectorSemiParm.h.

float SCEnergyCorrectorSemiParm::eThreshold_
private

Definition at line 72 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

const GBRForestD* SCEnergyCorrectorSemiParm::foresteb_
protected

Definition at line 44 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setEventSetup().

const GBRForestD* SCEnergyCorrectorSemiParm::forestee_
protected

Definition at line 45 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setEventSetup().

const GBRForestD* SCEnergyCorrectorSemiParm::forestsigmaeb_
protected

Definition at line 46 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setEventSetup().

const GBRForestD* SCEnergyCorrectorSemiParm::forestsigmaee_
protected

Definition at line 47 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setEventSetup().

bool SCEnergyCorrectorSemiParm::isHLT_
private

Definition at line 70 of file SCEnergyCorrectorSemiParm.h.

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

int SCEnergyCorrectorSemiParm::nHitsAboveThreshold_
private

Definition at line 71 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setEvent().

edm::Handle<EcalRecHitCollection> SCEnergyCorrectorSemiParm::rechitsEB_
protected

Definition at line 57 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setEvent().

edm::Handle<EcalRecHitCollection> SCEnergyCorrectorSemiParm::rechitsEE_
protected

Definition at line 58 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setEvent().

std::string SCEnergyCorrectorSemiParm::regressionKeyEB_
protected

Definition at line 64 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup(), and setTokens().

std::string SCEnergyCorrectorSemiParm::regressionKeyEE_
protected

Definition at line 66 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup(), and setTokens().

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

Definition at line 52 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

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

Definition at line 53 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

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

Definition at line 54 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

std::string SCEnergyCorrectorSemiParm::uncertaintyKeyEB_
protected

Definition at line 65 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup(), and setTokens().

std::string SCEnergyCorrectorSemiParm::uncertaintyKeyEE_
protected

Definition at line 67 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup(), and setTokens().

edm::InputTag SCEnergyCorrectorSemiParm::vertexInputTag_
protected

Definition at line 62 of file SCEnergyCorrectorSemiParm.h.

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

Definition at line 56 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setEvent().