CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
ESRecHitWorker Class Reference

#include <ESRecHitWorker.h>

Inheritance diagram for ESRecHitWorker:
ESRecHitWorkerBaseClass

Public Member Functions

 ESRecHitWorker (const edm::ParameterSet &ps, edm::ConsumesCollector cc)
 
bool run (const ESDigiCollection::const_iterator &digi, ESRecHitCollection &result) override
 
void set (const edm::EventSetup &es) override
 
 ~ESRecHitWorker () override
 
- Public Member Functions inherited from ESRecHitWorkerBaseClass
 ESRecHitWorkerBaseClass (const edm::ParameterSet &)
 
virtual ~ESRecHitWorkerBaseClass ()
 

Private Attributes

ESRecHitAnalyticAlgoalgoA_
 
ESRecHitFitAlgoalgoF_
 
ESRecHitSimAlgoalgoW_
 
edm::ESHandle< ESAngleCorrectionFactorsesAngleCorrFactors_
 
edm::ESGetToken< ESAngleCorrectionFactors, ESAngleCorrectionFactorsRcdesAngleCorrFactorsToken_
 
edm::ESHandle< ESChannelStatusesChannelStatus_
 
edm::ESGetToken< ESChannelStatus, ESChannelStatusRcdesChannelStatusToken_
 
edm::ESHandle< ESGainesgain_
 
edm::ESGetToken< ESGain, ESGainRcdesgainToken_
 
edm::ESHandle< ESIntercalibConstantsesMIPs_
 
edm::ESGetToken< ESIntercalibConstants, ESIntercalibConstantsRcdesMIPsToken_
 
edm::ESHandle< ESMIPToGeVConstantesMIPToGeV_
 
edm::ESGetToken< ESMIPToGeVConstant, ESMIPToGeVConstantRcdesMIPToGeVToken_
 
edm::ESHandle< ESPedestalsesPedestals_
 
edm::ESGetToken< ESPedestals, ESPedestalsRcdesPedestalsToken_
 
edm::ESHandle< ESRecHitRatioCutsesRatioCuts_
 
edm::ESGetToken< ESRecHitRatioCuts, ESRecHitRatioCutsRcdesRatioCutsToken_
 
edm::ESHandle< ESTimeSampleWeightsesWeights_
 
edm::ESGetToken< ESTimeSampleWeights, ESTimeSampleWeightsRcdesWeightsToken_
 
int recoAlgo_
 

Detailed Description

Definition at line 37 of file ESRecHitWorker.h.

Constructor & Destructor Documentation

◆ ESRecHitWorker()

ESRecHitWorker::ESRecHitWorker ( const edm::ParameterSet ps,
edm::ConsumesCollector  cc 
)

Definition at line 14 of file ESRecHitWorker.cc.

References algoA_, algoF_, algoW_, gpuPixelDoublets::cc, esAngleCorrFactorsToken_, esChannelStatusToken_, esgainToken_, esMIPsToken_, esMIPToGeVToken_, esPedestalsToken_, esRatioCutsToken_, esWeightsToken_, edm::ParameterSet::getParameter(), and recoAlgo_.

15  recoAlgo_ = ps.getParameter<int>("ESRecoAlgo");
16  esgainToken_ = cc.esConsumes<ESGain, ESGainRcd>();
24 
25  if (recoAlgo_ == 0)
26  algoW_ = new ESRecHitSimAlgo();
27  else if (recoAlgo_ == 1)
28  algoF_ = new ESRecHitFitAlgo();
29  else
31 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Definition: ESGain.h:7
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
ESRecHitSimAlgo * algoW_
ESPedestalsMap ESPedestals
Definition: ESPedestals.h:29
ESRecHitWorkerBaseClass(const edm::ParameterSet &)
ESChannelStatusMap ESChannelStatus
edm::ESGetToken< ESAngleCorrectionFactors, ESAngleCorrectionFactorsRcd > esAngleCorrFactorsToken_
ESRecHitFitAlgo * algoF_
ESAngleCorrectionFactorMap ESAngleCorrectionFactors
edm::ESGetToken< ESPedestals, ESPedestalsRcd > esPedestalsToken_
edm::ESGetToken< ESRecHitRatioCuts, ESRecHitRatioCutsRcd > esRatioCutsToken_
edm::ESGetToken< ESTimeSampleWeights, ESTimeSampleWeightsRcd > esWeightsToken_
ESIntercalibConstantMap ESIntercalibConstants
ESRecHitAnalyticAlgo * algoA_
edm::ESGetToken< ESGain, ESGainRcd > esgainToken_
edm::ESGetToken< ESMIPToGeVConstant, ESMIPToGeVConstantRcd > esMIPToGeVToken_
edm::ESGetToken< ESChannelStatus, ESChannelStatusRcd > esChannelStatusToken_
edm::ESGetToken< ESIntercalibConstants, ESIntercalibConstantsRcd > esMIPsToken_

