CMS 3D CMS Logo

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

#include <HcalHitReconstructor.h>

Inheritance diagram for HcalHitReconstructor:
edm::stream::EDProducer<>

Public Member Functions

void beginRun (edm::Run const &r, edm::EventSetup const &es) final
 
void endRun (edm::Run const &r, edm::EventSetup const &es) final
 
 HcalHitReconstructor (const edm::ParameterSet &ps)
 
void produce (edm::Event &e, const edm::EventSetup &c) override
 
 ~HcalHitReconstructor () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Types

typedef void(HcalSimpleRecAlgo::* SetCorrectionFcn) (std::shared_ptr< AbsOOTPileupCorrection >)
 

Private Attributes

std::string cat_
 
edm::ESGetToken< HcalDbService, HcalDbRecordconditionsToken_
 
bool correctTiming_
 
std::string corrName_
 
std::string dataOOTCorrectionCategory_
 
std::string dataOOTCorrectionName_
 
DetId::Detector det_
 
bool digiTimeFromDB_
 
edm::ESGetToken< HcalFlagHFDigiTimeParams, HcalFlagHFDigiTimeParamsRcddigiTimeToken_
 
bool dropZSmarkedPassed_
 
int firstAuxTS_
 
int firstSample_
 
HcalHFStatusBitFromDigishfdigibit_
 
std::unique_ptr< HcalFlagHFDigiTimeParamshFDigiTimeParams_
 
HcalHF_PETalgorithmhfPET_
 
HcalHF_S9S1algorithmhfS8S1_
 
HcalHF_S9S1algorithmhfS9S1_
 
HFTimingTrustFlagHFTimingTrustFlagSetter_
 
edm::ESGetToken< HcalTopology, HcalRecNumberingRecordhtopoToken_
 
edm::InputTag inputLabel_
 
std::string mcOOTCorrectionCategory_
 
std::string mcOOTCorrectionName_
 
edm::ESGetToken< HcalRecoParams, HcalRecoParamsRcdparamsToken_
 
std::unique_ptr< HcalRecoParamsparamTS_
 
edm::ESGetToken< HcalChannelQuality, HcalChannelQualityRcdqualToken_
 
HcalSimpleRecAlgo reco_
 
bool recoParamsFromDB_
 
int samplesToAdd_
 
HcalADCSaturationFlagsaturationFlagSetter_
 
bool setHSCPFlags_
 
bool setNegativeFlags_
 
bool setNoiseFlags_
 
SetCorrectionFcn setPileupCorrection_
 
bool setPulseShapeFlags_
 
bool setSaturationFlags_
 
bool setTimingTrustFlags_
 
edm::ESGetToken< HcalSeverityLevelComputer, HcalSeverityLevelComputerRcdsevToken_
 
int subdet_
 
HcalOtherSubdetector subdetOther_
 
edm::EDGetTokenT< HcalCalibDigiCollectiontok_calib_
 
edm::EDGetTokenT< HFDigiCollectiontok_hf_
 
edm::EDGetTokenT< HODigiCollectiontok_ho_
 
bool tsFromDB_
 
bool useLeakCorrection_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Author
J. Temple & E. Yazgan Based on HcalSimpleReconstructor.h by J. Mans

Definition at line 50 of file HcalHitReconstructor.h.

Member Typedef Documentation

◆ SetCorrectionFcn

typedef void(HcalSimpleRecAlgo::* HcalHitReconstructor::SetCorrectionFcn) (std::shared_ptr< AbsOOTPileupCorrection >)
private

Definition at line 60 of file HcalHitReconstructor.h.

Constructor & Destructor Documentation

◆ HcalHitReconstructor()

HcalHitReconstructor::HcalHitReconstructor ( const edm::ParameterSet ps)
explicit

Definition at line 26 of file HcalHitReconstructor.cc.

References DetId::Calo, conditionsToken_, dataOOTCorrectionCategory_, dataOOTCorrectionName_, det_, digiTimeFromDB_, digiTimeToken_, edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), HcalCalibration, HcalForward, HcalOther, HcalOuter, hfdigibit_, hfPET_, hfS8S1_, hfS9S1_, HFTimingTrustFlagSetter_, htopoToken_, inputLabel_, mcOOTCorrectionCategory_, mcOOTCorrectionName_, paramsToken_, qualToken_, recoParamsFromDB_, saturationFlagSetter_, setNegativeFlags_, setNoiseFlags_, setPileupCorrection_, setSaturationFlags_, setTimingTrustFlags_, sevToken_, AlCaHLTBitMon_QueryRunRegistry::string, subdet_, HcalZDCDetId::SubdetectorId, subdetOther_, tok_calib_, tok_hf_, tok_ho_, and tsFromDB_.

