22 using namespace particleFlowRecHitProducer;
24 template <
typename CAL>
28 auto cc = setWhatProduced(
this);
29 geomToken_ =
cc.consumes();
30 if constexpr (std::is_same_v<CAL, HCAL>)
31 hcalToken_ =
cc.consumes();
39 std::unique_ptr<typename CAL::TopologyTypeHost>
produce(
const typename CAL::TopologyRecordType& iRecord) {
40 const auto&
geom = iRecord.get(geomToken_);
42 auto view = product->view();
44 const int calEnums[2] = {CAL::kSubdetectorBarrelId, CAL::kSubdetectorEndcapId};
45 for (
const auto subdet : calEnums) {
51 std::variant<EcalBarrelTopology, EcalEndcapTopology> topoVar;
52 if constexpr (std::is_same_v<CAL, HCAL>)
53 topo = &iRecord.get(hcalToken_);
60 for (
auto const detId :
geom.getValidDetIds(CAL::kDetectorId, subdet)) {
61 const uint32_t denseId = CAL::detId2denseId(
detId);
62 assert(denseId < CAL::kSize);
65 view.positionX(denseId) =
pos.x();
66 view.positionY(denseId) =
pos.y();
67 view.positionZ(denseId) =
pos.z();
69 for (uint32_t
n = 0;
n < 8;
n++) {
70 uint32_t neighDetId = getNeighbourDetId(
detId,
n, *topo);
71 if (CAL::detIdInRange(neighDetId))
72 view.neighbours(denseId)(
n) = CAL::detId2denseId(neighDetId);
74 view.neighbours(denseId)(
n) = 0xffffffff;
80 if constexpr (std::is_same_v<CAL, HCAL>) {
81 for (
const auto subdet : calEnums)
82 for (
auto const detId :
geom.getValidDetIds(CAL::kDetectorId, subdet)) {
83 const uint32_t denseId = CAL::detId2denseId(
detId);
84 for (uint32_t
n = 0;
n < 8;
n++) {
86 view.neighbours(
view.neighbours(denseId)[
n]);
87 if (
std::find(neighboursOfNeighbour.begin(), neighboursOfNeighbour.end(), denseId) ==
88 neighboursOfNeighbour.end())
89 view.neighbours(denseId)[
n] = 0xffffffff;
95 LogDebug(
"PFRecHitTopologyESProducer").log([&](
auto&
log) {
96 for (
const auto subdet : calEnums)
97 for (
const auto detId :
geom.getValidDetIds(CAL::kDetectorId, subdet)) {
98 const uint32_t denseId = CAL::detId2denseId(
detId);
99 log.format(
"detId:{} denseId:{} pos:{},{},{} neighbours:{},{},{},{};{},{},{},{}\n",
102 view[denseId].positionX(),
103 view[denseId].positionY(),
104 view[denseId].positionZ(),
105 view[denseId].neighbours()(0),
106 view[denseId].neighbours()(1),
107 view[denseId].neighbours()(2),
108 view[denseId].neighbours()(3),
109 view[denseId].neighbours()(4),
110 view[denseId].neighbours()(5),
111 view[denseId].neighbours()(6),
112 view[denseId].neighbours()(7));
129 const uint32_t direction,
144 if (direction == 4) {
145 const uint32_t NE = getNeighbourDetId(getNeighbourDetId(
detId, 0, topo), 2, topo);
148 return getNeighbourDetId(getNeighbourDetId(
detId, 2, topo), 0, topo);
150 if (direction == 5) {
151 const uint32_t SW = getNeighbourDetId(getNeighbourDetId(
detId, 1, topo), 3, topo);
154 return getNeighbourDetId(getNeighbourDetId(
detId, 3, topo), 1, topo);
156 if (direction == 6) {
157 const uint32_t ES = getNeighbourDetId(getNeighbourDetId(
detId, 2, topo), 1, topo);
160 return getNeighbourDetId(getNeighbourDetId(
detId, 1, topo), 2, topo);
162 if (direction == 7) {
163 const uint32_t WN = getNeighbourDetId(getNeighbourDetId(
detId, 3, topo), 0, topo);
166 return getNeighbourDetId(getNeighbourDetId(
detId, 0, topo), 3, topo);
173 const uint32_t direction,
189 if (direction == 4) {
194 }
else if (direction == 5) {
199 }
else if (direction == 6) {
204 }
else if (direction == 7) {
212 const uint32_t nn2 = getNeighbourDetId(
detId,
directions.second, topo);
213 const uint32_t nnn = getNeighbourDetId(nn1,
directions.second, topo);
214 if (nnn == nn1 || nnn == nn2)
virtual DetId goNorth(const DetId &id) const
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
uint32_t cc[maxCellsPerHit]
PFRecHitTopologyESProducer< ECAL > PFRecHitECALTopologyESProducer
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
virtual DetId goSouth(const DetId &id) const
std::unique_ptr< typename CAL::TopologyTypeHost > produce(const typename CAL::TopologyRecordType &iRecord)
Eigen::Matrix< uint32_t, 8, 1 > PFRecHitsTopologyNeighbours
static constexpr int getZside(uint32_t detId)
PFRecHitTopologyESProducer(edm::ParameterSet const &iConfig)
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > hcalToken_
virtual DetId goEast(const DetId &id) const
PFRecHitTopologyESProducer< HCAL > PFRecHitHCALTopologyESProducer
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
uint32_t getNeighbourDetId(const uint32_t detId, const uint32_t direction, const CaloSubdetectorTopology &topo)
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geomToken_
#define DEFINE_FWK_EVENTSETUP_ALPAKA_MODULE(type)
virtual DetId goWest(const DetId &id) const