◆ ~ESRecHitWorker()

ESRecHitWorker::~ESRecHitWorker ( )
override

Definition at line 33 of file ESRecHitWorker.cc.

References algoA_, algoF_, algoW_, and recoAlgo_.

33  {
34  if (recoAlgo_ == 0)
35  delete algoW_;
36  else if (recoAlgo_ == 1)
37  delete algoF_;
38  else
39  delete algoA_;
40 }
ESRecHitSimAlgo * algoW_
ESRecHitFitAlgo * algoF_
ESRecHitAnalyticAlgo * algoA_

Member Function Documentation

◆ run()

bool ESRecHitWorker::run ( const ESDigiCollection::const_iterator digi,
ESRecHitCollection result 
)
overridevirtual

Implements ESRecHitWorkerBaseClass.

Definition at line 104 of file ESRecHitWorker.cc.

References algoA_, algoF_, algoW_, recoAlgo_, ESRecHitSimAlgo::reconstruct(), ESRecHitAnalyticAlgo::reconstruct(), ESRecHitFitAlgo::reconstruct(), and mps_fire::result.

104  {
105  if (recoAlgo_ == 0)
106  result.push_back(algoW_->reconstruct(*itdg));
107  else if (recoAlgo_ == 1)
108  result.push_back(algoF_->reconstruct(*itdg));
109  else
110  result.push_back(algoA_->reconstruct(*itdg));
111  return true;
112 }
EcalRecHit reconstruct(const ESDataFrame &digi) const
EcalRecHit reconstruct(const ESDataFrame &digi) const
ESRecHitSimAlgo * algoW_
ESRecHitFitAlgo * algoF_
EcalRecHit reconstruct(const ESDataFrame &digi) const
ESRecHitAnalyticAlgo * algoA_

◆ set()

void ESRecHitWorker::set ( const edm::EventSetup es)
overridevirtual

Implements ESRecHitWorkerBaseClass.

Definition at line 42 of file ESRecHitWorker.cc.

References algoA_, algoF_, algoW_, esAngleCorrFactors_, esAngleCorrFactorsToken_, esChannelStatus_, esChannelStatusToken_, esgain_, esgainToken_, esMIPs_, esMIPsToken_, esMIPToGeV_, esMIPToGeVToken_, esPedestals_, esPedestalsToken_, esRatioCuts_, esRatioCutsToken_, esWeights_, esWeightsToken_, PedestalClient_cfi::gain, ESMIPToGeVConstant::getESValueHigh(), ESMIPToGeVConstant::getESValueLow(), edm::EventSetup::getHandle(), ESTimeSampleWeights::getWeightForTS0(), ESTimeSampleWeights::getWeightForTS1(), ESTimeSampleWeights::getWeightForTS2(), edm::ESHandle< T >::product(), recoAlgo_, ESRecHitSimAlgo::setAngleCorrectionFactors(), ESRecHitAnalyticAlgo::setAngleCorrectionFactors(), ESRecHitFitAlgo::setAngleCorrectionFactors(), ESRecHitSimAlgo::setChannelStatus(), ESRecHitAnalyticAlgo::setChannelStatus(), ESRecHitFitAlgo::setChannelStatus(), ESRecHitSimAlgo::setESGain(), ESRecHitAnalyticAlgo::setESGain(), ESRecHitFitAlgo::setESGain(), ESRecHitSimAlgo::setIntercalibConstants(), ESRecHitAnalyticAlgo::setIntercalibConstants(), ESRecHitFitAlgo::setIntercalibConstants(), ESRecHitSimAlgo::setMIPGeV(), ESRecHitAnalyticAlgo::setMIPGeV(), ESRecHitFitAlgo::setMIPGeV(), ESRecHitSimAlgo::setPedestals(), ESRecHitAnalyticAlgo::setPedestals(), ESRecHitFitAlgo::setPedestals(), ESRecHitSimAlgo::setRatioCuts(), ESRecHitAnalyticAlgo::setRatioCuts(), ESRecHitFitAlgo::setRatioCuts(), ESRecHitSimAlgo::setW0(), ESRecHitSimAlgo::setW1(), ESRecHitSimAlgo::setW2(), w1, and w2.

