CMS 3D CMS Logo

HcalSimpleReconstructor.cc
Go to the documentation of this file.
13 
14 #include <iostream>
15 
17  reco_(conf.getParameter<bool>("correctForTimeslew"),
18  conf.getParameter<bool>("correctForPhaseContainment"),conf.getParameter<double>("correctionPhaseNS")),
19  det_(DetId::Hcal),
20  inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
21  dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
22  firstSample_(conf.getParameter<int>("firstSample")),
23  samplesToAdd_(conf.getParameter<int>("samplesToAdd")),
24  tsFromDB_(conf.getParameter<bool>("tsFromDB")),
25  paramTS(nullptr),
26  theTopology(nullptr)
27 {
28 
29  // register for data access
30  tok_hf_ = consumes<HFDigiCollection>(inputLabel_);
31  tok_ho_ = consumes<HODigiCollection>(inputLabel_);
32  tok_calib_ = consumes<HcalCalibDigiCollection>(inputLabel_);
33 
34  std::string subd=conf.getParameter<std::string>("Subdetector");
35  if (!strcasecmp(subd.c_str(),"HO")) {
37  produces<HORecHitCollection>();
38  }
39  else if (!strcasecmp(subd.c_str(),"HF")) {
41  produces<HFRecHitCollection>();
42  }
43  else {
44  std::cout << "HcalSimpleReconstructor is not associated with a specific subdetector!" << std::endl;
45  }
46 
47 }
48 
50  delete paramTS;
51  delete theTopology;
52 }
53 
55  if(tsFromDB_) {
57  es.get<HcalRecoParamsRcd>().get(p);
58  paramTS = new HcalRecoParams(*p.product());
59 
61  es.get<HcalRecNumberingRecord>().get(htopo);
62  theTopology=new HcalTopology(*htopo);
64 
65  }
66  reco_.beginRun(es);
67 }
68 
70  if(tsFromDB_ && paramTS) {
71  delete paramTS;
72  paramTS = nullptr;
73  reco_.endRun();
74  }
75 }
76 
77 
78 template<class DIGICOLL, class RECHITCOLL>
80 {
81  // get conditions
83  eventSetup.get<HcalDbRecord>().get(conditions);
84 
86  e.getByToken(tok,digi);
87 
88  // create empty output
89  auto rec = std::make_unique<RECHITCOLL>();
90  rec->reserve(digi->size());
91  // run the algorithm
92  int first = firstSample_;
93  int toadd = samplesToAdd_;
94  typename DIGICOLL::const_iterator i;
95  for (i=digi->begin(); i!=digi->end(); i++) {
96  HcalDetId cell = i->id();
97  DetId detcell=(DetId)cell;
98  // rof 27.03.09: drop ZS marked and passed digis:
100  if (i->zsMarkAndPass()) continue;
101 
102  const HcalCalibrations& calibrations=conditions->getHcalCalibrations(cell);
103  const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
104  const HcalQIEShape* shape = conditions->getHcalShape (channelCoder);
105  HcalCoderDb coder (*channelCoder, *shape);
106 
107  //>>> firstSample & samplesToAdd
108  if(tsFromDB_) {
109  const HcalRecoParam* param_ts = paramTS->getValues(detcell.rawId());
110  first = param_ts->firstSample();
111  toadd = param_ts->samplesToAdd();
112  }
113  rec->push_back(reco_.reconstruct(*i,first,toadd,coder,calibrations));
114  }
115  // return result
116  e.put(std::move(rec));
117 }
118 
120 {
121  if (det_==DetId::Hcal) {
122  if (subdet_==HcalForward) {
123  process<HFDigiCollection, HFRecHitCollection>(e, eventSetup, tok_hf_);
124  } else if (subdet_==HcalOuter) {
125  process<HODigiCollection, HORecHitCollection>(e, eventSetup, tok_ho_);
126  } else if (subdet_==HcalOther && subdetOther_==HcalCalibration) {
127  process<HcalCalibDigiCollection, HcalCalibRecHitCollection>(e, eventSetup, tok_calib_);
128  }
129  }
130 }
131 
132 void
134 
135 
136  // horeco
138  descHO.add<double>("correctionPhaseNS", 13.0);
139  descHO.add<edm::InputTag>("digiLabel", edm::InputTag("hcalDigis"));
140  descHO.add<bool>("tsFromDB", true);
141  descHO.add<int>("samplesToAdd", 4);
142  descHO.add<std::string>("Subdetector", "HO");
143  descHO.add<bool>("correctForTimeslew", true);
144  descHO.add<bool>("dropZSmarkedPassed", true);
145  descHO.add<bool>("correctForPhaseContainment", true);
146  descHO.add<int>("firstSample", 4);
147  descriptions.add("hosimplereco", descHO);
148 
149  // hfreco
151  descHF.add<double>("correctionPhaseNS", 0.0);
152  descHF.add<edm::InputTag>("digiLabel", edm::InputTag("hcalDigis"));
153  descHF.add<bool>("tsFromDB", true);
154  descHF.add<int>("samplesToAdd", 2);
155  descHF.add<std::string>("Subdetector", "HF");
156  descHF.add<bool>("correctForTimeslew", false);
157  descHF.add<bool>("dropZSmarkedPassed", true);
158  descHF.add<bool>("correctForPhaseContainment", false);
159  descHF.add<int>("firstSample", 4);
160  descriptions.add("hfsimplereco", descHF);
161 }
T getParameter(std::string const &) const
constexpr unsigned int samplesToAdd() const
Definition: HcalRecoParam.h:33
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
HcalSimpleReconstructor(const edm::ParameterSet &ps)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void beginRun(edm::EventSetup const &es)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:47
void process(edm::Event &e, const edm::EventSetup &c, const edm::EDGetTokenT< DIGICOLL > &tok)
const Item * getValues(DetId fId, bool throwOnFail=true) const
#define nullptr
HcalOtherSubdetector subdetOther_
void produce(edm::Event &e, const edm::EventSetup &c) final
void beginRun(edm::Run const &r, edm::EventSetup const &es) final
HFRecHit reconstruct(const HFDataFrame &digi, int first, int toadd, const HcalCoder &coder, const HcalCalibrations &calibs) const
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< HFDigiCollection > tok_hf_
edm::EDGetTokenT< HODigiCollection > tok_ho_
Definition: DetId.h:18
const HcalQIECoder * getHcalCoder(const HcalGenericDetId &fId) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const HcalQIEShape * getHcalShape(const HcalGenericDetId &fId) const
HLT enums.
T get() const
Definition: EventSetup.h:62
edm::EDGetTokenT< HcalCalibDigiCollection > tok_calib_
const HcalCalibrations & getHcalCalibrations(const HcalGenericDetId &fId) const
void endRun(edm::Run const &r, edm::EventSetup const &es) final
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511
void setTopo(const HcalTopology *topo)
constexpr unsigned int firstSample() const
Definition: HcalRecoParam.h:32
Definition: Run.h:44