27  : reco_(conf.getParameter<bool>("correctForTimeslew"),
28  conf.getParameter<bool>("correctForPhaseContainment"),
29  conf.getParameter<double>("correctionPhaseNS"),
30  consumesCollector()),
32  inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
33  correctTiming_(conf.getParameter<bool>("correctTiming")),
34  setNoiseFlags_(conf.getParameter<bool>("setNoiseFlags")),
35  setHSCPFlags_(conf.getParameter<bool>("setHSCPFlags")),
36  setSaturationFlags_(conf.getParameter<bool>("setSaturationFlags")),
37  setTimingTrustFlags_(conf.getParameter<bool>("setTimingTrustFlags")),
38  setPulseShapeFlags_(conf.getParameter<bool>("setPulseShapeFlags")),
39  setNegativeFlags_(false),
40  dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
41  firstAuxTS_(conf.getParameter<int>("firstAuxTS")),
42  firstSample_(conf.getParameter<int>("firstSample")),
43  samplesToAdd_(conf.getParameter<int>("samplesToAdd")),
44  tsFromDB_(conf.getParameter<bool>("tsFromDB")),
45  useLeakCorrection_(conf.getParameter<bool>("useLeakCorrection")),
50  setPileupCorrection_(nullptr) {
51  // register for data access
52  tok_ho_ = consumes<HODigiCollection>(inputLabel_);
53  tok_hf_ = consumes<HFDigiCollection>(inputLabel_);
54  tok_calib_ = consumes<HcalCalibDigiCollection>(inputLabel_);
55 
56  std::string subd = conf.getParameter<std::string>("Subdetector");
57  //Set all FlagSetters to 0
58  /* Important to do this! Otherwise, if the setters are turned off,
59  the "if (XSetter_) delete XSetter_;" commands can crash
60  */
61 
62  recoParamsFromDB_ = conf.getParameter<bool>("recoParamsFromDB");
63  // recoParamsFromDB_ = false ; // trun off for now.
64 
65  // std::cout<<" HcalHitReconstructor recoParamsFromDB_ "<<recoParamsFromDB_<<std::endl;
66 
67  if (conf.existsAs<bool>("setNegativeFlags"))
68  setNegativeFlags_ = conf.getParameter<bool>("setNegativeFlags");
69 
70  hfdigibit_ = nullptr;
71 
72  hfS9S1_ = nullptr;
73  hfS8S1_ = nullptr;
74  hfPET_ = nullptr;
75  saturationFlagSetter_ = nullptr;
76  HFTimingTrustFlagSetter_ = nullptr;
77  digiTimeFromDB_ = false; // only need for HF
78 
79  if (setSaturationFlags_) {
80  const edm::ParameterSet& pssat = conf.getParameter<edm::ParameterSet>("saturationParameters");
81  saturationFlagSetter_ = new HcalADCSaturationFlag(pssat.getParameter<int>("maxADCvalue"));
82  }
83 
84  if (!strcasecmp(subd.c_str(), "HO")) {
86  // setPileupCorrection_ = &HcalSimpleRecAlgo::setHOPileupCorrection;
87  setPileupCorrection_ = nullptr;
88  produces<HORecHitCollection>();
89  } else if (!strcasecmp(subd.c_str(), "HF")) {
91  // setPileupCorrection_ = &HcalSimpleRecAlgo::setHFPileupCorrection;
92  setPileupCorrection_ = nullptr;
93  digiTimeFromDB_ = conf.getParameter<bool>("digiTimeFromDB");
94 
96  const edm::ParameterSet& pstrust = conf.getParameter<edm::ParameterSet>("hfTimingTrustParameters");
97  HFTimingTrustFlagSetter_ = new HFTimingTrustFlag(pstrust.getParameter<int>("hfTimingTrustLevel1"),
98  pstrust.getParameter<int>("hfTimingTrustLevel2"));
99  }
100 
101  if (setNoiseFlags_) {
102  const edm::ParameterSet& psdigi = conf.getParameter<edm::ParameterSet>("digistat");
103  const edm::ParameterSet& psTimeWin = conf.getParameter<edm::ParameterSet>("HFInWindowStat");
104  hfdigibit_ = new HcalHFStatusBitFromDigis(psdigi, psTimeWin);
105 
106  const edm::ParameterSet& psS9S1 = conf.getParameter<edm::ParameterSet>("S9S1stat");
107  hfS9S1_ = new HcalHF_S9S1algorithm(psS9S1.getParameter<std::vector<double> >("short_optimumSlope"),
108  psS9S1.getParameter<std::vector<double> >("shortEnergyParams"),
109  psS9S1.getParameter<std::vector<double> >("shortETParams"),
110  psS9S1.getParameter<std::vector<double> >("long_optimumSlope"),
111  psS9S1.getParameter<std::vector<double> >("longEnergyParams"),
112  psS9S1.getParameter<std::vector<double> >("longETParams"),
113  psS9S1.getParameter<int>("HcalAcceptSeverityLevel"),
114  psS9S1.getParameter<bool>("isS8S1"));
115 
116  const edm::ParameterSet& psS8S1 = conf.getParameter<edm::ParameterSet>("S8S1stat");
117  hfS8S1_ = new HcalHF_S9S1algorithm(psS8S1.getParameter<std::vector<double> >("short_optimumSlope"),
118  psS8S1.getParameter<std::vector<double> >("shortEnergyParams"),
119  psS8S1.getParameter<std::vector<double> >("shortETParams"),
120  psS8S1.getParameter<std::vector<double> >("long_optimumSlope"),
121  psS8S1.getParameter<std::vector<double> >("longEnergyParams"),
122  psS8S1.getParameter<std::vector<double> >("longETParams"),
123  psS8S1.getParameter<int>("HcalAcceptSeverityLevel"),
124  psS8S1.getParameter<bool>("isS8S1"));
125 
126  const edm::ParameterSet& psPET = conf.getParameter<edm::ParameterSet>("PETstat");
127  hfPET_ = new HcalHF_PETalgorithm(psPET.getParameter<std::vector<double> >("short_R"),
128  psPET.getParameter<std::vector<double> >("shortEnergyParams"),
129  psPET.getParameter<std::vector<double> >("shortETParams"),
130  psPET.getParameter<std::vector<double> >("long_R"),
131  psPET.getParameter<std::vector<double> >("longEnergyParams"),
132  psPET.getParameter<std::vector<double> >("longETParams"),
133  psPET.getParameter<int>("HcalAcceptSeverityLevel"),
134  psPET.getParameter<std::vector<double> >("short_R_29"),
135  psPET.getParameter<std::vector<double> >("long_R_29"));
136  }
137  produces<HFRecHitCollection>();
138  } else if (!strcasecmp(subd.c_str(), "ZDC")) {
139  det_ = DetId::Calo;
141  produces<ZDCRecHitCollection>();
142  } else if (!strcasecmp(subd.c_str(), "CALIB")) {
143  subdet_ = HcalOther;
145  produces<HcalCalibRecHitCollection>();
146  } else {
147  edm::LogWarning("Configuration") << "HcalHitReconstructor is not associated with a specific subdetector!"
148  << std::endl;
149  }
150 
151  // If no valid OOT pileup correction name specified,
152  // disable the correction
153  if (conf.existsAs<std::string>("dataOOTCorrectionName"))
154  dataOOTCorrectionName_ = conf.getParameter<std::string>("dataOOTCorrectionName");
155  if (conf.existsAs<std::string>("dataOOTCorrectionCategory"))
156  dataOOTCorrectionCategory_ = conf.getParameter<std::string>("dataOOTCorrectionCategory");
157  if (conf.existsAs<std::string>("mcOOTCorrectionName"))
158  mcOOTCorrectionName_ = conf.getParameter<std::string>("mcOOTCorrectionName");
159  if (conf.existsAs<std::string>("mcOOTCorrectionCategory"))
160  mcOOTCorrectionCategory_ = conf.getParameter<std::string>("mcOOTCorrectionCategory");
161  if (dataOOTCorrectionName_.empty() && mcOOTCorrectionName_.empty())
162  setPileupCorrection_ = nullptr;
163 
164  // ES tokens
165  htopoToken_ = esConsumes<HcalTopology, HcalRecNumberingRecord, edm::Transition::BeginRun>();
167  paramsToken_ = esConsumes<HcalRecoParams, HcalRecoParamsRcd, edm::Transition::BeginRun>();
168  if (digiTimeFromDB_)
169  digiTimeToken_ = esConsumes<HcalFlagHFDigiTimeParams, HcalFlagHFDigiTimeParamsRcd, edm::Transition::BeginRun>();
170  conditionsToken_ = esConsumes<HcalDbService, HcalDbRecord>();
171  qualToken_ = esConsumes<HcalChannelQuality, HcalChannelQualityRcd>(edm::ESInputTag("", "withTopo"));
172  sevToken_ = esConsumes<HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd>();
173 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
HcalADCSaturationFlag * saturationFlagSetter_
edm::ESGetToken< HcalDbService, HcalDbRecord > conditionsToken_
SetCorrectionFcn setPileupCorrection_
edm::ESGetToken< HcalChannelQuality, HcalChannelQualityRcd > qualToken_
edm::ESGetToken< HcalFlagHFDigiTimeParams, HcalFlagHFDigiTimeParamsRcd > digiTimeToken_
std::string dataOOTCorrectionCategory_
edm::EDGetTokenT< HcalCalibDigiCollection > tok_calib_
HcalHFStatusBitFromDigis * hfdigibit_
edm::EDGetTokenT< HODigiCollection > tok_ho_
HcalHF_PETalgorithm * hfPET_
HcalHF_S9S1algorithm * hfS9S1_
edm::ESGetToken< HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd > sevToken_
static const int SubdetectorId
Definition: HcalZDCDetId.h:25
HcalHF_S9S1algorithm * hfS8S1_
HcalOtherSubdetector subdetOther_
Log< level::Warning, false > LogWarning
HFTimingTrustFlag * HFTimingTrustFlagSetter_
edm::ESGetToken< HcalRecoParams, HcalRecoParamsRcd > paramsToken_
std::string mcOOTCorrectionCategory_
edm::EDGetTokenT< HFDigiCollection > tok_hf_
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > htopoToken_
HcalSimpleRecAlgo reco_