42  {
44  const ESGain *gain = esgain_.product();
45 
47  const ESMIPToGeVConstant *mipToGeV = esMIPToGeV_.product();
48 
49  double ESGain = gain->getESGain();
50  double ESMIPToGeV = (ESGain == 1) ? mipToGeV->getESValueLow() : mipToGeV->getESValueHigh();
51 
53  const ESTimeSampleWeights *wgts = esWeights_.product();
54 
55  float w0 = wgts->getWeightForTS0();
56  float w1 = wgts->getWeightForTS1();
57  float w2 = wgts->getWeightForTS2();
58 
60  const ESPedestals *peds = esPedestals_.product();
61 
63  const ESIntercalibConstants *mips = esMIPs_.product();
64 
66  const ESAngleCorrectionFactors *ang = esAngleCorrFactors_.product();
67 
69  const ESChannelStatus *channelStatus = esChannelStatus_.product();
70 
72  const ESRecHitRatioCuts *ratioCuts = esRatioCuts_.product();
73 
74  if (recoAlgo_ == 0) {
76  algoW_->setMIPGeV(ESMIPToGeV);
77  algoW_->setW0(w0);
78  algoW_->setW1(w1);
79  algoW_->setW2(w2);
80  algoW_->setPedestals(peds);
82  algoW_->setChannelStatus(channelStatus);
83  algoW_->setRatioCuts(ratioCuts);
85  } else if (recoAlgo_ == 1) {
87  algoF_->setMIPGeV(ESMIPToGeV);
88  algoF_->setPedestals(peds);
90  algoF_->setChannelStatus(channelStatus);
91  algoF_->setRatioCuts(ratioCuts);
93  } else {
95  algoA_->setMIPGeV(ESMIPToGeV);
96  algoA_->setPedestals(peds);
98  algoA_->setChannelStatus(channelStatus);
99  algoA_->setRatioCuts(ratioCuts);
101  }
102 }
void setIntercalibConstants(const ESIntercalibConstants *mips)
edm::ESHandle< ESGain > esgain_
float getWeightForTS0() const
float getESValueHigh() 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)
void setW0(float value)
void setESGain(const double &value)
ESRecHitSimAlgo * algoW_
void setPedestals(const ESPedestals *peds)
void setMIPGeV(float value)
void setPedestals(const ESPedestals *peds)
edm::ESGetToken< ESAngleCorrectionFactors, ESAngleCorrectionFactorsRcd > esAngleCorrFactorsToken_
edm::ESHandle< ESRecHitRatioCuts > esRatioCuts_
void setW1(float value)
void setIntercalibConstants(const ESIntercalibConstants *mips)
T const * product() const
Definition: ESHandle.h:86
void setRatioCuts(const ESRecHitRatioCuts *ratioCuts)
ESRecHitFitAlgo * algoF_
edm::ESGetToken< ESPedestals, ESPedestalsRcd > esPedestalsToken_
edm::ESGetToken< ESRecHitRatioCuts, ESRecHitRatioCutsRcd > esRatioCutsToken_
void setMIPGeV(const double &value)
edm::ESHandle< ESPedestals > esPedestals_
void setMIPGeV(const double &value)
edm::ESHandle< ESChannelStatus > esChannelStatus_
float getWeightForTS2() const
float getESValueLow() const
void setRatioCuts(const ESRecHitRatioCuts *ratioCuts)
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
void setESGain(float value)
edm::ESHandle< ESTimeSampleWeights > esWeights_
edm::ESHandle< ESMIPToGeVConstant > esMIPToGeV_
edm::ESHandle< ESAngleCorrectionFactors > esAngleCorrFactors_
edm::ESGetToken< ESTimeSampleWeights, ESTimeSampleWeightsRcd > esWeightsToken_
weight_default_t w1[2000]
Definition: w1.h:9
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)
void setChannelStatus(const ESChannelStatus *status)
void setESGain(const double &value)
void setChannelStatus(const ESChannelStatus *status)
ESRecHitAnalyticAlgo * algoA_
edm::ESGetToken< ESGain, ESGainRcd > esgainToken_
edm::ESHandle< ESIntercalibConstants > esMIPs_
edm::ESGetToken< ESMIPToGeVConstant, ESMIPToGeVConstantRcd > esMIPToGeVToken_
float getWeightForTS1() const
edm::ESGetToken< ESChannelStatus, ESChannelStatusRcd > esChannelStatusToken_
edm::ESGetToken< ESIntercalibConstants, ESIntercalibConstantsRcd > esMIPsToken_

