CMS 3D CMS Logo

MTDThresholdClusterizer.cc
Go to the documentation of this file.
1 // Our own includes
4 
7 
10 
11 // STL
12 #include <stack>
13 #include <vector>
14 #include <iostream>
15 #include <atomic>
16 using namespace std;
17 
18 //#define DEBUG_ENABLED
19 #ifdef DEBUG_ENABLED
20 #define DEBUG(x) \
21  do { \
22  std::cout << x << std::endl; \
23  } while (0)
24 #else
25 #define DEBUG(x)
26 #endif
27 
28 //----------------------------------------------------------------------------
30 //----------------------------------------------------------------------------
32  : // Get energy thresholds
33  theHitThreshold(conf.getParameter<double>("HitThreshold")),
34  theSeedThreshold(conf.getParameter<double>("SeedThreshold")),
35  theClusterThreshold(conf.getParameter<double>("ClusterThreshold")),
36  theTimeThreshold(conf.getParameter<double>("TimeThreshold")),
37  theNumOfRows(0),
38  theNumOfCols(0),
39  theCurrentId(0),
40  theBuffer(theNumOfRows, theNumOfCols),
41  bufferAlreadySet(false) {}
44 
45 // Configuration descriptions
47  desc.add<double>("HitThreshold", 0.);
48  desc.add<double>("SeedThreshold", 0.);
49  desc.add<double>("ClusterThreshold", 0.);
50  desc.add<double>("TimeThreshold", 10.);
51 }
52 
53 //----------------------------------------------------------------------------
56 //----------------------------------------------------------------------------
57 bool MTDThresholdClusterizer::setup(const MTDGeometry* geom, const MTDTopology* topo, const DetId& id) {
58  theCurrentId = id;
59  //using geopraphicalId here
60  const auto& thedet = geom->idToDet(id);
61  if (thedet == nullptr) {
62  throw cms::Exception("MTDThresholdClusterizer")
63  << "GeographicalID: " << std::hex << id.rawId() << " is invalid!" << std::dec << std::endl;
64  return false;
65  }
66  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
67  const RectangularMTDTopology& topol = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
68 
69  // Get the new sizes.
70  unsigned int nrows = topol.nrows(); // rows in x
71  unsigned int ncols = topol.ncolumns(); // cols in y
72 
73  theNumOfRows = nrows; // Set new sizes
75 
76  DEBUG("Buffer size [" << theNumOfRows + 1 << "," << theNumOfCols + 1 << "]");
77 
78  if (nrows > theBuffer.rows() || ncols > theBuffer.columns()) { // change only when a larger is needed
79  // Resize the buffer
80  theBuffer.setSize(nrows + 1, ncols + 1); // +1 needed for MTD
81  bufferAlreadySet = true;
82  }
83 
84  return true;
85 }
86 //----------------------------------------------------------------------------
92 //----------------------------------------------------------------------------
94  const MTDGeometry* geom,
95  const MTDTopology* topo,
99 
100  // Do not bother for empty detectors
101  if (begin == end) {
102  edm::LogInfo("MTDThresholdClusterizer") << "No hits to clusterize";
103  return;
104  }
105 
106  DEBUG("Input collection " << input.size());
107  assert(output.empty());
108 
109  std::set<unsigned> geoIds;
110  std::multimap<unsigned, unsigned> geoIdToIdx;
111 
112  unsigned index = 0;
113  for (const auto& hit : input) {
114  MTDDetId mtdId = MTDDetId(hit.detid());
115  if (mtdId.subDetector() != MTDDetId::FastTime) {
116  throw cms::Exception("MTDThresholdClusterizer")
117  << "MTDDetId: " << std::hex << mtdId.rawId() << " is invalid!" << std::dec << std::endl;
118  }
119 
120  if (mtdId.mtdSubDetector() == MTDDetId::BTL) {
121  BTLDetId hitId(hit.detid());
122  DetId geoId = hitId.geographicalId(
123  (BTLDetId::CrysLayout)topo->getMTDTopologyMode()); //for BTL topology gives different layout id
124  geoIdToIdx.emplace(geoId, index);
125  geoIds.emplace(geoId);
126  ++index;
127  } else if (mtdId.mtdSubDetector() == MTDDetId::ETL) {
128  ETLDetId hitId(hit.detid());
129  DetId geoId = hitId.geographicalId();
130  geoIdToIdx.emplace(geoId, index);
131  geoIds.emplace(geoId);
132  ++index;
133  } else {
134  throw cms::Exception("MTDThresholdClusterizer")
135  << "MTDDetId: " << std::hex << mtdId.rawId() << " is invalid!" << std::dec << std::endl;
136  }
137  }
138 
139  //cluster hits within geoIds (modules)
140  for (unsigned id : geoIds) {
141  // Set up the clusterization on this DetId.
142  if (!setup(geom, topo, DetId(id)))
143  return;
144 
145  auto range = geoIdToIdx.equal_range(id);
146  DEBUG("Matching Ids for " << std::hex << id << std::dec << " [" << range.first->second << ","
147  << range.second->second << "]");
148 
149  // Copy MTDRecHits to the buffer array; select the seed hits
150  // on the way, and store them in theSeeds.
151  for (auto itr = range.first; itr != range.second; ++itr) {
152  const unsigned hitidx = itr->second;
153  copy_to_buffer(begin + hitidx);
154  }
155 
156  FTLClusterCollection::FastFiller clustersOnDet(output, id);
157 
158  for (unsigned int i = 0; i < theSeeds.size(); i++) {
159  if (theBuffer.energy(theSeeds[i]) > theSeedThreshold) { // Is this seed still valid?
160  // Make a cluster around this seed
161  const FTLCluster& cluster = make_cluster(theSeeds[i]);
162 
163  // Check if the cluster is above threshold
164  if (cluster.energy() > theClusterThreshold) {
165  DEBUG("putting in this cluster " << i << " #hits:" << cluster.size() << " E:" << cluster.energy()
166  << " T:" << cluster.time() << " X:" << cluster.x() << " Y:" << cluster.y());
167  clustersOnDet.push_back(cluster);
168  }
169  }
170  }
171 
172  // Erase the seeds.
173  theSeeds.clear();
174  // Need to clean unused hits from the buffer array.
175  for (auto itr = range.first; itr != range.second; ++itr) {
176  const unsigned hitidx = itr->second;
177  clear_buffer(begin + hitidx);
178  }
179  }
180 }
181 
182 //----------------------------------------------------------------------------
188 //----------------------------------------------------------------------------
189 void MTDThresholdClusterizer::clear_buffer(RecHitIterator itr) { theBuffer.clear(itr->row(), itr->column()); }
190 
191 //----------------------------------------------------------------------------
193 //----------------------------------------------------------------------------
195  int row = itr->row();
196  int col = itr->column();
197  float energy = itr->energy();
198  float time = itr->time();
199  float timeError = itr->timeError();
200 
201  DEBUG("ROW " << row << " COL " << col << " ENERGY " << energy << " TIME " << time);
202  if (energy > theHitThreshold) {
203  theBuffer.set(row, col, energy, time, timeError);
204  if (energy > theSeedThreshold)
205  theSeeds.push_back(FTLCluster::FTLHitPos(row, col));
206  //sort seeds?
207  }
208 }
209 
210 //----------------------------------------------------------------------------
212 //----------------------------------------------------------------------------
214  //First we acquire the seeds for the clusters
215  float seed_energy = theBuffer.energy(hit.row(), hit.col());
216  float seed_time = theBuffer.time(hit.row(), hit.col());
217  float seed_time_error = theBuffer.time_error(hit.row(), hit.col());
218  theBuffer.clear(hit);
219 
220  AccretionCluster acluster;
221  acluster.add(hit, seed_energy, seed_time, seed_time_error);
222 
223  bool stopClus = false;
224  //Here we search all hits adjacent to all hits in the cluster.
225  while (!acluster.empty() && !stopClus) {
226  //This is the standard algorithm to find and add a hit
227  auto curInd = acluster.top();
228  acluster.pop();
229  for (auto c = std::max(0, int(acluster.y[curInd] - 1));
230  c < std::min(int(acluster.y[curInd] + 2), int(theBuffer.columns())) && !stopClus;
231  ++c) {
232  for (auto r = std::max(0, int(acluster.x[curInd] - 1));
233  r < std::min(int(acluster.x[curInd] + 2), int(theBuffer.rows())) && !stopClus;
234  ++r) {
235  if (theBuffer.energy(r, c) > theHitThreshold) {
236  if (std::abs(theBuffer.time(r, c) - seed_time) >
238  sqrt(theBuffer.time_error(r, c) * theBuffer.time_error(r, c) + seed_time_error * seed_time_error))
239  continue;
240  FTLCluster::FTLHitPos newhit(r, c);
241  if (!acluster.add(newhit, theBuffer.energy(r, c), theBuffer.time(r, c), theBuffer.time_error(r, c))) {
242  stopClus = true;
243  break;
244  }
245  theBuffer.clear(newhit);
246  }
247  }
248  }
249  } // while accretion
250 
251  FTLCluster cluster(theCurrentId,
252  acluster.isize,
253  acluster.energy.data(),
254  acluster.time.data(),
255  acluster.timeError.data(),
256  acluster.x.data(),
257  acluster.y.data(),
258  acluster.xmin,
259  acluster.ymin);
260  return cluster;
261 }
float time_error(uint row, uint col) const
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...
void push_back(data_type const &d)
uint rows() const
float theHitThreshold
Clustering-related quantities:
float y() const
Definition: FTLCluster.h:126
CrysLayout
Definition: BTLDetId.h:70
FTLRecHitCollection::const_iterator RecHitIterator
uint columns() const
int getMTDTopologyMode() const
Definition: MTDTopology.h:73
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< T >::const_iterator const_iterator
float time() const
Definition: FTLCluster.h:131
int ncolumns() const override
int theNumOfRows
Geometry-related information.
MTDThresholdClusterizer(edm::ParameterSet const &conf)
Constructor:
bool add(FTLCluster::FTLHitPos const &p, float const ienergy, float const itime, float const itimeError)
int nrows() const override
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
static std::string const input
Definition: EdmProvDump.cc:48
static void fillDescriptions(edm::ParameterSetDescription &desc)
float time(uint row, uint col) const
const MTDGeomDet * idToDet(DetId) const override
Definition: MTDGeometry.cc:160
bool setup(const MTDGeometry *geometry, const MTDTopology *topo, const DetId &id)
std::array< float, MAXSIZE > timeError
float x() const
Definition: FTLCluster.h:121
T sqrt(T t)
Definition: SSEVec.h:19
virtual const PixelTopology & specificTopology() const
constexpr int row() const
Definition: FTLCluster.h:62
float energy(uint row, uint col) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define end
Definition: vmac.h:39
T min(T a, T b)
Definition: MathUtil.h:58
std::vector< FTLCluster::FTLHitPos > theSeeds
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void copy_to_buffer(RecHitIterator itr)
Copy FTLRecHit into the buffer, identify seeds.
int size() const
Definition: FTLCluster.h:142
std::array< UShort, MAXSIZE > x
std::array< float, MAXSIZE > time
const_iterator end() const
Definition: DetId.h:17
float energy() const
Definition: FTLCluster.h:150
#define DEBUG(x)
MTDArrayBuffer theBuffer
Data storage.
SubDetector subDetector() const
Definition: MTDDetId.h:53
FTLCluster make_cluster(const FTLCluster::FTLHitPos &hit)
The actual clustering algorithm: group the neighboring hits around the seed.
Detector identifier class for the Endcap Timing Layer.
Definition: ETLDetId.h:15
#define begin
Definition: vmac.h:32
size_type size() const
Detector identifier class for the Barrel Timing Layer. The crystal count must start from 0...
Definition: BTLDetId.h:18
int mtdSubDetector() const
Definition: MTDDetId.h:56
void set(uint row, uint col, float energy, float time, float time_error)
col
Definition: cuy.py:1010
void clear_buffer(RecHitIterator itr)
Clear the internal buffer array.
std::array< UShort, MAXSIZE > y
constexpr int col() const
Definition: FTLCluster.h:63
void clear(uint row, uint col)
std::array< float, MAXSIZE > energy
void setSize(uint rows, uint cols)
const_iterator begin() const