CMS 3D CMS Logo

HcalDigisSoAProducer.cc
Go to the documentation of this file.
5 
11 
17 
19 
21  public:
22  explicit HcalDigisSoAProducer(edm::ParameterSet const& ps);
23  ~HcalDigisSoAProducer() override = default;
25 
26  private:
27  void produce(device::Event&, device::EventSetup const&) override;
28 
29  private:
30  // input product tokens
33 
34  // type aliases
37 
40 
41  // output product tokens
45 
48  };
50  };
51 
54 
55  desc.add<edm::InputTag>("hbheDigisLabel", edm::InputTag("hcalDigis"));
56  desc.add<edm::InputTag>("qie11DigiLabel", edm::InputTag("hcalDigis"));
57  desc.add<std::string>("digisLabelF01HE", std::string{"f01HEDigis"});
58  desc.add<std::string>("digisLabelF5HB", std::string{"f5HBDigis"});
59  desc.add<std::string>("digisLabelF3HB", std::string{"f3HBDigis"});
60  desc.add<uint32_t>("maxChannelsF01HE", 10000u);
61  desc.add<uint32_t>("maxChannelsF5HB", 10000u);
62  desc.add<uint32_t>("maxChannelsF3HB", 10000u);
63 
64  confDesc.addWithDefaultLabel(desc);
65  }
66 
68  : hbheDigiToken_{consumes(ps.getParameter<edm::InputTag>("hbheDigisLabel"))},
69  qie11DigiToken_{consumes(ps.getParameter<edm::InputTag>("qie11DigiLabel"))},
70  digisF01HEToken_{produces(ps.getParameter<std::string>("digisLabelF01HE"))},
71  digisF5HBToken_{produces(ps.getParameter<std::string>("digisLabelF5HB"))},
72  digisF3HBToken_{produces(ps.getParameter<std::string>("digisLabelF3HB"))} {
73  config_.maxChannelsF01HE = ps.getParameter<uint32_t>("maxChannelsF01HE");
74  config_.maxChannelsF5HB = ps.getParameter<uint32_t>("maxChannelsF5HB");
75  config_.maxChannelsF3HB = ps.getParameter<uint32_t>("maxChannelsF3HB");
76  }
77 
79  const auto& hbheDigis = event.get(hbheDigiToken_);
80  const auto& qie11Digis = event.get(qie11DigiToken_);
81 
82  //Get the number of samples in data from the first digi
83  const int stride = HBHEDataFrame::MAXSAMPLES / 2 + 1;
84  const int size = hbheDigis.size() * stride; // number of channels * stride
85 
86  // stack host memory in the queue
87  HostCollectionPhase0 hf5_(size, event.queue());
88 
89  // device product
90  DeviceCollectionPhase0 df5_(size, event.queue());
91 
92  // set SoA_Scalar;
93  hf5_.view().stride() = stride;
94  hf5_.view().size() = hbheDigis.size();
95 
96  for (unsigned int i = 0; i < hbheDigis.size(); ++i) {
97  auto const hbhe = hbheDigis[i];
98  auto const id = hbhe.id().rawId();
99  auto const presamples = hbhe.presamples();
100  uint16_t header_word = (1 << 15) | (0x5 << 12) | (0 << 10) | ((hbhe.sample(0).capid() & 0x3) << 8);
101 
102  auto hf5_vi = hf5_.view()[i];
103  hf5_vi.ids() = id;
104  hf5_vi.npresamples() = presamples;
105  hf5_vi.data()[0] = header_word;
106  //TODO:: get HEADER_WORDS/WORDS_PER_SAMPLE from DataFormat
107  for (unsigned int i = 0; i < hf5_.view().stride() - 1; i++) {
108  uint16_t s0 = (0 << 7) | (static_cast<uint8_t>(hbhe.sample(2 * i).adc()) & 0x7f);
109  uint16_t s1 = (0 << 7) | (static_cast<uint8_t>(hbhe.sample(2 * i + 1).adc()) & 0x7f);
110  uint16_t sample = (s1 << 8) | s0;
111  hf5_vi.data()[i + 1] = sample;
112  }
113  }
114  // copy hf5 to df5
115  alpaka::memcpy(event.queue(), df5_.buffer(), hf5_.const_buffer());
116  event.emplace(digisF5HBToken_, std::move(df5_));
117 
118  if (qie11Digis.empty()) {
119  event.emplace(digisF01HEToken_);
120  event.emplace(digisF3HBToken_);
121 
122  } else {
123  auto size_f1 = 0;
124  auto size_f3 = 0;
125 
126  // count the size of the SOA;
127  for (unsigned int i = 0; i < qie11Digis.size(); i++) {
128  auto const digi = QIE11DataFrame{qie11Digis[i]};
129 
130  if (digi.flavor() == 0 or digi.flavor() == 1) {
131  if (digi.detid().subdetId() == HcalEndcap) {
132  size_f1++;
133  }
134  } else if (digi.flavor() == 3) {
135  if (digi.detid().subdetId() == HcalBarrel) {
136  size_f3++;
137  }
138  }
139  }
140  auto const nsamples = qie11Digis.samples();
141  auto const stride01 = nsamples * QIE11DataFrame::WORDS_PER_SAMPLE + QIE11DataFrame::HEADER_WORDS;
142 
143  // stack host memory in the queue
144  HostCollectionPhase1 hf1_(size_f1, event.queue());
145  HostCollectionPhase1 hf3_(size_f3, event.queue());
146 
147  // device product
148  DeviceCollectionPhase1 df1_(size_f1, event.queue());
149  DeviceCollectionPhase1 df3_(size_f3, event.queue());
150 
151  // set SoA_Scalar;
152  hf1_.view().stride() = stride01;
153  hf3_.view().stride() = stride01;
154 
155  unsigned int i_f1 = 0; //counters for f1 digis
156  unsigned int i_f3 = 0; //counters for f3 digis
157 
158  for (unsigned int i = 0; i < qie11Digis.size(); i++) {
159  auto const digi = QIE11DataFrame{qie11Digis[i]};
160  assert(digi.samples() == qie11Digis.samples() && "collection nsamples must equal per digi samples");
161 
162  if (digi.flavor() == 0 or digi.flavor() == 1) {
163  if (digi.detid().subdetId() != HcalEndcap)
164  continue;
165  auto hf01_vi = hf1_.view()[i_f1];
166 
167  hf01_vi.ids() = digi.detid().rawId();
168  for (int hw = 0; hw < QIE11DataFrame::HEADER_WORDS + digi.samples(); hw++) {
169  hf01_vi.data()[hw] = (qie11Digis[i][hw]);
170  }
171  i_f1++;
172  } else if (digi.flavor() == 3) {
173  if (digi.detid().subdetId() != HcalBarrel)
174  continue;
175  auto hf03_vi = hf3_.view()[i_f3];
176 
177  hf03_vi.ids() = digi.detid().rawId();
178 
179  for (int hw = 0; hw < QIE11DataFrame::HEADER_WORDS + digi.samples(); hw++) {
180  hf03_vi.data()[hw] = (qie11Digis[i][hw]);
181  }
182  i_f3++;
183  }
184  }
185 
186  hf1_.view().size() = size_f1;
187  hf3_.view().size() = size_f3;
188 
189  alpaka::memcpy(event.queue(), df1_.buffer(), hf1_.const_buffer());
190  alpaka::memcpy(event.queue(), df3_.buffer(), hf3_.const_buffer());
191 
192  event.emplace(digisF01HEToken_, std::move(df1_));
193  event.emplace(digisF3HBToken_, std::move(df3_));
194  }
195  }
196 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
198 DEFINE_FWK_ALPAKA_MODULE(HcalDigisSoAProducer);
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
static const int MAXSAMPLES
Definition: HBHEDataFrame.h:95
device::EDPutToken< DeviceCollectionPhase0 > digisF5HBToken_
PortableCollection< HcalPhase0DigiSoA > Phase0DigiDeviceCollection
void produce(device::Event &, device::EventSetup const &) override
static const int HEADER_WORDS
PortableHostCollection< HcalPhase0DigiSoA > Phase0DigiHostCollection
assert(be >=bs)
static const int WORDS_PER_SAMPLE
constexpr uint32_t stride
Definition: HelixFit.h:22
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
PortableHostCollection< HcalPhase1DigiSoA > Phase1DigiHostCollection
hcal::Phase0DigiDeviceCollection DeviceCollectionPhase0
static void fillDescriptions(edm::ConfigurationDescriptions &)
constexpr int samples() const
total number of samples in the digi
PortableCollection< HcalPhase1DigiSoA > Phase1DigiDeviceCollection
#define DEFINE_FWK_ALPAKA_MODULE(name)
Definition: MakerMacros.h:16
device::EDPutToken< DeviceCollectionPhase1 > digisF3HBToken_
edm::EDGetTokenT< QIE11DigiCollection > qie11DigiToken_
device::EDPutToken< DeviceCollectionPhase1 > digisF01HEToken_
edm::EDGetTokenT< HBHEDigiCollection > hbheDigiToken_
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1
hcal::Phase1DigiDeviceCollection DeviceCollectionPhase1