CMS 3D CMS Logo

HFPreReconstructor.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoLocalCalo/HcalRecProducers
4 // Class: HFPreReconstructor
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Thu, 05 May 2016 00:17:51 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 #include <utility>
22 #include <cassert>
23 
24 // user include files
30 
33 
36 
38 
40 
43 
47 
48 //
49 // class declaration
50 //
51 
53 public:
54  explicit HFPreReconstructor(const edm::ParameterSet&);
55  ~HFPreReconstructor() override;
56 
57  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
58 
59 private:
60  typedef std::pair<HcalDetId, int> PmtAnodeId;
61  typedef std::pair<PmtAnodeId, const HFQIE10Info*> QIE10InfoWithId;
62 
63  void produce(edm::Event&, const edm::EventSetup&) override;
64 
65  // Module configuration parameters
67  int forceSOI_;
68  int soiShift_;
70  bool tsFromDB_;
71 
72  // Other members
75  std::vector<HFQIE10Info> qie10Infos_;
76  std::vector<QIE10InfoWithId> sortedQIE10Infos_;
77 
78  // Fill qie10Infos_ from the event data
79  void fillInfos(const edm::Event& e, const edm::EventSetup& eventSetup);
80 
81  // Fill out sortedQIE10Infos_ from qie10Infos_ and return the PMT count
82  unsigned sortDataByPmt();
83 
84  // ES tokens
87 };
88 
89 //
90 // constructors and destructor
91 //
93  : inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
94  forceSOI_(conf.getParameter<int>("forceSOI")),
95  soiShift_(conf.getParameter<int>("soiShift")),
96  dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
97  tsFromDB_(conf.getParameter<bool>("tsFromDB")),
98  reco_(conf.getParameter<bool>("sumAllTimeSlices")) {
99  // Describe consumed data
100  tok_hfQIE10_ = consumes<QIE10DigiCollection>(inputLabel_);
101 
102  // Register the product
103  produces<HFPreRecHitCollection>();
104 
105  // ES tokens
106  htopoToken_ = esConsumes<HcalTopology, HcalRecNumberingRecord>();
107  propertiesToken_ = esConsumes<HcalChannelPropertiesVec, HcalChannelPropertiesRecord>();
108 }
109 
111  // do anything here that needs to be done at destruction time
112  // (e.g. close files, deallocate resources etc.)
113 }
114 
115 //
116 // member functions
117 //
119  sortedQIE10Infos_.clear();
120  unsigned pmtCount = 0;
121  const unsigned sz = qie10Infos_.size();
122  if (sz) {
123  // Perform sorting
124  sortedQIE10Infos_.reserve(sz);
125  const HFQIE10Info* info = &qie10Infos_[0];
126  for (unsigned i = 0; i < sz; ++i) {
127  const HcalDetId id(info[i].id());
128  sortedQIE10Infos_.push_back(QIE10InfoWithId(PmtAnodeId(id.baseDetId(), id.depth()), info + i));
129  }
131 
132  // Count the PMTs
133  HcalDetId previousBaseId(sortedQIE10Infos_[0].first.first);
134  pmtCount = 1;
135  for (unsigned i = 1; i < sz; ++i) {
136  const HcalDetId baseId(sortedQIE10Infos_[i].first.first);
137  if (baseId != previousBaseId) {
138  previousBaseId = baseId;
139  ++pmtCount;
140  }
141  }
142  }
143  return pmtCount;
144 }
145 
147  using namespace edm;
148 
149  // Clear the collection we want to fill in this method
150  qie10Infos_.clear();
151 
152  // Get the calibrations
153  const HcalTopology& htopo = eventSetup.getData(htopoToken_);
155 
156  // Get the input collection
158  e.getByToken(tok_hfQIE10_, digi);
159 
160  const unsigned inputSize = digi->size();
161  if (inputSize) {
162  // Process the digis and fill out the HFQIE10Info vector
163  qie10Infos_.reserve(inputSize);
164 
165  for (QIE10DigiCollection::const_iterator it = digi->begin(); it != digi->end(); ++it) {
166  const QIE10DataFrame& frame(*it);
167  const HcalDetId cell(frame.id());
168 
169  // Protection against calibration channels which are not
170  // in the database but can still come in the QIE10DataFrame
171  // in the laser calibs, etc.
172  if (cell.subdet() != HcalSubdetector::HcalForward)
173  continue;
174 
175  // Check zero suppression
177  if (frame.zsMarkAndPass())
178  continue;
179 
180  // Look up the channel properties. This lookup is O(1).
181  const HcalChannelProperties& properties(prop.at(htopo.detId2denseId(cell)));
182 
183  // ADC decoding tool
184  const HcalCoderDb coder(*properties.channelCoder, *properties.shape);
185 
186  int tsToUse = forceSOI_;
187  if (tsToUse < 0) {
188  if (tsFromDB_) {
189  tsToUse = properties.paramTs->firstSample();
190  } else {
191  // Get the "sample of interest" from the data frame itself
192  tsToUse = frame.presamples();
193  }
194  }
195 
196  // Reconstruct the charge, energy, etc
197  const HFQIE10Info& info = reco_.reconstruct(frame, tsToUse + soiShift_, coder, properties);
198  if (info.id().rawId())
199  qie10Infos_.push_back(info);
200  }
201  }
202 }
203 
204 // ------------ method called to produce the data ------------
206  // Process the input data
208 
209  // Create a new output collection
210  std::unique_ptr<HFPreRecHitCollection> out(std::make_unique<HFPreRecHitCollection>());
211 
212  // Fill the output collection
213  const unsigned pmtCount = sortDataByPmt();
214  if (pmtCount) {
215  out->reserve(pmtCount);
216  const unsigned sz = sortedQIE10Infos_.size();
217  HcalDetId previousBaseId(sortedQIE10Infos_[0].first.first);
218  unsigned nFound = 1;
219  for (unsigned i = 1; i <= sz; ++i) {
220  bool appendData = i == sz;
221  if (i < sz) {
222  const HcalDetId baseId(sortedQIE10Infos_[i].first.first);
223  if (baseId == previousBaseId)
224  ++nFound;
225  else {
226  appendData = true;
227  previousBaseId = baseId;
228  }
229  }
230 
231  if (appendData) {
232  // If we have found more than two QIE10 with the same base id,
233  // there is a bug somewhere in the dataframe. We can't do
234  // anything useful about it here. Once we make sure that
235  // everything works as expected, this assertion can be removed.
236  assert(nFound <= 2);
237 
238  const HFQIE10Info* first = nullptr;
239  const HFQIE10Info* second = sortedQIE10Infos_[i - 1].second;
240 
241  if (nFound >= 2)
242  first = sortedQIE10Infos_[i - 2].second;
243  else if (sortedQIE10Infos_[i - 1].first.second < 3) {
244  // Only one QIE10 readout found for this PMT.
245  // Arrange for depth 1 and 2 to be "first".
246  first = second;
247  second = nullptr;
248  }
249 
250  out->push_back(HFPreRecHit(sortedQIE10Infos_[i - nFound].first.first, first, second));
251 
252  // Reset the QIE find count for this base id
253  nFound = 1;
254  }
255  }
256 
257  assert(out->size() == pmtCount);
258  }
259 
260  // Add the output collection to the event record
261  e.put(std::move(out));
262 }
263 
264 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
267 
268  desc.add<edm::InputTag>("digiLabel");
269  desc.add<int>("forceSOI", -1);
270  desc.add<int>("soiShift", 0);
271  desc.add<bool>("dropZSmarkedPassed");
272  desc.add<bool>("tsFromDB");
273  desc.add<bool>("sumAllTimeSlices");
274 
275  descriptions.addDefault(desc);
276 }
277 
278 //define this as a plug-in
void produce(edm::Event &, const edm::EventSetup &) override
std::vector< HFQIE10Info > qie10Infos_
static const TGPicture * info(bool iBackgroundIsBlack)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const HcalQIEShape * shape
const HcalRecoParam * paramTs
constexpr unsigned int firstSample() const
Definition: HcalRecoParam.h:31
unsigned int detId2denseId(const DetId &id) const override
return a linear packed id
assert(be >=bs)
std::vector< QIE10InfoWithId > sortedQIE10Infos_
U second(std::pair< T, U > const &p)
edm::ESGetToken< HcalChannelPropertiesVec, HcalChannelPropertiesRecord > propertiesToken_
void addDefault(ParameterSetDescription const &psetDescription)
std::vector< HcalChannelProperties > HcalChannelPropertiesVec
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::pair< PmtAnodeId, const HFQIE10Info * > QIE10InfoWithId
void fillInfos(const edm::Event &e, const edm::EventSetup &eventSetup)
const_iterator end() const
const_iterator begin() const
The iterator returned can not safely be used across threads.
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
HLT enums.
edm::InputTag inputLabel_
edm::EDGetTokenT< QIE10DigiCollection > tok_hfQIE10_
HFQIE10Info reconstruct(const QIE10DataFrame &digi, int tsToUse, const HcalCoder &coder, const HcalChannelProperties &prop) const
Definition: HFPreRecAlgo.cc:12
def move(src, dest)
Definition: eostools.py:511
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > htopoToken_
HFPreReconstructor(const edm::ParameterSet &)
std::pair< HcalDetId, int > PmtAnodeId
const HcalQIECoder * channelCoder