CMS 3D CMS Logo

ZdcHitReconstructor_Run3.cc
Go to the documentation of this file.
14 
15 #include <iostream>
16 
17 #include <Eigen/Dense>
18 
19 #include <vector>
20 namespace zdchelper {
22  for (int i = 0; i < digi.samples(); i++) {
23  if (digi[i].adc() >= maxValue) {
25  break;
26  }
27  }
28  }
29 
30 } // namespace zdchelper
31 
33 
34  : reco_(conf.getParameter<int>("recoMethod")),
35  saturationFlagSetter_(nullptr),
36  det_(DetId::Hcal),
37  correctionMethodEM_(conf.getParameter<int>("correctionMethodEM")),
38  correctionMethodHAD_(conf.getParameter<int>("correctionMethodHAD")),
39  correctionMethodRPD_(conf.getParameter<int>("correctionMethodRPD")),
40  ootpuRatioEM_(conf.getParameter<double>("ootpuRatioEM")),
41  ootpuRatioHAD_(conf.getParameter<double>("ootpuRatioHAD")),
42  ootpuRatioRPD_(conf.getParameter<double>("ootpuRatioRPD")),
43  ootpuFracEM_(conf.getParameter<double>("ootpuFracEM")),
44  ootpuFracHAD_(conf.getParameter<double>("ootpuFracHAD")),
45  ootpuFracRPD_(conf.getParameter<double>("ootpuFracRPD")),
46  chargeRatiosEM_(conf.getParameter<std::vector<double>>("chargeRatiosEM")),
47  chargeRatiosHAD_(conf.getParameter<std::vector<double>>("chargeRatiosHAD")),
48  chargeRatiosRPD_(conf.getParameter<std::vector<double>>("chargeRatiosRPD")),
49  bxTs_(conf.getParameter<std::vector<unsigned int>>("bxTs")),
50  nTs_(conf.getParameter<int>("nTs")),
51  forceSOI_(conf.getParameter<bool>("forceSOI")),
52  signalSOI_(conf.getParameter<std::vector<unsigned int>>("signalSOI")),
53  noiseSOI_(conf.getParameter<std::vector<unsigned int>>("noiseSOI")),
54  setSaturationFlags_(conf.getParameter<bool>("setSaturationFlags")),
55  dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
56  skipRPD_(conf.getParameter<bool>("skipRPD")) {
57  tok_input_QIE10 = consumes<QIE10DigiCollection>(conf.getParameter<edm::InputTag>("digiLabelQIE10ZDC"));
58 
59  std::string subd = conf.getParameter<std::string>("Subdetector");
60 
61  if (setSaturationFlags_) {
62  const edm::ParameterSet& pssat = conf.getParameter<edm::ParameterSet>("saturationParameters");
63  maxADCvalue_ = pssat.getParameter<int>("maxADCvalue");
64  }
65  if (!strcasecmp(subd.c_str(), "ZDC")) {
66  det_ = DetId::Calo;
68  produces<ZDCRecHitCollection>();
69  } else if (!strcasecmp(subd.c_str(), "CALIB")) {
72  produces<HcalCalibRecHitCollection>();
73  } else {
74  std::cout << "ZdcHitReconstructor_Run3 is not associated with a specific subdetector!" << std::endl;
75  }
85  // ES tokens
86  htopoToken_ = esConsumes<HcalTopology, HcalRecNumberingRecord, edm::Transition::BeginRun>();
87  paramsToken_ = esConsumes<HcalLongRecoParams, HcalLongRecoParamsRcd, edm::Transition::BeginRun>();
88  conditionsToken_ = esConsumes<HcalDbService, HcalDbRecord>();
89  qualToken_ = esConsumes<HcalChannelQuality, HcalChannelQualityRcd>(edm::ESInputTag("", "withTopo"));
90  sevToken_ = esConsumes<HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd>();
91 }
92 
94 
96  const HcalTopology& htopo = es.getData(htopoToken_);
98  longRecoParams_ = std::make_unique<HcalLongRecoParams>(p);
99  longRecoParams_->setTopo(&htopo);
100 }
101 
103 
105  // get conditions
107  const HcalChannelQuality* myqual = &eventSetup.getData(qualToken_);
108  const HcalSeverityLevelComputer* mySeverity = &eventSetup.getData(sevToken_);
109 
110  // define vectors to pass noiseTS and signalTS
111  std::vector<unsigned int> mySignalTS;
112  std::vector<unsigned int> myNoiseTS;
113 
116  e.getByToken(tok_input_QIE10, digi);
117 
118  // create empty output
119  auto rec = std::make_unique<ZDCRecHitCollection>();
120  rec->reserve(digi->size());
121 
122  // testing QEI10 conditions
123  for (auto it = digi->begin(); it != digi->end(); it++) {
124  QIE10DataFrame QIE10_i = static_cast<QIE10DataFrame>(*it);
125  HcalZDCDetId cell = QIE10_i.id();
126  bool isRPD = cell.section() == 4;
127  if (isRPD && skipRPD_)
128  continue;
129  if (cell.section() == 1 && cell.channel() > 5)
130  continue; // ignore extra EM channels
131 
132  DetId detcell = (DetId)cell;
133 
134  // check on cells to be ignored and dropped: (rof,20.Feb.09)
135  const HcalChannelStatus* mydigistatus = myqual->getValues(detcell.rawId());
136  if (mySeverity->dropChannel(mydigistatus->getValue()))
137  continue;
139  if (QIE10_i.zsMarkAndPass())
140  continue;
141 
142  const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
143  const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
144  const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
145  HcalCoderDb coder(*channelCoder, *shape);
146 
147  // pass the effective pedestals to rec hit since both ped value and width used in subtraction of pedestals
148  const HcalPedestal* effPeds = conditions->getEffectivePedestal(cell);
149 
150  if (forceSOI_)
151  rec->push_back(reco_.reconstruct(QIE10_i, noiseSOI_, signalSOI_, coder, calibrations, *effPeds));
152 
153  else {
154  const HcalLongRecoParam* myParams = longRecoParams_->getValues(detcell);
155  mySignalTS.clear();
156  myNoiseTS.clear();
157  mySignalTS = myParams->signalTS();
158  myNoiseTS = myParams->noiseTS();
159 
160  rec->push_back(reco_.reconstruct(QIE10_i, myNoiseTS, mySignalTS, coder, calibrations, *effPeds));
161  }
162  // saturationFlagSetter_ doesn't work with QIE10
163  // created new function zdchelper::setZDCSaturation to work with QIE10
164  (rec->back()).setFlags(0);
166  zdchelper::setZDCSaturation(rec->back(), QIE10_i, maxADCvalue_);
167  }
168  // return result
169  e.put(std::move(rec));
170  } // else if (det_==DetId::Calo...)
171 
172 } // void HcalHitReconstructor::produce(...)
173 
175  // zdcreco
177  desc.add<edm::InputTag>("digiLabelQIE10ZDC", edm::InputTag("hcalDigis", "ZDC"));
178  desc.add<std::string>("Subdetector", "ZDC");
179  desc.add<bool>("dropZSmarkedPassed", true);
180  desc.add<bool>("skipRPD", true);
181  desc.add<int>("recoMethod", 1);
182  desc.add<int>("correctionMethodEM", 1);
183  desc.add<int>("correctionMethodHAD", 1);
184  desc.add<int>("correctionMethodRPD", 0);
185  desc.add<double>("ootpuRatioEM", 3.0);
186  desc.add<double>("ootpuRatioHAD", 3.0);
187  desc.add<double>("ootpuRatioRPD", -1.0);
188  desc.add<double>("ootpuFracEM", 1.0);
189  desc.add<double>("ootpuFracHAD", 1.0);
190  desc.add<double>("ootpuFracRPD", 0.0);
191  desc.add<std::vector<double>>("chargeRatiosEM",
192  {
193  1.0,
194  0.23157,
195  0.10477,
196  0.06312,
197  });
198  desc.add<std::vector<double>>("chargeRatiosHAD",
199  {
200  1.0,
201  0.23157,
202  0.10477,
203  0.06312,
204  });
205  desc.add<std::vector<double>>("chargeRatiosRPD",
206  {
207  1.0,
208  0.23157,
209  0.10477,
210  0.06312,
211  });
212  desc.add<std::vector<unsigned int>>("bxTs",
213  {
214  0,
215  2,
216  4,
217  });
218  desc.add<int>("nTs", 6);
219  desc.add<bool>("forceSOI", false);
220  desc.add<std::vector<unsigned int>>("signalSOI",
221  {
222  2,
223  });
224  desc.add<std::vector<unsigned int>>("noiseSOI",
225  {
226  1,
227  });
228  desc.add<bool>("setSaturationFlags", true);
229  {
231  psd0.add<int>("maxADCvalue", 255);
232  desc.add<edm::ParameterSetDescription>("saturationParameters", psd0);
233  }
234  descriptions.add("zdcrecoRun3", desc);
235  // or use the following to generate the label from the module's C++ type
236  //descriptions.addWithDefaultLabel(desc);
237 }
238 
239 //define this as a plug-in
edm::ESGetToken< HcalLongRecoParams, HcalLongRecoParamsRcd > paramsToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::vector< unsigned int > noiseSOI_
constexpr edm::DataFrame::id_type id() const
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
ZDCRecHit reconstruct(const QIE10DataFrame &digi, const std::vector< unsigned int > &myNoiseTS, const std::vector< unsigned int > &mySignalTS, const HcalCoder &coder, const HcalCalibrations &calibs, const HcalPedestal &effPeds) const
constexpr bool zsMarkAndPass() const
edm::ESGetToken< HcalChannelQuality, HcalChannelQualityRcd > qualToken_
std::vector< double > chargeRatiosEM_
edm::ESGetToken< HcalDbService, HcalDbRecord > conditionsToken_
std::vector< unsigned int > signalSOI_
void produce(edm::Event &e, const edm::EventSetup &c) final
constexpr void setFlagField(uint32_t value, int base, int width=1)
Definition: CaloRecHit.h:36
const Item * getValues(DetId fId, bool throwOnFail=true) const
std::vector< unsigned int > bxTs_
std::unique_ptr< HcalLongRecoParams > longRecoParams_
void beginRun(edm::Run const &r, edm::EventSetup const &es) final
void initRatioSubtraction(const float ratio, const float frac, const int ZdcSection)
void endRun(edm::Run const &r, edm::EventSetup const &es) final
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< double > chargeRatiosHAD_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
uint32_t getValue() const
edm::ESGetToken< HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd > sevToken_
const_iterator end() const
std::vector< double > chargeRatiosRPD_
HcalOtherSubdetector subdetOther_
Definition: DetId.h:17
void setZDCSaturation(ZDCRecHit &rh, QIE10DataFrame &digi, int maxValue)
const_iterator begin() const
The iterator returned can not safely be used across threads.
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > htopoToken_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< unsigned int > signalTS() const
constexpr Section section() const
get the section
Definition: HcalZDCDetId.h:92
ZdcHitReconstructor_Run3(const edm::ParameterSet &ps)
void initTemplateFit(const std::vector< unsigned int > &bxTs, const std::vector< double > &chargeRatios, const int nTs, const int ZdcSection)
static constexpr int32_t SubdetectorId
Definition: HcalZDCDetId.h:35
constexpr int32_t channel() const
get the channel
Definition: HcalZDCDetId.h:112
HcalADCSaturationFlag * saturationFlagSetter_
std::vector< unsigned int > noiseTS() const
void initCorrectionMethod(const int method, const int ZdcSection)
def move(src, dest)
Definition: eostools.py:511
constexpr int samples() const
total number of samples in the digi
edm::EDGetTokenT< QIE10DigiCollection > tok_input_QIE10
Definition: Run.h:45
uint16_t *__restrict__ uint16_t const *__restrict__ adc
bool dropChannel(const uint32_t &mystatus) const