◆ ~HcalHitReconstructor()

HcalHitReconstructor::~HcalHitReconstructor ( )
override

Definition at line 175 of file HcalHitReconstructor.cc.

References hfdigibit_, hfPET_, hfS8S1_, hfS9S1_, HFTimingTrustFlagSetter_, and saturationFlagSetter_.

175  {
176  delete hfdigibit_;
177 
178  delete hfS9S1_;
179  delete hfS8S1_;
180  delete hfPET_;
181  delete saturationFlagSetter_;
183 }
HcalADCSaturationFlag * saturationFlagSetter_
HcalHFStatusBitFromDigis * hfdigibit_
HcalHF_PETalgorithm * hfPET_
HcalHF_S9S1algorithm * hfS9S1_
HcalHF_S9S1algorithm * hfS8S1_
HFTimingTrustFlag * HFTimingTrustFlagSetter_

Member Function Documentation

◆ beginRun()

void HcalHitReconstructor::beginRun ( edm::Run const &  r,
edm::EventSetup const &  es 
)
final

Definition at line 185 of file HcalHitReconstructor.cc.

References HcalSimpleRecAlgo::beginRun(), digiTimeFromDB_, digiTimeToken_, edm::EventSetup::getData(), hFDigiTimeParams_, htopoToken_, AlCaHLTBitMon_ParallelJobs::p, paramsToken_, paramTS_, reco_, recoParamsFromDB_, and tsFromDB_.

