CMS 3D CMS Logo

MTDThresholdClusterizer.cc
Go to the documentation of this file.
1 // Our own includes
5 
8 
12 
14 
15 // STL
16 #include <stack>
17 #include <vector>
18 #include <iostream>
19 #include <atomic>
20 #include <cmath>
21 using namespace std;
22 
23 //----------------------------------------------------------------------------
25 //----------------------------------------------------------------------------
27  : // Get energy thresholds
28  theHitThreshold(conf.getParameter<double>("HitThreshold")),
29  theSeedThreshold(conf.getParameter<double>("SeedThreshold")),
30  theClusterThreshold(conf.getParameter<double>("ClusterThreshold")),
31  theTimeThreshold(conf.getParameter<double>("TimeThreshold")),
32  thePositionThreshold(conf.getParameter<double>("PositionThreshold")),
33  theNumOfRows(0),
34  theNumOfCols(0),
35  theCurrentId(0),
36  theBuffer(theNumOfRows, theNumOfCols),
37  bufferAlreadySet(false) {}
40 
41 // Configuration descriptions
43  desc.add<double>("HitThreshold", 0.);
44  desc.add<double>("SeedThreshold", 0.);
45  desc.add<double>("ClusterThreshold", 0.);
46  desc.add<double>("TimeThreshold", 10.);
47  desc.add<double>("PositionThreshold", -1.0);
48 }
49 
50 //----------------------------------------------------------------------------
53 //----------------------------------------------------------------------------
54 bool MTDThresholdClusterizer::setup(const MTDGeometry* geom, const MTDTopology* topo, const DetId& id) {
55  theCurrentId = id;
56  //using geopraphicalId here
57  const auto& thedet = geom->idToDet(id);
58  if (thedet == nullptr) {
59  throw cms::Exception("MTDThresholdClusterizer")
60  << "GeographicalID: " << std::hex << id.rawId() << " is invalid!" << std::dec << std::endl;
61  return false;
62  }
63  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
64  const RectangularMTDTopology& topol = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
65 
66  // Get the new sizes.
67  unsigned int nrows = topol.nrows(); // rows in x
68  unsigned int ncols = topol.ncolumns(); // cols in y
69 
70  theNumOfRows = nrows; // Set new sizes
72 
73  LogDebug("MTDThresholdClusterizer") << "Buffer size [" << theNumOfRows + 1 << "," << theNumOfCols + 1 << "]";
74 
75  if (nrows > theBuffer.rows() || ncols > theBuffer.columns()) { // change only when a larger is needed
76  // Resize the buffer
77  theBuffer.setSize(nrows + 1, ncols + 1); // +1 needed for MTD
78  bufferAlreadySet = true;
79  }
80 
81  return true;
82 }
83 //----------------------------------------------------------------------------
89 //----------------------------------------------------------------------------
91  const MTDGeometry* geom,
92  const MTDTopology* topo,
96 
97  // Do not bother for empty detectors
98  if (begin == end) {
99  edm::LogInfo("MTDThresholdClusterizer") << "No hits to clusterize";
100  return;
101  }
102 
103  LogDebug("MTDThresholdClusterizer") << "Input collection " << input.size();
104  assert(output.empty());
105 
106  std::set<unsigned> geoIds;
107  std::multimap<unsigned, unsigned> geoIdToIdx;
108 
109  unsigned index = 0;
110  for (const auto& hit : input) {
111  MTDDetId mtdId = MTDDetId(hit.detid());
112  if (mtdId.subDetector() != MTDDetId::FastTime) {
113  throw cms::Exception("MTDThresholdClusterizer")
114  << "MTDDetId: " << std::hex << mtdId.rawId() << " is invalid!" << std::dec << std::endl;
115  }
116 
117  if (mtdId.mtdSubDetector() == MTDDetId::BTL) {
118  BTLDetId hitId(hit.detid());
119  //for BTL topology gives different layout id
121  geoIdToIdx.emplace(geoId, index);
122  geoIds.emplace(geoId);
123  ++index;
124  } else if (mtdId.mtdSubDetector() == MTDDetId::ETL) {
125  ETLDetId hitId(hit.detid());
126  DetId geoId = hitId.geographicalId();
127  geoIdToIdx.emplace(geoId, index);
128  geoIds.emplace(geoId);
129  ++index;
130  } else {
131  throw cms::Exception("MTDThresholdClusterizer")
132  << "MTDDetId: " << std::hex << mtdId.rawId() << " is invalid!" << std::dec << std::endl;
133  }
134  }
135 
136  //cluster hits within geoIds (modules)
137  for (unsigned id : geoIds) {
138  // Set up the clusterization on this DetId.
139  if (!setup(geom, topo, DetId(id)))
140  return;
141 
142  auto range = geoIdToIdx.equal_range(id);
143  LogDebug("MTDThresholdClusterizer") << "Matching Ids for " << std::hex << id << std::dec << " ["
144  << range.first->second << "," << range.second->second << "]";
145 
146  // Copy MTDRecHits to the buffer array; select the seed hits
147  // on the way, and store them in theSeeds.
148  for (auto itr = range.first; itr != range.second; ++itr) {
149  const unsigned hitidx = itr->second;
150  copy_to_buffer(begin + hitidx, geom, topo);
151  }
152 
153  FTLClusterCollection::FastFiller clustersOnDet(output, id);
154 
155  for (unsigned int i = 0; i < theSeeds.size(); i++) {
156  if (theBuffer.energy(theSeeds[i]) > theSeedThreshold) { // Is this seed still valid?
157  // Make a cluster around this seed
158  const FTLCluster& cluster = make_cluster(theSeeds[i]);
159 
160  // Check if the cluster is above threshold
161  if (cluster.energy() > theClusterThreshold) {
162  LogDebug("MTDThresholdClusterizer")
163  << "putting in this cluster " << i << " #hits:" << cluster.size() << " E:" << cluster.energy()
164  << " T:" << cluster.time() << " X:" << cluster.x() << " Y:" << cluster.y();
165  clustersOnDet.push_back(cluster);
166  }
167  }
168  }
169 
170  // Erase the seeds.
171  theSeeds.clear();
172  // Need to clean unused hits from the buffer array.
173  for (auto itr = range.first; itr != range.second; ++itr) {
174  const unsigned hitidx = itr->second;
175  clear_buffer(begin + hitidx);
176  }
177  }
178 }
179 
180 //----------------------------------------------------------------------------
186 //----------------------------------------------------------------------------
187 void MTDThresholdClusterizer::clear_buffer(RecHitIterator itr) { theBuffer.clear(itr->row(), itr->column()); }
188 
189 //----------------------------------------------------------------------------
191 //----------------------------------------------------------------------------
193  MTDDetId mtdId = MTDDetId(itr->detid());
194  int row = itr->row();
195  int col = itr->column();
197  float energy = itr->energy();
198  float time = itr->time();
199  float timeError = itr->timeError();
200  float position = itr->position();
201  // position is the longitudinal offset that should be added into local x for bars in phi geometry
202  LocalError local_error(0, 0, 0);
203  GlobalPoint global_point(0, 0, 0);
204  if (mtdId.mtdSubDetector() == MTDDetId::BTL) {
206  BTLDetId id = itr->id();
208  const auto& det = geom->idToDet(geoId);
209  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(det->topology());
210  const RectangularMTDTopology& topol = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
211  MeasurementPoint mp(row, col);
212  LocalPoint lp_ctr = topol.localPosition(mp);
213  LocalPoint lp(lp_ctr.x() + position + topol.pitch().first * 0.5f, lp_ctr.y(), lp_ctr.z());
214  // local coordinates of BTL module locates RecHits on the left edge of the bar (-9.2, -3.067, 3.067)
215  // (position + topol.pitch().first/2.0) is the distance from the left edge to the Hit point
216  global_point = det->toGlobal(lp);
217  BTLRecHitsErrorEstimatorIM btlError(det, lp);
218  local_error = btlError.localError();
219  } else if (mtdId.mtdSubDetector() == MTDDetId::ETL) {
221  ETLDetId id = itr->id();
222  DetId geoId = id.geographicalId();
223  const auto& det = geom->idToDet(geoId);
224  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(det->topology());
225  const RectangularMTDTopology& topol = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
226 
227  MeasurementPoint mp(row, col);
228  LocalPoint lp = topol.localPosition(mp);
229  global_point = det->toGlobal(lp);
230  }
231 
232  LogDebug("MTDThresholdClusterizer") << "ROW " << row << " COL " << col << " ENERGY " << energy << " TIME " << time;
233  if (energy > theHitThreshold) {
234  theBuffer.set(row, col, subDet, energy, time, timeError, local_error, global_point);
235  if (energy > theSeedThreshold)
236  theSeeds.push_back(FTLCluster::FTLHitPos(row, col));
237  //sort seeds?
238  }
239 }
240 
241 //----------------------------------------------------------------------------
243 //----------------------------------------------------------------------------
245  //First we acquire the seeds for the clusters
246 
247  GeomDetEnumerators::Location seed_subdet = theBuffer.subDet(hit.row(), hit.col());
248  float seed_energy = theBuffer.energy(hit.row(), hit.col());
249  float seed_time = theBuffer.time(hit.row(), hit.col());
250  float seed_time_error = theBuffer.time_error(hit.row(), hit.col());
251  auto const seedPoint = theBuffer.global_point(hit.row(), hit.col());
252  double seed_error_xx = theBuffer.local_error(hit.row(), hit.col()).xx();
253  double seed_error_yy = theBuffer.local_error(hit.row(), hit.col()).yy();
255 
256  AccretionCluster acluster;
257  acluster.add(hit, seed_energy, seed_time, seed_time_error);
258 
259  bool stopClus = false;
260  //Here we search all hits adjacent to all hits in the cluster.
261  while (!acluster.empty() && !stopClus) {
262  //This is the standard algorithm to find and add a hit
263  auto curInd = acluster.top();
264  acluster.pop();
265  for (auto c = std::max(0, int(acluster.y[curInd] - 1));
266  c < std::min(int(acluster.y[curInd] + 2), int(theBuffer.columns())) && !stopClus;
267  ++c) {
268  for (auto r = std::max(0, int(acluster.x[curInd] - 1));
269  r < std::min(int(acluster.x[curInd] + 2), int(theBuffer.rows())) && !stopClus;
270  ++r) {
271  if (theBuffer.energy(r, c) > theHitThreshold) {
272  if (std::abs(theBuffer.time(r, c) - seed_time) >
274  sqrt(theBuffer.time_error(r, c) * theBuffer.time_error(r, c) + seed_time_error * seed_time_error))
275  continue;
276  if ((seed_subdet == GeomDetEnumerators::barrel) && (theBuffer.subDet(r, c) == GeomDetEnumerators::barrel)) {
277  double hit_error_xx = theBuffer.local_error(r, c).xx();
278  double hit_error_yy = theBuffer.local_error(r, c).yy();
279  if (thePositionThreshold > 0) {
280  if (((theBuffer.global_point(r, c) - seedPoint).mag2()) >
282  (hit_error_xx + seed_error_xx + hit_error_yy + seed_error_yy))
283  continue;
284  }
285  }
286  FTLCluster::FTLHitPos newhit(r, c);
287  if (!acluster.add(newhit, theBuffer.energy(r, c), theBuffer.time(r, c), theBuffer.time_error(r, c))) {
288  stopClus = true;
289  break;
290  }
291  theBuffer.clear(newhit);
292  }
293  }
294  }
295  } // while accretion
296 
297  FTLCluster cluster(theCurrentId,
298  acluster.isize,
299  acluster.energy.data(),
300  acluster.time.data(),
301  acluster.timeError.data(),
302  acluster.x.data(),
303  acluster.y.data(),
304  acluster.xmin,
305  acluster.ymin);
306  return cluster;
307 }
uint8_t geoId(const VFATFrame &frame)
retrieve the GEO information for this channel
void clusterize(const FTLRecHitCollection &input, const MTDGeometry *geom, const MTDTopology *topo, FTLClusterCollection &output) override
Cluster hits. This method operates on a matrix of hits and finds the largest contiguous cluster aroun...
int getMTDTopologyMode() const
Definition: MTDTopology.h:27
uint rows() const
void push_back(data_type const &d)
float energy(uint row, uint col) const
float x() const
Definition: FTLCluster.h:120
int size() const
Definition: FTLCluster.h:141
float time(uint row, uint col) const
LocalError local_error(uint row, uint col) const
static void fillPSetDescription(edm::ParameterSetDescription &desc)
T z() const
Definition: PV3DBase.h:61
uint columns() const
virtual const PixelTopology & specificTopology() const
FTLRecHitCollection::const_iterator RecHitIterator
std::vector< T >::const_iterator const_iterator
int mtdSubDetector() const
Definition: MTDDetId.h:56
const float theHitThreshold
Clustering-related quantities:
int theNumOfRows
Geometry-related information.
MTDThresholdClusterizer(edm::ParameterSet const &conf)
Constructor:
GeomDetEnumerators::Location subDet(uint row, uint col) const
Use subDet to identify whether the Hit is in BTL or ETL.
bool add(FTLCluster::FTLHitPos const &p, float const ienergy, float const itime, float const itimeError)
assert(be >=bs)
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
static std::string const input
Definition: EdmProvDump.cc:50
float time_error(uint row, uint col) const
void copy_to_buffer(RecHitIterator itr, const MTDGeometry *geom, const MTDTopology *topo)
Copy FTLRecHit into the buffer, identify seeds.
float yy() const
Definition: LocalError.h:24
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
float y() const
Definition: FTLCluster.h:125
GlobalPoint global_point(uint row, uint col) const
bool setup(const MTDGeometry *geometry, const MTDTopology *topo, const DetId &id)
std::array< float, MAXSIZE > timeError
int ncolumns() const override
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int nrows() const override
std::vector< FTLCluster::FTLHitPos > theSeeds
std::array< UShort, MAXSIZE > x
std::array< float, MAXSIZE > time
Log< level::Info, false > LogInfo
Definition: DetId.h:17
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
float energy() const
Definition: FTLCluster.h:149
std::pair< float, float > pitch() const override
MTDArrayBuffer theBuffer
Data storage.
FTLCluster make_cluster(const FTLCluster::FTLHitPos &hit)
The actual clustering algorithm: group the neighboring hits around the seed.
void set(uint row, uint col, GeomDetEnumerators::Location subDet, float energy, float time, float time_error, const LocalError &local_error, const GlobalPoint &global_point)
Detector identifier class for the Endcap Timing Layer.
Definition: ETLDetId.h:15
SubDetector subDetector() const
Definition: MTDDetId.h:53
static int position[264][3]
Definition: ReadPGInfo.cc:289
Detector identifier class for the Barrel Timing Layer. The crystal count must start from 0...
Definition: BTLDetId.h:19
col
Definition: cuy.py:1009
void clear_buffer(RecHitIterator itr)
Clear the internal buffer array.
float time() const
Definition: FTLCluster.h:130
BTLDetId::CrysLayout crysLayoutFromTopoMode(const int &topoMode)
std::array< UShort, MAXSIZE > y
LocalPoint localPosition(const MeasurementPoint &mp) const override
void clear(uint row, uint col)
std::array< float, MAXSIZE > energy
float xx() const
Definition: LocalError.h:22
void setSize(uint rows, uint cols)
#define LogDebug(id)