CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFHCALDenseIdNavigator.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFHCALDenseIdNavigator_h
2 #define RecoParticleFlow_PFClusterProducer_PFHCALDenseIdNavigator_h
3 
5 
10 
16 
21 
24 
25 template <typename DET, typename TOPO, bool ownsTopo = true>
27 public:
29  if (!ownsTopo) {
30  topology_.release();
31  }
32  }
33 
35  : vhcalEnum_(iConfig.getParameter<std::vector<int>>("hcalEnums")),
36  hcalToken_(cc.esConsumes<edm::Transition::BeginLuminosityBlock>()),
37  geomToken_(cc.esConsumes<edm::Transition::BeginLuminosityBlock>()) {}
38 
39  void init(const edm::EventSetup& iSetup) override {
40  bool check = theRecNumberWatcher_.check(iSetup);
41  if (!check)
42  return;
43 
44  edm::ESHandle<HcalTopology> hcalTopology = iSetup.getHandle(hcalToken_);
45  topology_.release();
46  topology_.reset(hcalTopology.product());
47 
48  // Fill a vector of valid denseid's
50  const CaloGeometry& caloGeom = *hGeom;
51 
52  std::vector<DetId> vecHcal;
53  std::vector<unsigned int> vDenseIdHcal;
54  neighboursHcal_.clear();
55  for (auto hcalSubdet : vhcalEnum_) {
56  std::vector<DetId> vecDetIds(caloGeom.getValidDetIds(DetId::Hcal, hcalSubdet));
57  vecHcal.insert(vecHcal.end(), vecDetIds.begin(), vecDetIds.end());
58  }
59  vDenseIdHcal.reserve(vecHcal.size());
60  for (auto hDetId : vecHcal) {
61  vDenseIdHcal.push_back(topology_.get()->detId2denseId(hDetId));
62  }
63  std::sort(vDenseIdHcal.begin(), vDenseIdHcal.end());
64 
65  // Fill a vector of cell neighbours
66  denseIdHcalMax_ = *max_element(vDenseIdHcal.begin(), vDenseIdHcal.end());
67  denseIdHcalMin_ = *min_element(vDenseIdHcal.begin(), vDenseIdHcal.end());
68  neighboursHcal_.resize(denseIdHcalMax_ - denseIdHcalMin_ + 1);
69 
70  for (auto denseid : vDenseIdHcal) {
71  DetId N(0);
72  DetId E(0);
73  DetId S(0);
74  DetId W(0);
75  DetId NW(0);
76  DetId NE(0);
77  DetId SW(0);
78  DetId SE(0);
79  std::vector<DetId> neighbours(9, DetId(0));
80 
81  // the centre
82  unsigned denseid_c = denseid;
83  DetId detid_c = topology_.get()->denseId2detId(denseid_c);
84  CaloNavigator<DET> navigator(detid_c, topology_.get());
85 
86  // Using enum in Geometry/CaloTopology/interface/CaloDirection.h
87  // Order: CENTER(NONE),SOUTH,SOUTHEAST,SOUTHWEST,EAST,WEST,NORTHEAST,NORTHWEST,NORTH
88  neighbours.at(NONE) = detid_c;
89 
90  navigator.home();
91  N = navigator.north();
92  neighbours.at(NORTH) = N;
93  if (N != DetId(0)) {
94  NE = navigator.east();
95  } else {
96  navigator.home();
97  E = navigator.east();
98  NE = navigator.north();
99  }
100  neighbours.at(NORTHEAST) = NE;
101 
102  navigator.home();
103  S = navigator.south();
104  neighbours.at(SOUTH) = S;
105  if (S != DetId(0)) {
106  SW = navigator.west();
107  } else {
108  navigator.home();
109  W = navigator.west();
110  SW = navigator.south();
111  }
112  neighbours.at(SOUTHWEST) = SW;
113 
114  navigator.home();
115  E = navigator.east();
116  neighbours.at(EAST) = E;
117  if (E != DetId(0)) {
118  SE = navigator.south();
119  } else {
120  navigator.home();
121  S = navigator.south();
122  SE = navigator.east();
123  }
124  neighbours.at(SOUTHEAST) = SE;
125 
126  navigator.home();
127  W = navigator.west();
128  neighbours.at(WEST) = W;
129  if (W != DetId(0)) {
130  NW = navigator.north();
131  } else {
132  navigator.home();
133  N = navigator.north();
134  NW = navigator.west();
135  }
136  neighbours.at(NORTHWEST) = NW;
137 
138  unsigned index = getIdx(denseid_c);
139  neighboursHcal_[index] = neighbours;
140  }
141  }
142 
144  std::unique_ptr<reco::PFRecHitCollection>& hits,
145  edm::RefProd<reco::PFRecHitCollection>& refProd) override {
146  DetId detid(hit.detId());
147  unsigned denseid = topology_.get()->detId2denseId(detid);
148 
149  std::vector<DetId> neighbours(9, DetId(0));
150 
151  if (denseid < denseIdHcalMin_ || denseid > denseIdHcalMax_) {
152  edm::LogWarning("PFRecHitHCALCachedNavigator") << " DenseId for this cell is out of the range." << std::endl;
153  } else if (!validNeighbours(denseid)) {
154  edm::LogWarning("PFRecHitHCALCachedNavigator")
155  << " DenseId for this cell does not have the neighbour information." << std::endl;
156  } else {
157  unsigned index = getIdx(denseid);
158  neighbours = neighboursHcal_.at(index);
159  }
160 
161  associateNeighbour(neighbours.at(NORTH), hit, hits, refProd, 0, 1, 0); // N
162  associateNeighbour(neighbours.at(NORTHEAST), hit, hits, refProd, 1, 1, 0); // NE
163  associateNeighbour(neighbours.at(SOUTH), hit, hits, refProd, 0, -1, 0); // S
164  associateNeighbour(neighbours.at(SOUTHWEST), hit, hits, refProd, -1, -1, 0); // SW
165  associateNeighbour(neighbours.at(EAST), hit, hits, refProd, 1, 0, 0); // E
166  associateNeighbour(neighbours.at(SOUTHEAST), hit, hits, refProd, 1, -1, 0); // SE
167  associateNeighbour(neighbours.at(WEST), hit, hits, refProd, -1, 0, 0); // W
168  associateNeighbour(neighbours.at(NORTHWEST), hit, hits, refProd, -1, 1, 0); // NW
169  }
170 
171  bool validNeighbours(const unsigned int denseid) const {
172  bool ok = true;
173  unsigned index = getIdx(denseid);
174  if (neighboursHcal_.at(index).size() != 9)
175  ok = false; // the neighbour vector size should be 3x3
176  return ok;
177  }
178 
179  unsigned int getIdx(const unsigned int denseid) const {
180  unsigned index = denseid - denseIdHcalMin_;
181  return index;
182  }
183 
184 protected:
186  std::unique_ptr<const TOPO> topology_;
187  std::vector<int> vhcalEnum_;
188  std::vector<std::vector<DetId>> neighboursHcal_;
189  unsigned int denseIdHcalMax_;
190  unsigned int denseIdHcalMin_;
191 
192 private:
195 };
196 
197 #endif
void associateNeighbour(const DetId &id, reco::PFRecHit &hit, std::unique_ptr< reco::PFRecHitCollection > &hits, edm::RefProd< reco::PFRecHitCollection > &refProd, short eta, short phi, short depth)
unsigned detId() const
rechit detId
Definition: PFRecHit.h:93
void associateNeighbours(reco::PFRecHit &hit, std::unique_ptr< reco::PFRecHitCollection > &hits, edm::RefProd< reco::PFRecHitCollection > &refProd) override
unsigned int getIdx(const unsigned int denseid) const
std::unique_ptr< const TOPO > topology_
edm::ESWatcher< HcalRecNumberingRecord > theRecNumberWatcher_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geomToken_
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
Transition
Definition: Transition.h:12
std::vector< std::vector< DetId > > neighboursHcal_
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > hcalToken_
Definition: DetId.h:17
#define N
Definition: blowfish.cc:9
PFHCALDenseIdNavigator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &cc)
T const * product() const
Definition: ESHandle.h:86
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:97
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
std::vector< DetId > getValidDetIds() const
Get the list of all valid detector ids.
Definition: CaloGeometry.cc:75
void init(const edm::EventSetup &iSetup) override
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
Log< level::Warning, false > LogWarning
bool validNeighbours(const unsigned int denseid) const
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283