CMS 3D CMS Logo

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< EcalRecHitCollectiontokenEBRecHits_
 
edm::EDGetTokenT< EcalRecHitCollectiontokenEERecHits_
 
edm::EDGetTokenT< reco::VertexCollectiontokenVertices_
 
std::string uncertaintyKeyEB_
 
std::string uncertaintyKeyEE_
 
edm::InputTag vertexInputTag_
 
edm::Handle< reco::VertexCollectionvertices_
 

Private Attributes

bool applySigmaIetaIphiBug_
 
float eThreshold_
 
bool isHLT_
 
int nHitsAboveThreshold_
 

Detailed Description

Definition at line 30 of file SCEnergyCorrectorSemiParm.h.

Constructor & Destructor Documentation

SCEnergyCorrectorSemiParm::SCEnergyCorrectorSemiParm ( )

Definition at line 17 of file SCEnergyCorrectorSemiParm.cc.

17  :
18 foresteb_(nullptr),
19 forestee_(nullptr),
20 forestsigmaeb_(nullptr),
21 forestsigmaee_(nullptr),
22 calotopo_(nullptr),
23 calogeom_(nullptr)
24 {}
edm::ESHandle< CaloTopology > calotopo_
edm::ESHandle< CaloGeometry > calogeom_
SCEnergyCorrectorSemiParm::~SCEnergyCorrectorSemiParm ( )

Definition at line 27 of file SCEnergyCorrectorSemiParm.cc.

28 {}

Member Function Documentation

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

Definition at line 102 of file SCEnergyCorrectorSemiParm.cc.

References applySigmaIetaIphiBug_, calotopo_, reco::SuperCluster::clusters(), reco::SuperCluster::clustersBegin(), reco::SuperCluster::clustersEnd(), constexpr, reco::deltaPhi(), reco::deltaR(), DetId::det(), EcalBarrel, photonPostprocessing_cfi::eMax, reco::CaloCluster::energy(), reco::CaloCluster::eta(), reco::SuperCluster::etaWidth(), foresteb_, forestee_, forestsigmaeb_, forestsigmaee_, GBRForestD::GetResponse(), reco::CaloCluster::hitsAndFractions(), EBDetId::ieta(), EcalTools::isHGCalDet(), isHLT_, edm::isNotFinite(), EEDetId::ix(), SiStripPI::max, jets_cff::maxDR, SiStripPI::mean, nHitsAboveThreshold_, AlCaHLTBitMon_ParallelJobs::p, reco::CaloCluster::phi(), reco::SuperCluster::phiWidth(), reco::SuperCluster::preshowerEnergy(), edm::ESHandle< T >::product(), reco::SuperCluster::rawEnergy(), rechitsEB_, rechitsEE_, reco::SuperCluster::seed(), reco::CaloCluster::seed(), edm::PtrVectorBase::size(), mathSSE::sqrt(), and vertices_.

Referenced by modifyObject().

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

Definition at line 315 of file SCEnergyCorrectorSemiParm.cc.

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

315  {
316 
317  std::pair<double, double> cor = getCorrections(sc);
318  if(cor.first<0) return;
319  sc.setEnergy(cor.first);
320  sc.setCorrectedEnergy(cor.first);
321  if(! isHLT_ && cor.second>=0.) sc.setCorrectedEnergyUncertainty(cor.second);
322 }
void setEnergy(double energy)
Definition: CaloCluster.h:113
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:114
void setCorrectedEnergyUncertainty(float energyerr)
Definition: CaloCluster.h:115
std::pair< double, double > getCorrections(const reco::SuperCluster &sc) const
void SCEnergyCorrectorSemiParm::setEvent ( const edm::Event e)

Definition at line 73 of file SCEnergyCorrectorSemiParm.cc.

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

73  {
74 
77 
78  if (not isHLT_)
80  else {
82  const EcalRecHitCollection *recHitsEB = ( rechitsEB_.isValid() ? rechitsEB_.product() : nullptr );
83  const EcalRecHitCollection *recHitsEE = ( rechitsEE_.isValid() ? rechitsEE_.product() : nullptr );
84 
85  if( nullptr != recHitsEB ) {
86  for (EcalRecHitCollection::const_iterator it=recHitsEB->begin(); it!=recHitsEB->end(); ++it) {
87  if (it->energy() > eThreshold_)
89  }
90  }
91 
92  if( nullptr != recHitsEE ) {
93  for (EcalRecHitCollection::const_iterator it=recHitsEE->begin(); it!=recHitsEE->end(); ++it) {
94  if (it->energy() > eThreshold_)
96  }
97  }
98  }
99 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
#define nullptr
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 51 of file SCEnergyCorrectorSemiParm.cc.

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

51  {
52 
53  es.get<CaloTopologyRecord>().get(calotopo_);
54  es.get<CaloGeometryRecord>().get(calogeom_);
55 
57  edm::ESHandle<GBRForestD> readerebvar;
59  edm::ESHandle<GBRForestD> readereevar;
60 
61  es.get<GBRDWrapperRcd>().get(regressionKeyEB_, readereb);
62  es.get<GBRDWrapperRcd>().get(uncertaintyKeyEB_, readerebvar);
63  es.get<GBRDWrapperRcd>().get(regressionKeyEE_, readeree);
64  es.get<GBRDWrapperRcd>().get(uncertaintyKeyEE_, readereevar);
65 
66  foresteb_ = readereb.product();
67  forestsigmaeb_ = readerebvar.product();
68  forestee_ = readeree.product();
69  forestsigmaee_ = readereevar.product();
70 }
edm::ESHandle< CaloTopology > calotopo_
T get() const
Definition: EventSetup.h:71
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 31 of file SCEnergyCorrectorSemiParm.cc.

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

31  {
32 
33  isHLT_ = iConfig.getParameter<bool>("isHLT");
34  applySigmaIetaIphiBug_ = iConfig.getParameter<bool>("applySigmaIetaIphiBug");
37 
38  regressionKeyEB_ = iConfig.getParameter<std::string>("regressionKeyEB");
39  uncertaintyKeyEB_ = iConfig.getParameter<std::string>("uncertaintyKeyEB");
40  regressionKeyEE_ = iConfig.getParameter<std::string>("regressionKeyEE");
41  uncertaintyKeyEE_ = iConfig.getParameter<std::string>("uncertaintyKeyEE");
42 
43  if (not isHLT_){
44  tokenVertices_ = cc.consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexCollection"));
45  }else{
46  eThreshold_ = iConfig.getParameter<double>("eRecHitThreshold");
47  }
48 }
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

bool SCEnergyCorrectorSemiParm::applySigmaIetaIphiBug_
private

Definition at line 71 of file SCEnergyCorrectorSemiParm.h.

Referenced by getCorrections(), and setTokens().

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 73 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 72 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().