All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MergeClusterAlg.cxx
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////
2 // Merge Cluster algorithm
3 //
4 // Runs on the output of previous clustering algorithms to merge
5 // clusters together which lie on a straight line and are within
6 // some separation threshold.
7 // Runs recursively over all clusters, including new ones formed
8 // in the algorithm.
9 //
10 // M Wallbank (m.wallbank@sheffield.ac.uk), July 2015
11 ////////////////////////////////////////////////////////////////////
12 
14 
15 #include "TPrincipal.h"
16 #include "TTree.h"
17 
18 cluster::MergeClusterAlg::MergeClusterAlg(fhicl::ParameterSet const& pset) {
19  this->reconfigure(pset);
20  fTree = tfs->make<TTree>("MatchingVariables","MatchingVariables");
21  fTree->Branch("Angle",&fAngle);
22  fTree->Branch("Eigenvalue",&fEigenvalue);
23  fTree->Branch("Cluster1Size",&fCluster1Size);
24  fTree->Branch("Cluster2Size",&fCluster2Size);
25  fTree->Branch("Length1",&fLength1);
26  fTree->Branch("Length2",&fLength2);
27  fTree->Branch("Separation",&fSeparation);
28  fTree->Branch("CrossingDistance",&fCrossingDistance);
29  fTree->Branch("ProjectedWidth",&fProjectedWidth);
30  fTree->Branch("Overlap",&fOverlap);
31  fTree->Branch("TrueMerge",&fTrueMerge);
32 }
33 
34 void cluster::MergeClusterAlg::FindClusterEndPoints(art::PtrVector<recob::Hit> const& cluster, TVector2 const& centre, TVector2 const& direction, TVector2& start, TVector2& end) const {
35 
36  /// Find estimates of cluster start/end points
37 
38  TVector2 pos;
39  std::map<double,TVector2> hitProjection;
40 
41  // Project all hits onto line to determine end points
42  for (auto &hit : cluster) {
43  pos = HitCoordinates(hit) - centre;
44  hitProjection[direction*pos] = pos;
45  }
46 
47  // Project end points onto line which passes through centre of cluster
48  start = hitProjection.begin()->second.Proj(direction) + centre;
49  end = hitProjection.rbegin()->second.Proj(direction) + centre;
50 
51  return;
52 
53 }
54 
55 double cluster::MergeClusterAlg::FindClusterOverlap(TVector2 const& direction, TVector2 const& centre, TVector2 const& start1, TVector2 const& end1, TVector2 const& start2, TVector2 const& end2) const {
56 
57  /// Calculates the overlap of the clusters on the line projected between them
58 
59  double clusterOverlap = 0;
60 
61  // Project onto the average direction through both clusters
62  double s1 = (start1-centre)*direction;
63  double e1 = (end1-centre)*direction;
64  double s2 = (start2-centre)*direction;
65  double e2 = (end2-centre)*direction;
66 
67  // Make sure end > start
68  if (s1 > e1) {
69  std::cout << "s1>e1: " << s1 << " and " << e1 << std::endl;
70  double tmp = e1;
71  e1 = s1;
72  s1 = tmp;
73  }
74  if (s2 > e2) {
75  std::cout << "s1>e1: " << s1 << " and " << e1 << std::endl;
76  double tmp = e2;
77  e2 = s2;
78  s2 = tmp;
79  }
80 
81  // Find the overlap of the clusters on the centre line
82  if ((e1 > s2) && (e2 > s1))
83  clusterOverlap = std::min((e1 - s2), (e2 - s1));
84 
85  return clusterOverlap;
86 
87 }
88 
89 double cluster::MergeClusterAlg::FindCrossingDistance(TVector2 const &direction1, TVector2 const &centre1, TVector2 const &direction2, TVector2 const &centre2) const {
90 
91  /// Finds the distance between the crossing point of the lines and the closest line centre
92 
93  // Find intersection point of two lines drawn through the centre of the clusters
94  double dcross = (direction1.X() * direction2.Y()) - (direction1.Y() * direction2.X());
95  TVector2 p = centre2 - centre1;
96  double pcrossd = (p.X() * direction2.Y()) - (p.Y() * direction2.X());
97  TVector2 crossing = centre1 + ((pcrossd/dcross) * direction1);
98 
99  // Get distance from this point to the clusters
100  double crossingDistance = std::min((centre1-crossing).Mod(),(centre2-crossing).Mod());
101 
102  return crossingDistance;
103 
104 }
105 
106 double cluster::MergeClusterAlg::FindMinSeparation(art::PtrVector<recob::Hit> const& cluster1, art::PtrVector<recob::Hit> const& cluster2) const {
107 
108  /// Calculates the minimum separation between two clusters
109 
110  double minDistance = 99999.;
111 
112  // Loop over the two clusters to find the smallest distance
113  for (auto const& hit1 : cluster1) {
114  for (auto const& hit2 : cluster2) {
115 
116  TVector2 pos1 = HitCoordinates(hit1);
117  TVector2 pos2 = HitCoordinates(hit2);
118 
119  double distance = (pos1 - pos2).Mod();
120 
121  if (distance < minDistance) minDistance = distance;
122 
123  }
124  }
125 
126  return minDistance;
127 
128 }
129 
130 double cluster::MergeClusterAlg::FindProjectedWidth(TVector2 const& centre1, TVector2 const& start1, TVector2 const& end1, TVector2 const& centre2, TVector2 const& start2, TVector2 const& end2) const {
131 
132  /// Projects clusters parallel to the line which runs through their centres and finds the minimum containing width
133 
134  // Get the line running through the centre of the two clusters
135  TVector2 parallel = (centre2 - centre1).Unit();
136  TVector2 perpendicular = parallel.Rotate(TMath::Pi()/2);
137 
138  // Project the cluster vector onto this perpendicular line
139  double s1 = (start1-centre1)*perpendicular;
140  double e1 = (end1-centre1)*perpendicular;
141  double s2 = (start2-centre2)*perpendicular;
142  double e2 = (end2-centre2)*perpendicular;
143 
144  // Find the width in each direction
145  double projectionStart = std::max(TMath::Abs(s1), TMath::Abs(s2));
146  double projectionEnd = std::max(TMath::Abs(e1), TMath::Abs(e2));
147 
148  double projectionWidth = projectionStart + projectionEnd;
149 
150  return projectionWidth;
151 
152 }
153 
155 
156  /// Find the global wire position
157 
158  double wireCentre[3];
159  fGeom->WireIDToWireGeo(wireID).GetCenter(wireCentre);
160 
161  double globalWire;
162  if (fGeom->SignalType(wireID) == geo::kInduction) {
163  if (wireID.TPC % 2 == 0) globalWire = fGeom->WireCoordinate(wireCentre[1], wireCentre[2], wireID.Plane, 0, wireID.Cryostat);
164  else globalWire = fGeom->WireCoordinate(wireCentre[1], wireCentre[2], wireID.Plane, 1, wireID.Cryostat);
165  }
166  else {
167  if (wireID.TPC % 2 == 0) globalWire = wireID.Wire + ((wireID.TPC/2) * fGeom->Nwires(wireID.Plane, 0, wireID.Cryostat));
168  else globalWire = wireID.Wire + ((int)(wireID.TPC/2) * fGeom->Nwires(wireID.Plane, 1, wireID.Cryostat));
169  }
170 
171  return globalWire;
172 
173 }
174 
175 TVector2 cluster::MergeClusterAlg::HitCoordinates(art::Ptr<recob::Hit> const& hit) const {
176 
177  /// Return the coordinates of this hit in global wire/tick space
178 
179  return TVector2(GlobalWire(hit->WireID()), hit->PeakTime());
180 
181 }
182 
183 int cluster::MergeClusterAlg::MergeClusters(std::vector<art::PtrVector<recob::Hit> > const &planeClusters, std::vector<art::PtrVector<recob::Hit> > &clusters) const {
184 
185  /// Merges clusters which lie along a straight line
186 
187  // // ////// MAKE SOME MESSY CODE! CHECK FEATURES OF THESE CLUSTERS
188 
189  // // Truth matching
190  // for (unsigned int cluster = 0; cluster < planeClusters.size(); ++cluster) {
191  // std::map<int,double> trackMap;
192  // art::PtrVector<recob::Hit> hits = planeClusters.at(cluster);
193  // for (auto &hit : hits) {
194  // std::vector<sim::TrackIDE> ides = backtracker->HitToTrackID(hit);
195  // for (auto &ide : ides)
196  // trackMap[ide.trackID] += ide.energy;
197  // }
198  // // Find the true particle associated with this track
199  // double highEnergy = 0;
200  // int bestTrack = 0;
201  // for (auto &track : trackMap) {
202  // if (track.second > highEnergy) {
203  // highEnergy = track.second;
204  // bestTrack = track.first;
205  // }
206  // }
207  // trueClusterMap[cluster] = bestTrack;
208  // }
209 
210  // for (unsigned int cluster1It = 0; cluster1It < planeClusters.size(); ++cluster1It) {
211  // for (unsigned int cluster2It = cluster1It+1; cluster2It < planeClusters.size(); ++cluster2It) {
212 
213  // const art::PtrVector<recob::Hit> cluster1 = planeClusters.at(cluster1It);
214  // const art::PtrVector<recob::Hit> cluster2 = planeClusters.at(cluster2It);
215 
216  // // true merge
217  // if (trueClusterMap[cluster1It] == trueClusterMap[cluster2It])
218  // fTrueMerge = true;
219  // else fTrueMerge = false;
220 
221  // // geometry
222  // fCluster1Size = cluster1.size();
223  // fCluster2Size = cluster2.size();
224  // fSeparation = this->FindMinSeparation(cluster1, cluster2);
225 
226  // // PCA
227  // TPrincipal *pca = new TPrincipal(2,"");
228  // TPrincipal *pca1 = new TPrincipal(2,"");
229  // TPrincipal *pca2 = new TPrincipal(2,"");
230  // double hits[2];
231  // TVector2 pos;
232 
233  // // Cluster centre
234  // TVector2 chargePoint1 = TVector2(0,0), chargePoint2 = TVector2(0,0);
235  // double totalCharge1 = 0, totalCharge2 = 0;
236 
237  // for (auto &hit1 : cluster1) {
238  // pos = HitCoordinates(hit1);
239  // hits[0] = pos.X();
240  // hits[1] = pos.Y();
241  // pca->AddRow(hits);
242  // pca1->AddRow(hits);
243  // chargePoint1 += hit1->Integral() * pos;
244  // totalCharge1 += hit1->Integral();
245  // }
246  // for (auto &hit2 : cluster2) {
247  // pos = HitCoordinates(hit2);
248  // hits[0] = pos.X();
249  // hits[1] = pos.Y();
250  // pca->AddRow(hits);
251  // pca2->AddRow(hits);
252  // chargePoint2 += hit2->Integral() * pos;
253  // totalCharge2 += hit2->Integral();
254  // }
255 
256  // pca->MakePrincipals();
257  // pca1->MakePrincipals();
258  // pca2->MakePrincipals();
259 
260  // // Properties of these clusters
261  // TVector2 direction1 = TVector2( (*pca1->GetEigenVectors())[0][0], (*pca1->GetEigenVectors())[1][0] ).Unit();
262  // TVector2 direction2 = TVector2( (*pca2->GetEigenVectors())[0][0], (*pca2->GetEigenVectors())[1][0] ).Unit();
263  // TVector2 directionAv = ((direction1+direction2)/2).Unit();
264  // TVector2 centre1 = chargePoint1 / totalCharge1;
265  // TVector2 centre2 = chargePoint2 / totalCharge2;
266  // TVector2 centre = (centre1+centre2)/2;
267  // TVector2 start1, end1;
268  // TVector2 start2, end2;
269  // FindClusterEndPoints(cluster1, centre1, direction1, start1, end1);
270  // FindClusterEndPoints(cluster2, centre2, direction2, start2, end2);
271  // fLength1 = (end1-start1).Mod();
272  // fLength2 = (end2-start2).Mod();
273 
274  // // Properties of the pair of clusters
275  // fCrossingDistance = FindCrossingDistance(direction1, centre1, direction2, centre2);
276  // fProjectedWidth = FindProjectedWidth(centre1, start1, end1, centre2, start2, end2);
277  // fAngle = direction1.DeltaPhi(direction2);
278  // if (fAngle > 1.57) fAngle = 3.14159 - fAngle;
279  // fOverlap = FindClusterOverlap(directionAv, centre, start1, end1, start2, end2);
280  // fSeparation = FindMinSeparation(cluster1, cluster2);
281  // fEigenvalue = (*pca->GetEigenValues())[0];
282 
283  // fTree->Fill();
284 
285  // // std::cout << std::endl << "Plane " << fPlane << ": Clusters " << cluster1It << " and " << cluster2It << " have overlap " << fOverlap << " and start and end ... " << std::endl;
286  // // start1.Print();
287  // // end1.Print();
288  // // start2.Print();
289  // // end2.Print();
290 
291  // // // Find if this is merged!
292  // // if (fCrossingDistance < 6 + (5 / (fAngle - 0.05)))
293  // // fMerge = true;
294  // // else fMerge = false;
295 
296  // // if (fCluster1Size >= 10 && fCluster2Size >= 10) std::cout << "Merge " << fMerge << " and true merge " << fTrueMerge << std::endl;
297 
298  // }
299  // }
300 
301  // ----------------------------- END OF MESSY CODE! --------------------------------------------------------------------------------------------------------------
302 
303  std::vector<unsigned int> mergedClusters;
304 
305  std::vector<art::PtrVector<recob::Hit> > oldClusters = planeClusters;
306 
307  // Sort the clusters by size
308  std::sort(oldClusters.begin(), oldClusters.end(), [](const art::PtrVector<recob::Hit> &a, const art::PtrVector<recob::Hit> &b) {return a.size() > b.size();} );
309 
310  // Find the numbers of clusters above size threshold
311  unsigned int nclusters = 0;
312  for (auto &cluster : oldClusters)
313  if (cluster.size() >= fMinMergeClusterSize) ++nclusters;
314 
315  // Until all clusters are merged, create new clusters
316  bool mergedAllClusters = false;
317  while (!mergedAllClusters) {
318 
319  // New cluster
320  art::PtrVector<recob::Hit> cluster;
321 
322  // Put the largest unmerged cluster in this new cluster
323  for (unsigned int initCluster = 0; initCluster < oldClusters.size(); ++initCluster) {
324  if (oldClusters.at(initCluster).size() < fMinMergeClusterSize or std::find(mergedClusters.begin(), mergedClusters.end(), initCluster) != mergedClusters.end()) continue;
325  cluster = oldClusters.at(initCluster);
326  mergedClusters.push_back(initCluster);
327  break;
328  }
329 
330  // Merge all aligned clusters to this
331  bool mergedAllToThisCluster = false;
332  while (!mergedAllToThisCluster) {
333 
334  // Look at all clusters and merge
335  int nadded = 0;
336  for (unsigned int trialCluster = 0; trialCluster < oldClusters.size(); ++trialCluster) {
337 
338  if (oldClusters.at(trialCluster).size() < fMinMergeClusterSize or std::find(mergedClusters.begin(), mergedClusters.end(), trialCluster) != mergedClusters.end()) continue;
339 
340  // Calculate the PCA for each
341  TPrincipal *pca1 = new TPrincipal(2,""), *pca2 = new TPrincipal(2,"");
342  double hits[2];
343  TVector2 pos;
344 
345  // Cluster centre
346  TVector2 chargePoint1 = TVector2(0,0), chargePoint2 = TVector2(0,0);
347  double totalCharge1 = 0, totalCharge2 = 0;
348 
349  for (auto &hit1 : cluster) {
350  pos = HitCoordinates(hit1);
351  hits[0] = pos.X();
352  hits[1] = pos.Y();
353  pca1->AddRow(hits);
354  chargePoint1 += hit1->Integral() * pos;
355  totalCharge1 += hit1->Integral();
356  }
357  for (auto &hit2 : oldClusters.at(trialCluster)) {
358  pos = HitCoordinates(hit2);
359  hits[0] = pos.X();
360  hits[1] = pos.Y();
361  pca2->AddRow(hits);
362  chargePoint2 += hit2->Integral() * pos;
363  totalCharge2 += hit2->Integral();
364  }
365 
366  pca1->MakePrincipals();
367  pca2->MakePrincipals();
368 
369  // Properties of these clusters
370  TVector2 direction1 = TVector2( (*pca1->GetEigenVectors())[0][0], (*pca1->GetEigenVectors())[1][0] ).Unit();
371  TVector2 direction2 = TVector2( (*pca2->GetEigenVectors())[0][0], (*pca2->GetEigenVectors())[1][0] ).Unit();
372  TVector2 direction = ((direction1+direction2)/2).Unit();
373  TVector2 centre1 = chargePoint1 / totalCharge1;
374  TVector2 centre2 = chargePoint2 / totalCharge2;
375  TVector2 centre = (centre1+centre2)/2;
376  TVector2 start1, end1;
377  TVector2 start2, end2;
378  FindClusterEndPoints(cluster, centre1, direction1, start1, end1);
379  FindClusterEndPoints(oldClusters.at(trialCluster), centre2, direction2, start2, end2);
380  double length1 = (end1-start1).Mod();
381  double length2 = (end2-start2).Mod();
382 
383  // Properties of the pair of clusters
384  double crossingDistance = FindCrossingDistance(direction1, centre1, direction2, centre2);
385  double projectedWidth = FindProjectedWidth(centre1, start1, end1, centre2, start2, end2);
386  double angle = direction1.DeltaPhi(direction2);
387  if (angle > 1.57) angle = 3.14159 - angle;
388  double overlap = FindClusterOverlap(direction, centre, start1, end1, start2, end2);
389  double separation = FindMinSeparation(cluster, oldClusters.at(trialCluster));
390 
391  if (separation > fMaxMergeSeparation)
392  continue;
393  if (PassCuts(angle, crossingDistance, projectedWidth, separation, overlap, TMath::Max(length1, length2))) {
394 
395  for (auto &hit : oldClusters.at(trialCluster))
396  cluster.push_back(hit);
397 
398  mergedClusters.push_back(trialCluster);
399  ++nadded;
400 
401  }
402 
403  delete pca1;
404  delete pca2;
405 
406  } // loop over clusters to add
407 
408  if (nadded == 0) mergedAllToThisCluster = true;
409 
410  } // while loop
411 
412  clusters.push_back(cluster);
413  if (mergedClusters.size() == nclusters) mergedAllClusters = true;
414 
415  }
416 
417  return clusters.size();
418 
419 }
420 
421 bool cluster::MergeClusterAlg::PassCuts(double const& angle, double const& crossingDistance, double const& projectedWidth, double const& separation, double const& overlap, double const& longLength) const {
422 
423  /// Boolean function which decides whether or not two clusters should be merged, depending on their properties
424 
425  bool passCrossingDistanceAngle = false;
426  if (crossingDistance < (-2 + (5 / (1 * TMath::Abs(angle)) - 0) ) )
427  passCrossingDistanceAngle = true;
428 
429  bool passSeparationAngle = false;
430  if (separation < (200 * TMath::Abs(angle) + 40))
431  passSeparationAngle = true;
432 
433  bool passProjectedWidth = false;
434  if (((double)projectedWidth/(double)longLength) < fProjWidthThreshold)
435  passProjectedWidth = true;
436 
437  return passCrossingDistanceAngle and passSeparationAngle and passProjectedWidth;
438 
439 }
440 
441 void cluster::MergeClusterAlg::reconfigure(fhicl::ParameterSet const& p) {
442  fMinMergeClusterSize = p.get<int> ("MinMergeClusterSize");
443  fMaxMergeSeparation = p.get<double>("MaxMergeSeparation");
444  fProjWidthThreshold = p.get<double>("ProjWidthThreshold");
445 }
MergeClusterAlg(fhicl::ParameterSet const &pset)
process_name cluster
Definition: cheaterreco.fcl:51
pdgs p
Definition: selectors.fcl:22
void FindClusterEndPoints(art::PtrVector< recob::Hit > const &cluster, TVector2 const &centre, TVector2 const &direction, TVector2 &start, TVector2 &end) const
art::ServiceHandle< art::TFileService const > tfs
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:212
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:580
process_name hit
Definition: cheaterreco.fcl:51
void reconfigure(fhicl::ParameterSet const &p)
double FindClusterOverlap(TVector2 const &direction, TVector2 const &centre, TVector2 const &start1, TVector2 const &end1, TVector2 const &start2, TVector2 const &end2) const
process_name gaushit a
int MergeClusters(std::vector< art::PtrVector< recob::Hit > > const &planeClusters, std::vector< art::PtrVector< recob::Hit > > &clusters) const
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
Signal from induction planes.
Definition: geo_types.h:145
bool PassCuts(double const &angle, double const &crossingDistance, double const &projectedWidth, double const &separation, double const &overlap, double const &longLength) const
auto end(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:585
return match has_match and(match.match_pdg==11 or match.match_pdg==-11)
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:493
TVector2 HitCoordinates(art::Ptr< recob::Hit > const &hit) const
double FindMinSeparation(art::PtrVector< recob::Hit > const &cluster1, art::PtrVector< recob::Hit > const &cluster2) const
double FindCrossingDistance(TVector2 const &direction1, TVector2 const &centre1, TVector2 const &direction2, TVector2 const &centre2) const
finds tracks best matching by angle
double GlobalWire(geo::WireID const &wireID) const
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:406
double FindProjectedWidth(TVector2 const &centre1, TVector2 const &start1, TVector2 const &end1, TVector2 const &centre2, TVector2 const &start2, TVector2 const &end2) const
physics pm2 e1
BEGIN_PROLOG could also be cout