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

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

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

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_, GBRForestD::GetResponse(), reco::CaloCluster::hitsAndFractions(), EBDetId::ieta(), isHLT_, edm::isNotFinite(), EEDetId::ix(), bookConverter::max, timingPdfMaker::mean, nHitsAboveThreshold_, 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(), reco::CaloCluster::setCorrectedEnergy(), reco::CaloCluster::setCorrectedEnergyUncertainty(), reco::CaloCluster::setEnergy(), edm::PtrVectorBase::size(), mathSSE::sqrt(), and vertices_.

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

Referenced by setEventSetup().

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

Definition at line 48 of file SCEnergyCorrectorSemiParm.h.

Referenced by modifyObject(), and setEventSetup().

edm::InputTag SCEnergyCorrectorSemiParm::ecalHitsEBInputTag_
protected

Definition at line 59 of file SCEnergyCorrectorSemiParm.h.

edm::InputTag SCEnergyCorrectorSemiParm::ecalHitsEEInputTag_
protected

Definition at line 60 of file SCEnergyCorrectorSemiParm.h.

float SCEnergyCorrectorSemiParm::eThreshold_
private

Definition at line 71 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

const GBRForestD* SCEnergyCorrectorSemiParm::foresteb_
protected

Definition at line 43 of file SCEnergyCorrectorSemiParm.h.

Referenced by modifyObject(), and setEventSetup().

const GBRForestD* SCEnergyCorrectorSemiParm::forestee_
protected

Definition at line 44 of file SCEnergyCorrectorSemiParm.h.

Referenced by modifyObject(), and setEventSetup().

const GBRForestD* SCEnergyCorrectorSemiParm::forestsigmaeb_
protected

Definition at line 45 of file SCEnergyCorrectorSemiParm.h.

Referenced by modifyObject(), and setEventSetup().

const GBRForestD* SCEnergyCorrectorSemiParm::forestsigmaee_
protected

Definition at line 46 of file SCEnergyCorrectorSemiParm.h.

Referenced by modifyObject(), and setEventSetup().

bool SCEnergyCorrectorSemiParm::isHLT_
private

Definition at line 69 of file SCEnergyCorrectorSemiParm.h.

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

int SCEnergyCorrectorSemiParm::nHitsAboveThreshold_
private

Definition at line 70 of file SCEnergyCorrectorSemiParm.h.

Referenced by modifyObject(), and setEvent().

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

Definition at line 56 of file SCEnergyCorrectorSemiParm.h.

Referenced by modifyObject(), and setEvent().

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

Definition at line 57 of file SCEnergyCorrectorSemiParm.h.

Referenced by modifyObject(), and setEvent().

std::string SCEnergyCorrectorSemiParm::regressionKeyEB_
protected

Definition at line 63 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup(), and setTokens().

std::string SCEnergyCorrectorSemiParm::regressionKeyEE_
protected

Definition at line 65 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup(), and setTokens().

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

Definition at line 51 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

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

Definition at line 52 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

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

Definition at line 53 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEvent(), and setTokens().

std::string SCEnergyCorrectorSemiParm::uncertaintyKeyEB_
protected

Definition at line 64 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup(), and setTokens().

std::string SCEnergyCorrectorSemiParm::uncertaintyKeyEE_
protected

Definition at line 66 of file SCEnergyCorrectorSemiParm.h.

Referenced by setEventSetup(), and setTokens().

edm::InputTag SCEnergyCorrectorSemiParm::vertexInputTag_
protected

Definition at line 61 of file SCEnergyCorrectorSemiParm.h.

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

Definition at line 55 of file SCEnergyCorrectorSemiParm.h.

Referenced by modifyObject(), and setEvent().