185  {
186  const HcalTopology& htopo = es.getData(htopoToken_);
187 
188  if (tsFromDB_ || recoParamsFromDB_) {
189  const HcalRecoParams& p = es.getData(paramsToken_);
190  paramTS_ = std::make_unique<HcalRecoParams>(p);
191  paramTS_->setTopo(&htopo);
192 
193  // std::cout<<" skdump in HcalHitReconstructor::beginRun dupm RecoParams "<<std::endl;
194  // std::ofstream skfile("skdumpRecoParamsNewFormat.txt");
195  // HcalDbASCIIIO::dumpObject(skfile, (*paramTS_) );
196  }
197 
198  if (digiTimeFromDB_) {
199  const HcalFlagHFDigiTimeParams& p = es.getData(digiTimeToken_);
200  hFDigiTimeParams_ = std::make_unique<HcalFlagHFDigiTimeParams>(p);
201  hFDigiTimeParams_->setTopo(&htopo);
202  }
203 
204  reco_.beginRun(es);
205 }
void beginRun(edm::EventSetup const &es)
std::unique_ptr< HcalRecoParams > paramTS_
edm::ESGetToken< HcalFlagHFDigiTimeParams, HcalFlagHFDigiTimeParamsRcd > digiTimeToken_
std::unique_ptr< HcalFlagHFDigiTimeParams > hFDigiTimeParams_
edm::ESGetToken< HcalRecoParams, HcalRecoParamsRcd > paramsToken_
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > htopoToken_
HcalSimpleRecAlgo reco_

◆ endRun()

void HcalHitReconstructor::endRun ( edm::Run const &  r,
edm::EventSetup const &  es 
)
final

Definition at line 207 of file HcalHitReconstructor.cc.

References HcalSimpleRecAlgo::endRun(), and reco_.

207 { reco_.endRun(); }
HcalSimpleRecAlgo reco_

◆ produce()

void HcalHitReconstructor::produce ( edm::Event e,
const edm::EventSetup c 
)
override

Definition at line 209 of file HcalHitReconstructor.cc.

References funct::abs(), edm::SortedCollection< T, SORT >::begin(), AlignmentProducer_cff::calibrations, submitPVValidationJobs::conditions, conditionsToken_, HcalTimingCorrector::Correct(), HcalRecoParam::correctForPhaseContainment(), castor_dqm_sourceclient-live_cfg::correctForPhaseContainment, HcalRecoParam::correctForTimeslew(), castor_dqm_sourceclient-live_cfg::correctForTimeslew, HcalRecoParam::correctionPhaseNS(), HcalRecoParam::correctTiming(), correctTiming_, hcalRecHitTable_cff::depth, det_, digiTimeFromDB_, HcalSeverityLevelComputer::dropChannel(), dropZSmarkedPassed_, MillePedeFileConverter_cfg::e, edm::SortedCollection< T, SORT >::end(), options_cfi::eventSetup, dqmdumpme::first, HcalRecoParam::firstAuxTS(), firstAuxTS_, HcalRecoParam::firstSample(), firstSample_, HcalChannelStatus::getValue(), HcalCondObjectContainer< Item >::getValues(), DetId::Hcal, HcalCalibration, HcalForward, HcalOther, HcalOuter, hfdigibit_, HcalFlagHFDigiTimeParam::HFdigiflagCoefficients(), HcalFlagHFDigiTimeParam::HFdigiflagExpectedPeak(), HcalFlagHFDigiTimeParam::HFdigiflagFirstSample(), HcalFlagHFDigiTimeParam::HFdigiflagMinEThreshold(), HcalFlagHFDigiTimeParam::HFdigiflagSamplesToAdd(), hFDigiTimeParams_, hfPET_, hfS8S1_, hfS9S1_, HcalHFStatusBitFromDigis::hfSetFlagFromDigi(), HcalHF_PETalgorithm::HFSetFlagFromPET(), HcalHF_S9S1algorithm::HFSetFlagFromS9S1(), HFTimingTrustFlagSetter_, mps_fire::i, hcalRecHitTable_cff::ieta, dqmdumpme::k, eostools::move(), paramTS_, HcalRecoParam::pileupCleaningID(), HcalCaloFlagLabels::PresampleADC, qualToken_, DetId::rawId(), reco_, HcalSimpleRecAlgo::reconstruct(), recoParamsFromDB_, HcalHFStatusBitFromDigis::resetParamsFromDB(), HcalRecoParam::samplesToAdd(), samplesToAdd_, saturationFlagSetter_, HFTimingTrustFlag::setHFTimingTrustFlag(), HcalSimpleRecAlgo::setLeakCorrection(), setNoiseFlags_, HcalSimpleRecAlgo::setRecoParams(), HcalADCSaturationFlag::setSaturationFlag(), setSaturationFlags_, setTimingTrustFlags_, sevToken_, l1trig_cff::shape, edm::SortedCollection< T, SORT >::size(), subdet_, subdetOther_, tok_calib_, tok_hf_, tok_ho_, tsFromDB_, HcalRecoParam::useLeakCorrection(), useLeakCorrection_, and geometryCSVtoXML::xx.

