CMS 3D CMS Logo

PFRecHitProducerKernel.dev.cc
Go to the documentation of this file.
1 #include <type_traits>
2 
3 #include <alpaka/alpaka.hpp>
4 
7 
9 
11  using namespace particleFlowRecHitProducer;
12 
13  // Kernel to apply cuts to calorimeter hits and construct PFRecHits
14  template <typename CAL>
16  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
17  ALPAKA_FN_ACC void operator()(const TAcc& acc,
18  const typename CAL::ParameterType::ConstView params,
19  const typename CAL::TopologyTypeDevice::ConstView topology,
20  const typename CAL::CaloRecHitSoATypeDevice::ConstView recHits,
22  uint32_t* __restrict__ denseId2pfRecHit,
23  uint32_t* __restrict__ num_pfRecHits) const {
24  // Strided loop over CaloRecHits
25  for (int32_t i : cms::alpakatools::elements_with_stride(acc, recHits.metadata().size())) {
26  // Check energy thresholds/quality cuts (specialised for HCAL/ECAL)
27  if (!applyCuts(recHits[i], params, topology))
28  continue;
29 
30  // Use atomic operation to determine index of the PFRecHit to be constructed
31  // The index needs to be unique and consequtive across all threads in all blocks.
32  // This is achieved using the alpaka::hierarchy::Blocks argument.
33  const uint32_t j = alpaka::atomicInc(acc, num_pfRecHits, 0xffffffff, alpaka::hierarchy::Blocks{});
34 
35  // Construct PFRecHit from CAL recHit (specialised for HCAL/ECAL)
36  constructPFRecHit(pfRecHits[j], recHits[i]);
37 
38  // Fill denseId -> pfRecHit index map
39  denseId2pfRecHit[CAL::detId2denseId(pfRecHits.detId(j))] = j;
40  }
41  }
42 
43  ALPAKA_FN_ACC static bool applyCuts(const typename CAL::CaloRecHitSoATypeDevice::ConstView::const_element rh,
44  const typename CAL::ParameterType::ConstView params,
45  const typename CAL::TopologyTypeDevice::ConstView topology);
46 
47  ALPAKA_FN_ACC static void constructPFRecHit(
48  reco::PFRecHitDeviceCollection::View::element pfrh,
49  const typename CAL::CaloRecHitSoATypeDevice::ConstView::const_element rh);
50  };
51 
52  template <>
54  const typename HCAL::CaloRecHitSoATypeDevice::ConstView::const_element rh,
57  // Reject HCAL recHits below enery threshold
58  float threshold = 9999.;
59  const uint32_t detId = rh.detId();
60  const uint32_t depth = HCAL::getDepth(detId);
61  const uint32_t subdet = getSubdet(detId);
62  if (topology.cutsFromDB()) {
63  threshold = topology.noiseThreshold()[HCAL::detId2denseId(detId)];
64  } else {
65  if (subdet == HcalBarrel) {
66  threshold = params.energyThresholds()[depth - 1];
67  } else if (subdet == HcalEndcap) {
68  threshold = params.energyThresholds()[depth - 1 + HCAL::kMaxDepthHB];
69  } else {
70  printf("Rechit with detId %u has invalid subdetector %u!\n", detId, subdet);
71  return false;
72  }
73  }
74  return rh.energy() >= threshold;
75  }
76 
77  template <>
79  const ECAL::CaloRecHitSoATypeDevice::ConstView::const_element rh,
82  // Reject ECAL recHits below energy threshold
83  if (rh.energy() < params.energyThresholds()[ECAL::detId2denseId(rh.detId())])
84  return false;
85 
86  // Reject ECAL recHits of bad quality
87  if ((ECAL::checkFlag(rh.flags(), ECAL::Flags::kOutOfTime) && rh.energy() > params.cleaningThreshold()) ||
89  ECAL::checkFlag(rh.flags(), ECAL::Flags::kDiWeird))
90  return false;
91 
92  return true;
93  }
94 
95  template <>
97  reco::PFRecHitDeviceCollection::View::element pfrh,
98  const HCAL::CaloRecHitSoATypeDevice::ConstView::const_element rh) {
99  pfrh.detId() = rh.detId();
100  pfrh.denseId() = HCAL::detId2denseId(rh.detId());
101  pfrh.energy() = rh.energy();
102  pfrh.time() = rh.time();
103  pfrh.depth() = HCAL::getDepth(pfrh.detId());
104  const uint32_t subdet = getSubdet(pfrh.detId());
105  if (subdet == HcalBarrel)
106  pfrh.layer() = PFLayer::HCAL_BARREL1;
107  else if (subdet == HcalEndcap)
108  pfrh.layer() = PFLayer::HCAL_ENDCAP;
109  else
110  pfrh.layer() = PFLayer::NONE;
111  }
112 
113  template <>
115  reco::PFRecHitDeviceCollection::View::element pfrh,
116  const ECAL::CaloRecHitSoATypeDevice::ConstView::const_element rh) {
117  pfrh.detId() = rh.detId();
118  pfrh.denseId() = ECAL::detId2denseId(rh.detId());
119  pfrh.energy() = rh.energy();
120  pfrh.time() = rh.time();
121  pfrh.depth() = 1;
122  const uint32_t subdet = getSubdet(pfrh.detId());
123  if (subdet == EcalBarrel)
124  pfrh.layer() = PFLayer::ECAL_BARREL;
125  else if (subdet == EcalEndcap)
126  pfrh.layer() = PFLayer::ECAL_ENDCAP;
127  else
128  pfrh.layer() = PFLayer::NONE;
129  }
130 
131  // Kernel to associate topology information of PFRecHits
132  template <typename CAL>
134  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
135  ALPAKA_FN_ACC void operator()(const TAcc& acc,
136  const typename CAL::TopologyTypeDevice::ConstView topology,
138  const uint32_t* __restrict__ denseId2pfRecHit,
139  uint32_t* __restrict__ num_pfRecHits) const {
140  // First thread updates size field pfRecHits SoA
141  if (const int32_t thread = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc)[0u]; thread == 0)
142  pfRecHits.size() = *num_pfRecHits;
143 
144  // Assign position information and associate neighbours
145  for (int32_t i : cms::alpakatools::elements_with_stride(acc, *num_pfRecHits)) {
146  const uint32_t denseId = CAL::detId2denseId(pfRecHits.detId(i));
147 
148  pfRecHits.x(i) = topology.positionX(denseId);
149  pfRecHits.y(i) = topology.positionY(denseId);
150  pfRecHits.z(i) = topology.positionZ(denseId);
151 
152  for (uint32_t n = 0; n < 8; n++) {
153  pfRecHits.neighbours(i)(n) = -1;
154  const uint32_t denseId_neighbour = topology.neighbours(denseId)(n);
155  if (denseId_neighbour != 0xffffffff) {
156  const uint32_t pfRecHit_neighbour = denseId2pfRecHit[denseId_neighbour];
157  if (pfRecHit_neighbour != 0xffffffff)
158  pfRecHits.neighbours(i)(n) = (int32_t)pfRecHit_neighbour;
159  }
160  }
161  }
162  }
163  };
164 
165  template <typename CAL>
166  PFRecHitProducerKernel<CAL>::PFRecHitProducerKernel(Queue& queue, const uint32_t num_recHits)
167  : denseId2pfRecHit_(cms::alpakatools::make_device_buffer<uint32_t[]>(queue, CAL::kSize)),
168  num_pfRecHits_(cms::alpakatools::make_device_buffer<uint32_t>(queue)),
169  work_div_(cms::alpakatools::make_workdiv<Acc1D>(1, 1)) {
170  alpaka::memset(queue, denseId2pfRecHit_, 0xff); // Reset denseId -> pfRecHit index map
171  alpaka::memset(queue, num_pfRecHits_, 0x00); // Reset global pfRecHit counter
172 
173  const uint32_t items = 64;
174  const uint32_t groups = cms::alpakatools::divide_up_by(num_recHits, items);
175  work_div_ = cms::alpakatools::make_workdiv<Acc1D>(groups, items);
176  }
177 
178  template <typename CAL>
180  const typename CAL::CaloRecHitSoATypeDevice& recHits,
181  const typename CAL::ParameterType& params,
182  const typename CAL::TopologyTypeDevice& topology,
184  alpaka::exec<Acc1D>(queue,
185  work_div_,
187  params.view(),
188  topology.view(),
189  recHits.view(),
190  pfRecHits.view(),
191  denseId2pfRecHit_.data(),
192  num_pfRecHits_.data());
193  }
194 
195  template <typename CAL>
197  const typename CAL::TopologyTypeDevice& topology,
199  alpaka::exec<Acc1D>(queue,
200  work_div_,
202  topology.view(),
203  pfRecHits.view(),
204  denseId2pfRecHit_.data(),
205  num_pfRecHits_.data());
206  }
207 
208  // Instantiate templates
209  template class PFRecHitProducerKernel<HCAL>;
210  template class PFRecHitProducerKernel<ECAL>;
211 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
cms::alpakatools::device_buffer< Device, uint32_t > num_pfRecHits_
WorkDiv< Dim1D > make_workdiv(Idx blocks, Idx elements)
Definition: workdivision.h:46
void associateTopologyInfo(Queue &queue, const typename CAL::TopologyTypeDevice &topology, reco::PFRecHitDeviceCollection &pfRecHits)
constexpr Idx divide_up_by(Idx value, Idx divisor)
Definition: workdivision.h:19
PFRecHitProducerKernel(Queue &queue, const uint32_t num_recHits)
ALPAKA_FN_ACC void operator()(const TAcc &acc, const typename CAL::TopologyTypeDevice::ConstView topology, reco::PFRecHitDeviceCollection::View pfRecHits, const uint32_t *__restrict__ denseId2pfRecHit, uint32_t *__restrict__ num_pfRecHits) const
std::enable_if_t< alpaka::isDevice< TDev > and not std::is_array_v< T >, device_buffer< TDev, T > > make_device_buffer(TDev const &device)
Definition: memory.h:177
static ALPAKA_FN_ACC void constructPFRecHit(reco::PFRecHitDeviceCollection::View::element pfrh, const typename CAL::CaloRecHitSoATypeDevice::ConstView::const_element rh)
typename Layout::ConstView ConstView
static constexpr bool checkFlag(uint32_t flagBits, int flag)
T1 atomicInc(T1 *a, T2 b)
Definition: cudaCompat.h:48
std::vector< Block > Blocks
Definition: Block.h:99
Namespace of DDCMS conversion namespace.
static ALPAKA_FN_ACC bool applyCuts(const typename CAL::CaloRecHitSoATypeDevice::ConstView::const_element rh, const typename CAL::ParameterType::ConstView params, const typename CAL::TopologyTypeDevice::ConstView topology)
void processRecHits(Queue &queue, const typename CAL::CaloRecHitSoATypeDevice &recHits, const typename CAL::ParameterType &params, const typename CAL::TopologyTypeDevice &topology, reco::PFRecHitDeviceCollection &pfRecHits)
cms::alpakatools::device_buffer< Device, uint32_t[]> denseId2pfRecHit_
ALPAKA_FN_ACC void operator()(const TAcc &acc, const typename CAL::ParameterType::ConstView params, const typename CAL::TopologyTypeDevice::ConstView topology, const typename CAL::CaloRecHitSoATypeDevice::ConstView recHits, reco::PFRecHitDeviceCollection::View pfRecHits, uint32_t *__restrict__ denseId2pfRecHit, uint32_t *__restrict__ num_pfRecHits) const