CMS 3D CMS Logo

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

#include <HcalSimpleReconstructor.h>

Inheritance diagram for HcalSimpleReconstructor:
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
 
 HcalSimpleReconstructor (const edm::ParameterSet &ps)
 
void produce (edm::Event &e, const edm::EventSetup &c) final
 
 ~HcalSimpleReconstructor () 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
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

template<class DIGICOLL , class RECHITCOLL >
void process (edm::Event &e, const edm::EventSetup &c, const edm::EDGetTokenT< DIGICOLL > &tok)
 

Private Attributes

edm::ESGetToken< HcalDbService, HcalDbRecordconditionsToken_
 
DetId::Detector det_
 
bool dropZSmarkedPassed_
 
int firstSample_
 
edm::ESGetToken< HcalTopology, HcalRecNumberingRecordhtopoToken_
 
edm::InputTag inputLabel_
 
edm::ESGetToken< HcalRecoParams, HcalRecoParamsRcdparamsToken_
 
std::unique_ptr< HcalRecoParamsparamTS_
 
HcalSimpleRecAlgo reco_
 
int samplesToAdd_
 
int subdet_
 
HcalOtherSubdetector subdetOther_
 
edm::EDGetTokenT< HcalCalibDigiCollectiontok_calib_
 
edm::EDGetTokenT< HFDigiCollectiontok_hf_
 
edm::EDGetTokenT< HODigiCollectiontok_ho_
 
bool tsFromDB_
 

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. Mans - Minnesota
E. Garcia - CSU Based on HcalSimpleReconstructor.h by J. Mans

Definition at line 33 of file HcalSimpleReconstructor.h.

Constructor & Destructor Documentation

◆ HcalSimpleReconstructor()

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

Definition at line 16 of file HcalSimpleReconstructor.cc.

References conditionsToken_, gather_cfg::cout, edm::ParameterSet::getParameter(), HcalForward, HcalOuter, htopoToken_, inputLabel_, paramsToken_, AlCaHLTBitMon_QueryRunRegistry::string, subdet_, tok_calib_, tok_hf_, tok_ho_, and tsFromDB_.

17  : reco_(conf.getParameter<bool>("correctForTimeslew"),
18  conf.getParameter<bool>("correctForPhaseContainment"),
19  conf.getParameter<double>("correctionPhaseNS"),
20  consumesCollector()),
22  inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
23  dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
24  firstSample_(conf.getParameter<int>("firstSample")),
25  samplesToAdd_(conf.getParameter<int>("samplesToAdd")),
26  tsFromDB_(conf.getParameter<bool>("tsFromDB")) {
27  // register for data access
28  tok_hf_ = consumes<HFDigiCollection>(inputLabel_);
29  tok_ho_ = consumes<HODigiCollection>(inputLabel_);
30  tok_calib_ = consumes<HcalCalibDigiCollection>(inputLabel_);
31 
32  std::string subd = conf.getParameter<std::string>("Subdetector");
33  if (!strcasecmp(subd.c_str(), "HO")) {
35  produces<HORecHitCollection>();
36  } else if (!strcasecmp(subd.c_str(), "HF")) {
38  produces<HFRecHitCollection>();
39  } else {
40  std::cout << "HcalSimpleReconstructor is not associated with a specific subdetector!" << std::endl;
41  }
42 
43  // ES tokens
44  if (tsFromDB_) {
45  htopoToken_ = esConsumes<HcalTopology, HcalRecNumberingRecord, edm::Transition::BeginRun>();
46  paramsToken_ = esConsumes<HcalRecoParams, HcalRecoParamsRcd, edm::Transition::BeginRun>();
47  }
48  conditionsToken_ = esConsumes<HcalDbService, HcalDbRecord>();
49 }
edm::ESGetToken< HcalRecoParams, HcalRecoParamsRcd > paramsToken_
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > htopoToken_
edm::ESGetToken< HcalDbService, HcalDbRecord > conditionsToken_
edm::EDGetTokenT< HFDigiCollection > tok_hf_
edm::EDGetTokenT< HODigiCollection > tok_ho_
edm::EDGetTokenT< HcalCalibDigiCollection > tok_calib_

◆ ~HcalSimpleReconstructor()

HcalSimpleReconstructor::~HcalSimpleReconstructor ( )
override

Definition at line 51 of file HcalSimpleReconstructor.cc.

51 {}

Member Function Documentation

◆ beginRun()

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

