CMS 3D CMS Logo

ESRecHitWorker.cc
Go to the documentation of this file.
2 
17 
18 #include <cmath>
19 #include <iomanip>
20 #include <iostream>
21 
24 {
25  recoAlgo_ = ps.getParameter<int>("ESRecoAlgo");
26 
27  if (recoAlgo_ == 0)
28  algoW_ = new ESRecHitSimAlgo();
29  else if (recoAlgo_ == 1)
30  algoF_ = new ESRecHitFitAlgo();
31  else
33 }
34 
36  if (recoAlgo_ == 0)
37  delete algoW_;
38  else if (recoAlgo_ == 1)
39  delete algoF_;
40  else
41  delete algoA_;
42 }
43 
45 
46  es.get<ESGainRcd>().get(esgain_);
47  const ESGain *gain = esgain_.product();
48 
50  const ESMIPToGeVConstant *mipToGeV = esMIPToGeV_.product();
51 
52  double ESGain = gain->getESGain();
53  double ESMIPToGeV = (ESGain == 1) ? mipToGeV->getESValueLow() : mipToGeV->getESValueHigh();
54 
56  const ESTimeSampleWeights *wgts = esWeights_.product();
57 
58  float w0 = wgts->getWeightForTS0();
59  float w1 = wgts->getWeightForTS1();
60  float w2 = wgts->getWeightForTS2();
61 
62  es.get<ESPedestalsRcd>().get(esPedestals_);
63  const ESPedestals *peds = esPedestals_.product();
64 
66  const ESIntercalibConstants *mips = esMIPs_.product();
67 
69  const ESAngleCorrectionFactors *ang = esAngleCorrFactors_.product();
70 
72  const ESChannelStatus *channelStatus = esChannelStatus_.product();
73 
75  const ESRecHitRatioCuts *ratioCuts = esRatioCuts_.product();
76 
77  if (recoAlgo_ == 0) {
78  algoW_->setESGain(ESGain);
79  algoW_->setMIPGeV(ESMIPToGeV);
80  algoW_->setW0(w0);
81  algoW_->setW1(w1);
82  algoW_->setW2(w2);
83  algoW_->setPedestals(peds);
85  algoW_->setChannelStatus(channelStatus);
86  algoW_->setRatioCuts(ratioCuts);
88  } else if (recoAlgo_ == 1) {
89  algoF_->setESGain(ESGain);
90  algoF_->setMIPGeV(ESMIPToGeV);
91  algoF_->setPedestals(peds);
93  algoF_->setChannelStatus(channelStatus);
94  algoF_->setRatioCuts(ratioCuts);
96  } else {
97  algoA_->setESGain(ESGain);
98  algoA_->setMIPGeV(ESMIPToGeV);
99  algoA_->setPedestals(peds);
101  algoA_->setChannelStatus(channelStatus);
102  algoA_->setRatioCuts(ratioCuts);
104  }
105 }
106 
107 bool
110 {
111  if (recoAlgo_ == 0)
112  result.push_back( algoW_->reconstruct(*itdg) );
113  else if (recoAlgo_ == 1)
114  result.push_back( algoF_->reconstruct(*itdg) );
115  else
116  result.push_back( algoA_->reconstruct(*itdg) );
117  return true;
118 }
119 
T getParameter(std::string const &) const
void setIntercalibConstants(const ESIntercalibConstants *mips)
edm::ESHandle< ESGain > esgain_
ESRecHitWorker(const edm::ParameterSet &ps)
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
EcalRecHit reconstruct(const ESDataFrame &digi) const
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
Definition: ESGain.h:7
void setPedestals(const ESPedestals *peds)
void setRatioCuts(const ESRecHitRatioCuts *ratioCuts)
float getWeightForTS1() const
float getESValueLow() const
void setW0(float value)
void push_back(T const &t)
void setESGain(const double &value)
ESRecHitSimAlgo * algoW_
void setPedestals(const ESPedestals *peds)
void setMIPGeV(float value)
void setPedestals(const ESPedestals *peds)
edm::ESHandle< ESRecHitRatioCuts > esRatioCuts_
void setW1(float value)
void setIntercalibConstants(const ESIntercalibConstants *mips)
float getWeightForTS2() const
void setRatioCuts(const ESRecHitRatioCuts *ratioCuts)
ESRecHitFitAlgo * algoF_
float getWeightForTS0() const
EcalRecHit reconstruct(const ESDataFrame &digi) const
void setMIPGeV(const double &value)
edm::ESHandle< ESPedestals > esPedestals_
void setMIPGeV(const double &value)
edm::ESHandle< ESChannelStatus > esChannelStatus_
bool run(const ESDigiCollection::const_iterator &digi, ESRecHitCollection &result) override
void setRatioCuts(const ESRecHitRatioCuts *ratioCuts)
void setESGain(float value)
edm::ESHandle< ESTimeSampleWeights > esWeights_
edm::ESHandle< ESMIPToGeVConstant > esMIPToGeV_
edm::ESHandle< ESAngleCorrectionFactors > esAngleCorrFactors_
EcalRecHit reconstruct(const ESDataFrame &digi) const
float getESGain() const
Definition: ESGain.h:13
void setIntercalibConstants(const ESIntercalibConstants *mips)
void setAngleCorrectionFactors(const ESAngleCorrectionFactors *ang)
void setChannelStatus(const ESChannelStatus *status)
void setAngleCorrectionFactors(const ESAngleCorrectionFactors *ang)
void setAngleCorrectionFactors(const ESAngleCorrectionFactors *ang)
void setW2(float value)
T get() const
Definition: EventSetup.h:63
void setChannelStatus(const ESChannelStatus *status)
void setESGain(const double &value)
void setChannelStatus(const ESChannelStatus *status)
void set(const edm::EventSetup &es) override
#define DEFINE_EDM_PLUGIN(factory, type, name)
ESRecHitAnalyticAlgo * algoA_
edm::ESHandle< ESIntercalibConstants > esMIPs_
float getESValueHigh() const
~ESRecHitWorker() override
T const * product() const
Definition: ESHandle.h:86