209  {
210  // get conditions
212  const HcalChannelQuality* myqual = &eventSetup.getData(qualToken_);
213  const HcalSeverityLevelComputer* mySeverity = &eventSetup.getData(sevToken_);
214 
215  if (useLeakCorrection_)
217 
218  // GET THE BEAM CROSSING INFO HERE, WHEN WE UNDERSTAND HOW THINGS WORK.
219  // Then, call "setBXInfo" method of the reco_ object.
220  // Also remember to call SetBXInfo in the negative energy flag setter.
221 
222  if (det_ == DetId::Hcal) {
223  // HO ------------------------------------------------------------------
224  if (subdet_ == HcalOuter) {
226  e.getByToken(tok_ho_, digi);
227 
228  // create empty output
229  auto rec = std::make_unique<HORecHitCollection>();
230  rec->reserve(digi->size());
231  // run the algorithm
233 
234  // Vote on majority TS0 CapId
235  int favorite_capid = 0;
236  if (correctTiming_) {
237  long capid_votes[4] = {0, 0, 0, 0};
238  for (i = digi->begin(); i != digi->end(); i++) {
239  capid_votes[(*i)[0].capid()]++;
240  }
241  for (int k = 0; k < 4; k++)
242  if (capid_votes[k] > capid_votes[favorite_capid])
243  favorite_capid = k;
244  }
245 
246  for (i = digi->begin(); i != digi->end(); i++) {
247  HcalDetId cell = i->id();
248  DetId detcell = (DetId)cell;
249  // firstSample & samplesToAdd
250  if (tsFromDB_ || recoParamsFromDB_) {
251  const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
252  if (tsFromDB_) {
253  firstSample_ = param_ts->firstSample();
254  samplesToAdd_ = param_ts->samplesToAdd();
255  }
256  if (recoParamsFromDB_) {
257  bool correctForTimeslew = param_ts->correctForTimeslew();
259  float phaseNS = param_ts->correctionPhaseNS();
261  correctTiming_ = param_ts->correctTiming();
262  firstAuxTS_ = param_ts->firstAuxTS();
263  int pileupCleaningID = param_ts->pileupCleaningID();
266  }
267  }
268 
269  int first = firstSample_;
270  int toadd = samplesToAdd_;
271 
272  if (first >= i->size() || first < 0)
273  edm::LogWarning("Configuration")
274  << "HcalHitReconstructor: illegal firstSample" << first << " in subdet " << subdet_ << std::endl;
275 
276  // check on cells to be ignored and dropped: (rof,20.Feb.09)
277  const HcalChannelStatus* mydigistatus = myqual->getValues(detcell.rawId());
278  if (mySeverity->dropChannel(mydigistatus->getValue()))
279  continue;
281  if (i->zsMarkAndPass())
282  continue;
283 
284  const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
285  const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
286  const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
287  HcalCoderDb coder(*channelCoder, *shape);
288 
289  rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
290 
291  // Set auxiliary flag
292  int auxflag = 0;
293  int fTS = firstAuxTS_;
294  if (fTS < 0)
295  fTS = 0; //silly protection against negative time slice values
296  for (int xx = fTS; xx < fTS + 4 && xx < i->size(); ++xx)
297  auxflag += (i->sample(xx).adc())
298  << (7 *
299  (xx - fTS)); // store the time slices in the first 28 bits of aux, a set of 4 7-bit adc values
300  // bits 28 and 29 are reserved for capid of the first time slice saved in aux
301  auxflag += ((i->sample(fTS).capid()) << 28);
302  (rec->back()).setAux(auxflag);
303  // (rec->back()).setFlags(0); Don't want to do this because the algorithm
304  // can already set some flags
305  // Fill Presample ADC flag
306  if (fTS > 0)
307  (rec->back()).setFlagField((i->sample(fTS - 1).adc()), HcalCaloFlagLabels::PresampleADC, 7);
308 
310  saturationFlagSetter_->setSaturationFlag(rec->back(), *i);
311  if (correctTiming_)
312  HcalTimingCorrector::Correct(rec->back(), *i, favorite_capid);
313  }
314  // return result
315  e.put(std::move(rec));
316 
317  // HF -------------------------------------------------------------------
318  } else if (subdet_ == HcalForward) {
320  e.getByToken(tok_hf_, digi);
321 
323  // create empty output
324  auto rec = std::make_unique<HFRecHitCollection>();
325  rec->reserve(digi->size());
326  // run the algorithm
328 
329  // Vote on majority TS0 CapId
330  int favorite_capid = 0;
331  if (correctTiming_) {
332  long capid_votes[4] = {0, 0, 0, 0};
333  for (i = digi->begin(); i != digi->end(); i++) {
334  capid_votes[(*i)[0].capid()]++;
335  }
336  for (int k = 0; k < 4; k++)
337  if (capid_votes[k] > capid_votes[favorite_capid])
338  favorite_capid = k;
339  }
340 
341  for (i = digi->begin(); i != digi->end(); i++) {
342  HcalDetId cell = i->id();
343  DetId detcell = (DetId)cell;
344 
345  if (tsFromDB_ || recoParamsFromDB_) {
346  const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
347  if (tsFromDB_) {
348  firstSample_ = param_ts->firstSample();
349  samplesToAdd_ = param_ts->samplesToAdd();
350  }
351  if (recoParamsFromDB_) {
352  bool correctForTimeslew = param_ts->correctForTimeslew();
354  float phaseNS = param_ts->correctionPhaseNS();
356  correctTiming_ = param_ts->correctTiming();
357  firstAuxTS_ = param_ts->firstAuxTS();
358  int pileupCleaningID = param_ts->pileupCleaningID();
361  }
362  }
363 
364  int first = firstSample_;
365  int toadd = samplesToAdd_;
366 
367  if (first >= i->size() || first < 0)
368  edm::LogWarning("Configuration")
369  << "HcalHitReconstructor: illegal firstSample" << first << " in subdet " << subdet_ << std::endl;
370 
371  // check on cells to be ignored and dropped: (rof,20.Feb.09)
372  const HcalChannelStatus* mydigistatus = myqual->getValues(detcell.rawId());
373  if (mySeverity->dropChannel(mydigistatus->getValue()))
374  continue;
376  if (i->zsMarkAndPass())
377  continue;
378 
379  const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
380  const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
381  const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
382  HcalCoderDb coder(*channelCoder, *shape);
383 
384  // Set HFDigiTime flag values from digiTimeFromDB_
385  if (digiTimeFromDB_ && hfdigibit_) {
386  const HcalFlagHFDigiTimeParam* hfDTparam = hFDigiTimeParams_->getValues(detcell.rawId());
388  hfDTparam->HFdigiflagSamplesToAdd(),
389  hfDTparam->HFdigiflagExpectedPeak(),
390  hfDTparam->HFdigiflagMinEThreshold(),
391  hfDTparam->HFdigiflagCoefficients());
392  }
393 
394  //std::cout << "TOADDHF " << toadd << " " << first << " " << std::endl;
395  rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
396 
397  // Set auxiliary flag
398  int auxflag = 0;
399  int fTS = firstAuxTS_;
400  if (fTS < 0)
401  fTS = 0; // silly protection against negative time slice values
402  for (int xx = fTS; xx < fTS + 4 && xx < i->size(); ++xx)
403  auxflag += (i->sample(xx).adc())
404  << (7 *
405  (xx - fTS)); // store the time slices in the first 28 bits of aux, a set of 4 7-bit adc values
406  // bits 28 and 29 are reserved for capid of the first time slice saved in aux
407  auxflag += ((i->sample(fTS).capid()) << 28);
408  (rec->back()).setAux(auxflag);
409 
410  // (rec->back()).setFlags(0); Don't want to do this because the algorithm
411  // can already set some flags
412 
413  // Fill Presample ADC flag
414  if (fTS > 0)
415  (rec->back()).setFlagField((i->sample(fTS - 1).adc()), HcalCaloFlagLabels::PresampleADC, 7);
416 
417  // This calls the code for setting the HF noise bit determined from digi shape
418  if (setNoiseFlags_)
419  hfdigibit_->hfSetFlagFromDigi(rec->back(), *i, coder, calibrations);
421  saturationFlagSetter_->setSaturationFlag(rec->back(), *i);
424  if (correctTiming_)
425  HcalTimingCorrector::Correct(rec->back(), *i, favorite_capid);
426  } // for (i=digi->begin(); i!=digi->end(); i++) -- loop on all HF digis
427 
428  // The following flags require the full set of rechits
429  // These need to be set consecutively, so an energy check should be the first
430  // test performed on these hits (to minimize the loop time)
431  if (setNoiseFlags_) {
432  // Step 1: Set PET flag (short fibers of |ieta|==29)
433  // Neighbor/partner channels that are flagged by Pulse Shape algorithm (HFDigiTime)
434  // won't be considered in these calculations
435  for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
436  int depth = i->id().depth();
437  int ieta = i->id().ieta();
438  // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
439  if (depth == 2 || abs(ieta) == 29)
440  hfPET_->HFSetFlagFromPET(*i, *rec, myqual, mySeverity);
441  }
442 
443  // Step 2: Set S8S1 flag (short fibers or |ieta|==29)
444  for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
445  int depth = i->id().depth();
446  int ieta = i->id().ieta();
447  // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
448  if (depth == 2 || abs(ieta) == 29)
449  hfS8S1_->HFSetFlagFromS9S1(*i, *rec, myqual, mySeverity);
450  }
451 
452  // Set 3: Set S9S1 flag (long fibers)
453  for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
454  int depth = i->id().depth();
455  int ieta = i->id().ieta();
456  // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
457  if (depth == 1 && abs(ieta) != 29)
458  hfS9S1_->HFSetFlagFromS9S1(*i, *rec, myqual, mySeverity);
459  }
460  }
461 
462  // return result
463  e.put(std::move(rec));
464  } else if (subdet_ == HcalOther && subdetOther_ == HcalCalibration) {
466  e.getByToken(tok_calib_, digi);
467 
468  // create empty output
469  auto rec = std::make_unique<HcalCalibRecHitCollection>();
470  rec->reserve(digi->size());
471  // run the algorithm
472  int first = firstSample_;
473  int toadd = samplesToAdd_;
474 
476  for (i = digi->begin(); i != digi->end(); i++) {
477  HcalCalibDetId cell = i->id();
478  // HcalDetId cellh = i->id();
479  DetId detcell = (DetId)cell;
480  // check on cells to be ignored and dropped: (rof,20.Feb.09)
481  const HcalChannelStatus* mydigistatus = myqual->getValues(detcell.rawId());
482  if (mySeverity->dropChannel(mydigistatus->getValue()))
483  continue;
485  if (i->zsMarkAndPass())
486  continue;
487 
488  const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
489  const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
490  const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
491  HcalCoderDb coder(*channelCoder, *shape);
492 
493  // firstSample & samplesToAdd
494  if (tsFromDB_) {
495  const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
496  first = param_ts->firstSample();
497  toadd = param_ts->samplesToAdd();
498  }
499  rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
500 
501  /*
502  // Flag setting not available for calibration rechits
503  // Set auxiliary flag
504  int auxflag=0;
505  int fTS = firstAuxTS_;
506  for (int xx=fTS; xx<fTS+4 && xx<i->size();++xx)
507  auxflag+=(i->sample(xx).adc())<<(7*(xx-fTS)); // store the time slices in the first 28 bits of aux, a set of 4 7-bit adc values
508  // bits 28 and 29 are reserved for capid of the first time slice saved in aux
509  auxflag+=((i->sample(fTS).capid())<<28);
510  (rec->back()).setAux(auxflag);
511 
512  (rec->back()).setFlags(0); // Not yet implemented for HcalCalibRecHit
513  */
514  }
515  // return result
516  e.put(std::move(rec));
517  }
518  }
519  //DL delete myqual;
520 } // void HcalHitReconstructor::produce(...)
std::vector< double > HFdigiflagCoefficients() const
void HFSetFlagFromS9S1(HFRecHit &hf, HFRecHitCollection &rec, const HcalChannelQuality *myqual, const HcalSeverityLevelComputer *mySeverity)
constexpr unsigned int firstAuxTS() const
Definition: HcalRecoParam.h:40
void setHFTimingTrustFlag(HFRecHit &rechit, const HFDataFrame &digi)
HcalADCSaturationFlag * saturationFlagSetter_
size_type size() const
std::unique_ptr< HcalRecoParams > paramTS_
void hfSetFlagFromDigi(HFRecHit &hf, const HFDataFrame &digi, const HcalCoder &coder, const HcalCalibrations &calib)
constexpr unsigned int firstSample() const
Definition: HcalRecoParam.h:31
std::vector< T >::const_iterator const_iterator
uint32_t HFdigiflagExpectedPeak() const
void resetParamsFromDB(int firstSample, int samplesToAdd, int expectedPeak, double minthreshold, const std::vector< double > &coef)
edm::ESGetToken< HcalDbService, HcalDbRecord > conditionsToken_
constexpr bool correctForPhaseContainment() const
Definition: HcalRecoParam.h:28
const Item * getValues(DetId fId, bool throwOnFail=true) const
edm::ESGetToken< HcalChannelQuality, HcalChannelQualityRcd > qualToken_
void HFSetFlagFromPET(HFRecHit &hf, HFRecHitCollection &rec, const HcalChannelQuality *myqual, const HcalSeverityLevelComputer *mySeverity)
constexpr bool correctForTimeslew() const
Definition: HcalRecoParam.h:37
std::unique_ptr< HcalFlagHFDigiTimeParams > hFDigiTimeParams_
edm::EDGetTokenT< HcalCalibDigiCollection > tok_calib_
constexpr bool correctTiming() const
Definition: HcalRecoParam.h:39
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HcalHFStatusBitFromDigis * hfdigibit_
const_iterator begin() const
void setSaturationFlag(HBHERecHit &rechit, const HBHEDataFrame &digi)
constexpr float correctionPhaseNS() const
Definition: HcalRecoParam.h:30
edm::EDGetTokenT< HODigiCollection > tok_ho_
uint32_t getValue() const
std::vector< T >::iterator iterator
const_iterator end() const
HcalHF_PETalgorithm * hfPET_
Definition: DetId.h:17
HcalHF_S9S1algorithm * hfS9S1_
edm::ESGetToken< HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd > sevToken_
void setRecoParams(bool correctForTimeslew, bool correctForPulse, bool setLeakCorrection, int pileupCleaningID, float phaseNS)
constexpr bool useLeakCorrection() const
Definition: HcalRecoParam.h:35
uint32_t HFdigiflagFirstSample() const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
constexpr unsigned int samplesToAdd() const
Definition: HcalRecoParam.h:32
HcalHF_S9S1algorithm * hfS8S1_
HcalOtherSubdetector subdetOther_
static void Correct(HBHERecHit &rechit, const HBHEDataFrame &digi, int favorite_capid)
constexpr unsigned int pileupCleaningID() const
Definition: HcalRecoParam.h:43
HFRecHit reconstruct(const HFDataFrame &digi, int first, int toadd, const HcalCoder &coder, const HcalCalibrations &calibs) const
HFTimingTrustFlag * HFTimingTrustFlagSetter_
edm::EDGetTokenT< HFDigiCollection > tok_hf_
def move(src, dest)
Definition: eostools.py:511
HcalSimpleRecAlgo reco_
uint32_t HFdigiflagSamplesToAdd() const
bool dropChannel(const uint32_t &mystatus) const