Definition at line 53 of file HcalSimpleReconstructor.cc.

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

53  {
54  if (tsFromDB_) {
55  const HcalTopology& htopo = es.getData(htopoToken_);
56  const HcalRecoParams& p = es.getData(paramsToken_);
57  paramTS_ = std::make_unique<HcalRecoParams>(p);
58  paramTS_->setTopo(&htopo);
59  }
60  reco_.beginRun(es);
61 }
std::unique_ptr< HcalRecoParams > paramTS_
edm::ESGetToken< HcalRecoParams, HcalRecoParamsRcd > paramsToken_
void beginRun(edm::EventSetup const &es)
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > htopoToken_

◆ endRun()

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

Definition at line 63 of file HcalSimpleReconstructor.cc.

References HcalSimpleRecAlgo::endRun(), and reco_.

◆ fillDescriptions()

void HcalSimpleReconstructor::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 119 of file HcalSimpleReconstructor.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), HLT_2022v15_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

119  {
120  // horeco
122  descHO.add<double>("correctionPhaseNS", 13.0);
123  descHO.add<edm::InputTag>("digiLabel", edm::InputTag("hcalDigis"));
124  descHO.add<bool>("tsFromDB", true);
125  descHO.add<int>("samplesToAdd", 4);
126  descHO.add<std::string>("Subdetector", "HO");
127  descHO.add<bool>("correctForTimeslew", true);
128  descHO.add<bool>("dropZSmarkedPassed", true);
129  descHO.add<bool>("correctForPhaseContainment", true);
130  descHO.add<int>("firstSample", 4);
131  descriptions.add("hosimplereco", descHO);
132 
133  // hfreco
135  descHF.add<double>("correctionPhaseNS", 0.0);
136  descHF.add<edm::InputTag>("digiLabel", edm::InputTag("hcalDigis"));
137  descHF.add<bool>("tsFromDB", true);
138  descHF.add<int>("samplesToAdd", 2);
139  descHF.add<std::string>("Subdetector", "HF");
140  descHF.add<bool>("correctForTimeslew", false);
141  descHF.add<bool>("dropZSmarkedPassed", true);
142  descHF.add<bool>("correctForPhaseContainment", false);
143  descHF.add<int>("firstSample", 4);
144  descriptions.add("hfsimplereco", descHF);
145 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ process()

template<class DIGICOLL , class RECHITCOLL >
void HcalSimpleReconstructor::process ( edm::Event e,
const edm::EventSetup c,
const edm::EDGetTokenT< DIGICOLL > &  tok 
)
private

Definition at line 66 of file HcalSimpleReconstructor.cc.

References AlignmentProducer_cff::calibrations, submitPVValidationJobs::conditions, conditionsToken_, dropZSmarkedPassed_, MillePedeFileConverter_cfg::e, options_cfi::eventSetup, dqmdumpme::first, HcalRecoParam::firstSample(), firstSample_, mps_fire::i, eostools::move(), paramTS_, DetId::rawId(), reco_, HcalSimpleRecAlgo::reconstruct(), HcalRecoParam::samplesToAdd(), samplesToAdd_, and tsFromDB_.

68  {
69  // get conditions
71 
73  e.getByToken(tok, digi);
74 
75  // create empty output
76  auto rec = std::make_unique<RECHITCOLL>();
77  rec->reserve(digi->size());
78  // run the algorithm
79  int first = firstSample_;
80  int toadd = samplesToAdd_;
81  typename DIGICOLL::const_iterator i;
82  for (i = digi->begin(); i != digi->end(); i++) {
83  HcalDetId cell = i->id();
84  DetId detcell = (DetId)cell;
85  // rof 27.03.09: drop ZS marked and passed digis:
87  if (i->zsMarkAndPass())
88  continue;
89 
90  const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
91  const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
92  const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
93  HcalCoderDb coder(*channelCoder, *shape);
94 
95  //>>> firstSample & samplesToAdd
96  if (tsFromDB_) {
97  const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
98  first = param_ts->firstSample();
99  toadd = param_ts->samplesToAdd();
100  }
101  rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
102  }
103  // return result
104  e.put(std::move(rec));
105 }
std::unique_ptr< HcalRecoParams > paramTS_
constexpr unsigned int firstSample() const
Definition: HcalRecoParam.h:31
edm::ESGetToken< HcalDbService, HcalDbRecord > conditionsToken_
Definition: DetId.h:17
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
constexpr unsigned int samplesToAdd() const
Definition: HcalRecoParam.h:32
HFRecHit reconstruct(const HFDataFrame &digi, int first, int toadd, const HcalCoder &coder, const HcalCalibrations &calibs) const
def move(src, dest)
Definition: eostools.py:511

◆ produce()

void HcalSimpleReconstructor::produce ( edm::Event e,
const edm::EventSetup c 
)
final

Definition at line 107 of file HcalSimpleReconstructor.cc.

References det_, MillePedeFileConverter_cfg::e, options_cfi::eventSetup, DetId::Hcal, HcalCalibration, HcalForward, HcalOther, HcalOuter, subdet_, subdetOther_, tok_calib_, tok_hf_, and tok_ho_.

107  {
108  if (det_ == DetId::Hcal) {
109  if (subdet_ == HcalForward) {
110  process<HFDigiCollection, HFRecHitCollection>(e, eventSetup, tok_hf_);
111  } else if (subdet_ == HcalOuter) {
112  process<HODigiCollection, HORecHitCollection>(e, eventSetup, tok_ho_);
113  } else if (subdet_ == HcalOther && subdetOther_ == HcalCalibration) {
114  process<HcalCalibDigiCollection, HcalCalibRecHitCollection>(e, eventSetup, tok_calib_);
115  }
116  }
117 }
HcalOtherSubdetector subdetOther_
edm::EDGetTokenT< HFDigiCollection > tok_hf_
edm::EDGetTokenT< HODigiCollection > tok_ho_
edm::EDGetTokenT< HcalCalibDigiCollection > tok_calib_

Member Data Documentation

◆ conditionsToken_

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

Definition at line 68 of file HcalSimpleReconstructor.h.

Referenced by HcalSimpleReconstructor(), and process().

◆ det_

DetId::Detector HcalSimpleReconstructor::det_
private

Definition at line 46 of file HcalSimpleReconstructor.h.

Referenced by produce().

◆ dropZSmarkedPassed_

bool HcalSimpleReconstructor::dropZSmarkedPassed_
private

Definition at line 55 of file HcalSimpleReconstructor.h.

Referenced by process().

◆ firstSample_

int HcalSimpleReconstructor::firstSample_
private

Definition at line 59 of file HcalSimpleReconstructor.h.

Referenced by process().

◆ htopoToken_

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

Definition at line 66 of file HcalSimpleReconstructor.h.

Referenced by beginRun(), and HcalSimpleReconstructor().

◆ inputLabel_

edm::InputTag HcalSimpleReconstructor::inputLabel_
private

Definition at line 49 of file HcalSimpleReconstructor.h.

Referenced by HcalSimpleReconstructor().

◆ paramsToken_

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

Definition at line 67 of file HcalSimpleReconstructor.h.

Referenced by beginRun(), and HcalSimpleReconstructor().

◆ paramTS_

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

Definition at line 63 of file HcalSimpleReconstructor.h.

Referenced by beginRun(), and process().

◆ reco_

HcalSimpleRecAlgo HcalSimpleReconstructor::reco_
private

Definition at line 45 of file HcalSimpleReconstructor.h.

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

◆ samplesToAdd_

int HcalSimpleReconstructor::samplesToAdd_
private

Definition at line 60 of file HcalSimpleReconstructor.h.

Referenced by process().

◆ subdet_

int HcalSimpleReconstructor::subdet_
private

Definition at line 47 of file HcalSimpleReconstructor.h.

Referenced by HcalSimpleReconstructor(), and produce().

◆ subdetOther_

HcalOtherSubdetector HcalSimpleReconstructor::subdetOther_
private

Definition at line 48 of file HcalSimpleReconstructor.h.

Referenced by produce().

◆ tok_calib_

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

Definition at line 53 of file HcalSimpleReconstructor.h.

Referenced by HcalSimpleReconstructor(), and produce().

◆ tok_hf_

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

Definition at line 51 of file HcalSimpleReconstructor.h.

Referenced by HcalSimpleReconstructor(), and produce().

◆ tok_ho_

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

Definition at line 52 of file HcalSimpleReconstructor.h.

Referenced by HcalSimpleReconstructor(), and produce().

◆ tsFromDB_

bool HcalSimpleReconstructor::tsFromDB_
private

Definition at line 61 of file HcalSimpleReconstructor.h.

Referenced by beginRun(), HcalSimpleReconstructor(), and process().