Member Data Documentation

◆ algoA_

ESRecHitAnalyticAlgo* ESRecHitWorker::algoA_
private

Definition at line 49 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), run(), set(), and ~ESRecHitWorker().

◆ algoF_

ESRecHitFitAlgo* ESRecHitWorker::algoF_
private

Definition at line 48 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), run(), set(), and ~ESRecHitWorker().

◆ algoW_

ESRecHitSimAlgo* ESRecHitWorker::algoW_
private

Definition at line 47 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), run(), set(), and ~ESRecHitWorker().

◆ esAngleCorrFactors_

edm::ESHandle<ESAngleCorrectionFactors> ESRecHitWorker::esAngleCorrFactors_
private

Definition at line 58 of file ESRecHitWorker.h.

Referenced by set().

◆ esAngleCorrFactorsToken_

edm::ESGetToken<ESAngleCorrectionFactors, ESAngleCorrectionFactorsRcd> ESRecHitWorker::esAngleCorrFactorsToken_
private

Definition at line 66 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), and set().

◆ esChannelStatus_

edm::ESHandle<ESChannelStatus> ESRecHitWorker::esChannelStatus_
private

Definition at line 56 of file ESRecHitWorker.h.

Referenced by set().

◆ esChannelStatusToken_

edm::ESGetToken<ESChannelStatus, ESChannelStatusRcd> ESRecHitWorker::esChannelStatusToken_
private

Definition at line 64 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), and set().

◆ esgain_

edm::ESHandle<ESGain> ESRecHitWorker::esgain_
private

Definition at line 51 of file ESRecHitWorker.h.

Referenced by set().

◆ esgainToken_

edm::ESGetToken<ESGain, ESGainRcd> ESRecHitWorker::esgainToken_
private

Definition at line 59 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), and set().

◆ esMIPs_

edm::ESHandle<ESIntercalibConstants> ESRecHitWorker::esMIPs_
private

Definition at line 55 of file ESRecHitWorker.h.

Referenced by set().

◆ esMIPsToken_

edm::ESGetToken<ESIntercalibConstants, ESIntercalibConstantsRcd> ESRecHitWorker::esMIPsToken_
private

Definition at line 63 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), and set().

◆ esMIPToGeV_

edm::ESHandle<ESMIPToGeVConstant> ESRecHitWorker::esMIPToGeV_
private

Definition at line 52 of file ESRecHitWorker.h.

Referenced by set().

◆ esMIPToGeVToken_

edm::ESGetToken<ESMIPToGeVConstant, ESMIPToGeVConstantRcd> ESRecHitWorker::esMIPToGeVToken_
private

Definition at line 60 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), and set().

◆ esPedestals_

edm::ESHandle<ESPedestals> ESRecHitWorker::esPedestals_
private

Definition at line 54 of file ESRecHitWorker.h.

Referenced by set().

◆ esPedestalsToken_

edm::ESGetToken<ESPedestals, ESPedestalsRcd> ESRecHitWorker::esPedestalsToken_
private

Definition at line 62 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), and set().

◆ esRatioCuts_

edm::ESHandle<ESRecHitRatioCuts> ESRecHitWorker::esRatioCuts_
private

Definition at line 57 of file ESRecHitWorker.h.

Referenced by set().

◆ esRatioCutsToken_

edm::ESGetToken<ESRecHitRatioCuts, ESRecHitRatioCutsRcd> ESRecHitWorker::esRatioCutsToken_
private

Definition at line 65 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), and set().

◆ esWeights_

edm::ESHandle<ESTimeSampleWeights> ESRecHitWorker::esWeights_
private

Definition at line 53 of file ESRecHitWorker.h.

Referenced by set().

◆ esWeightsToken_

edm::ESGetToken<ESTimeSampleWeights, ESTimeSampleWeightsRcd> ESRecHitWorker::esWeightsToken_
private

Definition at line 61 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), and set().

◆ recoAlgo_

int ESRecHitWorker::recoAlgo_
private

Definition at line 46 of file ESRecHitWorker.h.

Referenced by ESRecHitWorker(), run(), set(), and ~ESRecHitWorker().