Member Data Documentation

◆ cat_

std::string HcalHitReconstructor::cat_
private

Definition at line 109 of file HcalHitReconstructor.h.

◆ conditionsToken_

edm::ESGetToken<HcalDbService, HcalDbRecord> HcalHitReconstructor::conditionsToken_
private

Definition at line 115 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ correctTiming_

bool HcalHitReconstructor::correctTiming_
private

Definition at line 78 of file HcalHitReconstructor.h.

Referenced by produce().

◆ corrName_

std::string HcalHitReconstructor::corrName_
private

Definition at line 109 of file HcalHitReconstructor.h.

◆ dataOOTCorrectionCategory_

std::string HcalHitReconstructor::dataOOTCorrectionCategory_
private

Definition at line 101 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor().

◆ dataOOTCorrectionName_

std::string HcalHitReconstructor::dataOOTCorrectionName_
private

Definition at line 100 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor().

◆ det_

DetId::Detector HcalHitReconstructor::det_
private

Definition at line 70 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ digiTimeFromDB_

bool HcalHitReconstructor::digiTimeFromDB_
private

Definition at line 94 of file HcalHitReconstructor.h.

Referenced by beginRun(), HcalHitReconstructor(), and produce().

◆ digiTimeToken_

edm::ESGetToken<HcalFlagHFDigiTimeParams, HcalFlagHFDigiTimeParamsRcd> HcalHitReconstructor::digiTimeToken_
private

