24 using namespace particleFlowRecHitProducer;
26 template <
typename CAL>
30 :
ESProducer(iConfig), cutsFromDB_(iConfig.getParameter<
bool>(
"usePFThresholdsFromDB")) {
31 auto cc = setWhatProduced(
this);
32 geomToken_ =
cc.consumes();
33 if constexpr (std::is_same_v<CAL, HCAL>) {
34 hcalToken_ =
cc.consumes();
35 hcalCutsToken_ =
cc.consumes();
42 desc.add<
bool>(
"usePFThresholdsFromDB",
true);
44 desc.add<
bool>(
"usePFThresholdsFromDB",
false);
48 std::unique_ptr<typename CAL::TopologyTypeHost>
produce(
const typename CAL::TopologyRecordType& iRecord) {
49 const auto&
geom = iRecord.get(geomToken_);
51 auto view = product->view();
52 const int calEnums[2] = {CAL::kSubdetectorBarrelId, CAL::kSubdetectorEndcapId};
53 for (
const auto subdet : calEnums) {
59 std::variant<EcalBarrelTopology, EcalEndcapTopology> topoVar;
61 topo = &iRecord.get(hcalToken_);
68 for (
auto const detId :
geom.getValidDetIds(CAL::kDetectorId, subdet)) {
69 const uint32_t denseId = CAL::detId2denseId(
detId);
70 assert(denseId < CAL::kSize);
73 if constexpr (std::is_same_v<CAL, HCAL>) {
74 view.cutsFromDB() =
false;
76 view.cutsFromDB() =
true;
77 const HcalPFCuts& pfCuts = iRecord.get(hcalCutsToken_);
79 std::unique_ptr<HcalPFCuts>
prod = std::make_unique<HcalPFCuts>(pfCuts);
80 prod->setTopo(&htopo);
82 view.seedThreshold(denseId) =
prod->getValues(
detId.rawId())->seedThreshold();
87 view.positionX(denseId) =
pos.x();
88 view.positionY(denseId) =
pos.y();
89 view.positionZ(denseId) =
pos.z();
91 for (uint32_t
n = 0;
n < 8;
n++) {
92 uint32_t neighDetId = getNeighbourDetId(
detId,
n, *topo);
93 if (CAL::detIdInRange(neighDetId))
94 view.neighbours(denseId)(
n) = CAL::detId2denseId(neighDetId);
96 view.neighbours(denseId)(
n) = 0xffffffff;
102 if constexpr (std::is_same_v<CAL, HCAL>) {
103 for (
const auto subdet : calEnums)
104 for (
auto const detId :
geom.getValidDetIds(CAL::kDetectorId, subdet)) {
105 const uint32_t denseId = CAL::detId2denseId(
detId);
106 for (uint32_t
n = 0;
n < 8;
n++) {
108 view.neighbours(
view.neighbours(denseId)[
n]);
109 if (
std::find(neighboursOfNeighbour.begin(), neighboursOfNeighbour.end(), denseId) ==
110 neighboursOfNeighbour.end())
111 view.neighbours(denseId)[
n] = 0xffffffff;
117 LogDebug(
"PFRecHitTopologyESProducer").log([&](
auto&
log) {
118 for (
const auto subdet : calEnums)
119 for (
const auto detId :
geom.getValidDetIds(CAL::kDetectorId, subdet)) {
120 const uint32_t denseId = CAL::detId2denseId(
detId);
121 log.format(
"detId:{} denseId:{} pos:{},{},{} neighbours:{},{},{},{};{},{},{},{}\n",
124 view[denseId].positionX(),
125 view[denseId].positionY(),
126 view[denseId].positionZ(),
127 view[denseId].neighbours()(0),
128 view[denseId].neighbours()(1),
129 view[denseId].neighbours()(2),
130 view[denseId].neighbours()(3),
131 view[denseId].neighbours()(4),
132 view[denseId].neighbours()(5),
133 view[denseId].neighbours()(6),
134 view[denseId].neighbours()(7));
153 const uint32_t direction,
168 if (direction == 4) {
169 const uint32_t NE = getNeighbourDetId(getNeighbourDetId(
detId, 0, topo), 2, topo);
172 return getNeighbourDetId(getNeighbourDetId(
detId, 2, topo), 0, topo);
174 if (direction == 5) {
175 const uint32_t SW = getNeighbourDetId(getNeighbourDetId(
detId, 1, topo), 3, topo);
178 return getNeighbourDetId(getNeighbourDetId(
detId, 3, topo), 1, topo);
180 if (direction == 6) {
181 const uint32_t ES = getNeighbourDetId(getNeighbourDetId(
detId, 2, topo), 1, topo);
184 return getNeighbourDetId(getNeighbourDetId(
detId, 1, topo), 2, topo);
186 if (direction == 7) {
187 const uint32_t WN = getNeighbourDetId(getNeighbourDetId(
detId, 3, topo), 0, topo);
190 return getNeighbourDetId(getNeighbourDetId(
detId, 0, topo), 3, topo);
197 const uint32_t direction,
213 if (direction == 4) {
218 }
else if (direction == 5) {
223 }
else if (direction == 6) {
228 }
else if (direction == 7) {
236 const uint32_t nn2 = getNeighbourDetId(
detId,
directions.second, topo);
237 const uint32_t nnn = getNeighbourDetId(nn1,
directions.second, topo);
238 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
edm::ESGetToken< HcalPFCuts, HcalPFCutsRcd > hcalCutsToken_