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::CaloRecHitSoATypeDevice::ConstView recHits,
21  uint32_t* __restrict__ denseId2pfRecHit,
22  uint32_t* __restrict__ num_pfRecHits) const {
23  // Strided loop over CaloRecHits
24  for (int32_t i : cms::alpakatools::elements_with_stride(acc, recHits.metadata().size())) {
25  // Check energy thresholds/quality cuts (specialised for HCAL/ECAL)
26  if (!applyCuts(recHits[i], params))
27  continue;
28 
29  // Use atomic operation to determine index of the PFRecHit to be constructed
30  // The index needs to be unique and consequtive across all threads in all blocks.
31  // This is achieved using the alpaka::hierarchy::Blocks argument.
32  const uint32_t j = alpaka::atomicInc(acc, num_pfRecHits, 0xffffffff, alpaka::hierarchy::Blocks{});
33 
34  // Construct PFRecHit from CAL recHit (specialised for HCAL/ECAL)
35  constructPFRecHit(pfRecHits[j], recHits[i]);
36 
37  // Fill denseId -> pfRecHit index map
38  denseId2pfRecHit[CAL::detId2denseId(pfRecHits.detId(j))] = j;
39  }
40  }
41 
42  ALPAKA_FN_ACC static bool applyCuts(const typename CAL::CaloRecHitSoATypeDevice::ConstView::const_element rh,
43  const typename CAL::ParameterType::ConstView params);
44 
45  ALPAKA_FN_ACC static void constructPFRecHit(
46  reco::PFRecHitDeviceCollection::View::element pfrh,
47  const typename CAL::CaloRecHitSoATypeDevice::ConstView::const_element rh);
48  };
49 
50  template <>
52  const typename HCAL::CaloRecHitSoATypeDevice::ConstView::const_element rh,
54  // Reject HCAL recHits below enery threshold
55  float threshold = 9999.;
56  const uint32_t detId = rh.detId();
57  const uint32_t depth = HCAL::getDepth(detId);
58  const uint32_t subdet = getSubdet(detId);
59  if (subdet == HcalBarrel) {
60  threshold = params.energyThresholds()[depth - 1];
61  } else if (subdet == HcalEndcap) {
62  threshold = params.energyThresholds()[depth - 1 + HCAL::kMaxDepthHB];
63  } else {
64  printf("Rechit with detId %u has invalid subdetector %u!\n", detId, subdet);
65  return false;
66  }
67  return rh.energy() >= threshold;
68  }
69 
70  template <>
72  const ECAL::CaloRecHitSoATypeDevice::ConstView::const_element rh, const ECAL::ParameterType::ConstView params) {
73  // Reject ECAL recHits below energy threshold
74  if (rh.energy() < params.energyThresholds()[ECAL::detId2denseId(rh.detId())])
75  return false;
76 
77  // Reject ECAL recHits of bad quality
78  if ((ECAL::checkFlag(rh.flags(), ECAL::Flags::kOutOfTime) && rh.energy() > params.cleaningThreshold()) ||
80  ECAL::checkFlag(rh.flags(), ECAL::Flags::kDiWeird))
81  return false;
82 
83  return true;
84  }
85 
86  template <>
88  reco::PFRecHitDeviceCollection::View::element pfrh,
89  const HCAL::CaloRecHitSoATypeDevice::ConstView::const_element rh) {
90  pfrh.detId() = rh.detId();
91  pfrh.energy() = rh.energy();
92  pfrh.time() = rh.time();
93  pfrh.depth() = HCAL::getDepth(pfrh.detId());
94  const uint32_t subdet = getSubdet(pfrh.detId());
95  if (subdet == HcalBarrel)
96  pfrh.layer() = PFLayer::HCAL_BARREL1;
97  else if (subdet == HcalEndcap)
98  pfrh.layer() = PFLayer::HCAL_ENDCAP;
99  else
100  pfrh.layer() = PFLayer::NONE;
101  }
102 
103  template <>
105  reco::PFRecHitDeviceCollection::View::element pfrh,
106  const ECAL::CaloRecHitSoATypeDevice::ConstView::const_element rh) {
107  pfrh.detId() = rh.detId();
108  pfrh.energy() = rh.energy();
109  pfrh.time() = rh.time();
110  pfrh.depth() = 1;
111  const uint32_t subdet = getSubdet(pfrh.detId());
112  if (subdet == EcalBarrel)
113  pfrh.layer() = PFLayer::ECAL_BARREL;
114  else if (subdet == EcalEndcap)
115  pfrh.layer() = PFLayer::ECAL_ENDCAP;
116  else
117  pfrh.layer() = PFLayer::NONE;
118  }
119 
120  // Kernel to associate topology information of PFRecHits
121  template <typename CAL>
123  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
124  ALPAKA_FN_ACC void operator()(const TAcc& acc,
125  const typename CAL::TopologyTypeDevice::ConstView topology,
127  const uint32_t* __restrict__ denseId2pfRecHit,
128  uint32_t* __restrict__ num_pfRecHits) const {
129  // First thread updates size field pfRecHits SoA
130  if (const int32_t thread = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc)[0u]; thread == 0)
131  pfRecHits.size() = *num_pfRecHits;
132 
133  // Assign position information and associate neighbours
134  for (int32_t i : cms::alpakatools::elements_with_stride(acc, *num_pfRecHits)) {
135  const uint32_t denseId = CAL::detId2denseId(pfRecHits.detId(i));
136 
137  pfRecHits.x(i) = topology.positionX(denseId);
138  pfRecHits.y(i) = topology.positionY(denseId);
139  pfRecHits.z(i) = topology.positionZ(denseId);
140 
141  for (uint32_t n = 0; n < 8; n++) {
142  pfRecHits.neighbours(i)(n) = -1;
143  const uint32_t denseId_neighbour = topology.neighbours(denseId)(n);
144  if (denseId_neighbour != 0xffffffff) {
145  const uint32_t pfRecHit_neighbour = denseId2pfRecHit[denseId_neighbour];
146  if (pfRecHit_neighbour != 0xffffffff)
147  pfRecHits.neighbours(i)(n) = (int32_t)pfRecHit_neighbour;
148  }
149  }
150  }
151  }
152  };
153 
154  template <typename CAL>
155  PFRecHitProducerKernel<CAL>::PFRecHitProducerKernel(Queue& queue, const uint32_t num_recHits)
156  : denseId2pfRecHit_(cms::alpakatools::make_device_buffer<uint32_t[]>(queue, CAL::kSize)),
157  num_pfRecHits_(cms::alpakatools::make_device_buffer<uint32_t>(queue)),
158  work_div_(cms::alpakatools::make_workdiv<Acc1D>(1, 1)) {
159  alpaka::memset(queue, denseId2pfRecHit_, 0xff); // Reset denseId -> pfRecHit index map
160  alpaka::memset(queue, num_pfRecHits_, 0x00); // Reset global pfRecHit counter
161 
162  const uint32_t items = 64;
163  const uint32_t groups = cms::alpakatools::divide_up_by(num_recHits, items);
164  work_div_ = cms::alpakatools::make_workdiv<Acc1D>(groups, items);
165  }
166 
167  template <typename CAL>
169  const typename CAL::CaloRecHitSoATypeDevice& recHits,
170  const typename CAL::ParameterType& params,
171  reco::PFRecHitDeviceCollection& pfRecHits) {
172  alpaka::exec<Acc1D>(queue,
173  work_div_,
175  params.view(),
176  recHits.view(),
177  pfRecHits.view(),
178  denseId2pfRecHit_.data(),
179  num_pfRecHits_.data());
180  }
181 
182  template <typename CAL>
184  const typename CAL::TopologyTypeDevice& topology,
185  reco::PFRecHitDeviceCollection& pfRecHits) {
186  alpaka::exec<Acc1D>(queue,
187  work_div_,
189  topology.view(),
190  pfRecHits.view(),
191  denseId2pfRecHit_.data(),
192  num_pfRecHits_.data());
193  }
194 
195  // Instantiate templates
196  template class PFRecHitProducerKernel<HCAL>;
197  template class PFRecHitProducerKernel<ECAL>;
198 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
void processRecHits(Queue &queue, const typename CAL::CaloRecHitSoATypeDevice &recHits, const typename CAL::ParameterType &params, reco::PFRecHitDeviceCollection &pfRecHits)
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
ALPAKA_FN_ACC void operator()(const TAcc &acc, const typename CAL::ParameterType::ConstView params, const typename CAL::CaloRecHitSoATypeDevice::ConstView recHits, reco::PFRecHitDeviceCollection::View pfRecHits, uint32_t *__restrict__ denseId2pfRecHit, uint32_t *__restrict__ num_pfRecHits) const
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
static ALPAKA_FN_ACC bool applyCuts(const typename CAL::CaloRecHitSoATypeDevice::ConstView::const_element rh, const typename CAL::ParameterType::ConstView params)
Namespace of DDCMS conversion namespace.
cms::alpakatools::device_buffer< Device, uint32_t[]> denseId2pfRecHit_