Definition at line 114 of file HcalHitReconstructor.h.

Referenced by beginRun(), and HcalHitReconstructor().

◆ dropZSmarkedPassed_

bool HcalHitReconstructor::dropZSmarkedPassed_
private

Definition at line 85 of file HcalHitReconstructor.h.

Referenced by produce().

◆ firstAuxTS_

int HcalHitReconstructor::firstAuxTS_
private

Definition at line 87 of file HcalHitReconstructor.h.

Referenced by produce().

◆ firstSample_

int HcalHitReconstructor::firstSample_
private

Definition at line 90 of file HcalHitReconstructor.h.

Referenced by produce().

◆ hfdigibit_

HcalHFStatusBitFromDigis* HcalHitReconstructor::hfdigibit_
private

Definition at line 65 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), produce(), and ~HcalHitReconstructor().

◆ hFDigiTimeParams_

std::unique_ptr<HcalFlagHFDigiTimeParams> HcalHitReconstructor::hFDigiTimeParams_
private

Definition at line 107 of file HcalHitReconstructor.h.

Referenced by beginRun(), and produce().

◆ hfPET_

HcalHF_PETalgorithm* HcalHitReconstructor::hfPET_
private

Definition at line 68 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), produce(), and ~HcalHitReconstructor().

◆ hfS8S1_

HcalHF_S9S1algorithm* HcalHitReconstructor::hfS8S1_
private

Definition at line 67 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), produce(), and ~HcalHitReconstructor().

◆ hfS9S1_

HcalHF_S9S1algorithm* HcalHitReconstructor::hfS9S1_
private

Definition at line 66 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), produce(), and ~HcalHitReconstructor().

◆ HFTimingTrustFlagSetter_

HFTimingTrustFlag* HcalHitReconstructor::HFTimingTrustFlagSetter_
private

Definition at line 64 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), produce(), and ~HcalHitReconstructor().

◆ htopoToken_

edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> HcalHitReconstructor::htopoToken_
private

Definition at line 112 of file HcalHitReconstructor.h.

Referenced by beginRun(), and HcalHitReconstructor().

◆ inputLabel_

edm::InputTag HcalHitReconstructor::inputLabel_
private

Definition at line 73 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor().

◆ mcOOTCorrectionCategory_

std::string HcalHitReconstructor::mcOOTCorrectionCategory_
private

Definition at line 103 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor().

◆ mcOOTCorrectionName_

std::string HcalHitReconstructor::mcOOTCorrectionName_
private

Definition at line 102 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor().

◆ paramsToken_

edm::ESGetToken<HcalRecoParams, HcalRecoParamsRcd> HcalHitReconstructor::paramsToken_
private

Definition at line 113 of file HcalHitReconstructor.h.

Referenced by beginRun(), and HcalHitReconstructor().

◆ paramTS_

std::unique_ptr<HcalRecoParams> HcalHitReconstructor::paramTS_
private

Definition at line 106 of file HcalHitReconstructor.h.

Referenced by beginRun(), and produce().

◆ qualToken_

edm::ESGetToken<HcalChannelQuality, HcalChannelQualityRcd> HcalHitReconstructor::qualToken_
private

Definition at line 116 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ reco_

HcalSimpleRecAlgo HcalHitReconstructor::reco_
private

Definition at line 62 of file HcalHitReconstructor.h.

Referenced by beginRun(), endRun(), and produce().

◆ recoParamsFromDB_

bool HcalHitReconstructor::recoParamsFromDB_
private

Definition at line 93 of file HcalHitReconstructor.h.

Referenced by beginRun(), HcalHitReconstructor(), and produce().

◆ samplesToAdd_

int HcalHitReconstructor::samplesToAdd_
private

Definition at line 91 of file HcalHitReconstructor.h.

Referenced by produce().

◆ saturationFlagSetter_

HcalADCSaturationFlag* HcalHitReconstructor::saturationFlagSetter_
private

Definition at line 63 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), produce(), and ~HcalHitReconstructor().

◆ setHSCPFlags_

bool HcalHitReconstructor::setHSCPFlags_
private

Definition at line 80 of file HcalHitReconstructor.h.

◆ setNegativeFlags_

bool HcalHitReconstructor::setNegativeFlags_
private

Definition at line 84 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor().

◆ setNoiseFlags_

bool HcalHitReconstructor::setNoiseFlags_
private

Definition at line 79 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ setPileupCorrection_

SetCorrectionFcn HcalHitReconstructor::setPileupCorrection_
private

Definition at line 104 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor().

◆ setPulseShapeFlags_

bool HcalHitReconstructor::setPulseShapeFlags_
private

Definition at line 83 of file HcalHitReconstructor.h.

◆ setSaturationFlags_

bool HcalHitReconstructor::setSaturationFlags_
private

Definition at line 81 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ setTimingTrustFlags_

bool HcalHitReconstructor::setTimingTrustFlags_
private

Definition at line 82 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ sevToken_

edm::ESGetToken<HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd> HcalHitReconstructor::sevToken_
private

Definition at line 117 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ subdet_

int HcalHitReconstructor::subdet_
private

Definition at line 71 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ subdetOther_

HcalOtherSubdetector HcalHitReconstructor::subdetOther_
private

Definition at line 72 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ tok_calib_

edm::EDGetTokenT<HcalCalibDigiCollection> HcalHitReconstructor::tok_calib_
private

Definition at line 76 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ tok_hf_

edm::EDGetTokenT<HFDigiCollection> HcalHitReconstructor::tok_hf_
private

Definition at line 75 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ tok_ho_

edm::EDGetTokenT<HODigiCollection> HcalHitReconstructor::tok_ho_
private

Definition at line 74 of file HcalHitReconstructor.h.

Referenced by HcalHitReconstructor(), and produce().

◆ tsFromDB_

bool HcalHitReconstructor::tsFromDB_
private

Definition at line 92 of file HcalHitReconstructor.h.

Referenced by beginRun(), HcalHitReconstructor(), and produce().

◆ useLeakCorrection_

bool HcalHitReconstructor::useLeakCorrection_
private

Definition at line 97 of file HcalHitReconstructor.h.

Referenced by produce().