9 #include "art/Utilities/ToolMacros.h"
10 #include "art/Utilities/make_tool.h"
11 #include "cetlib/cpu_timer.h"
12 #include "fhiclcpp/ParameterSet.h"
13 #include "messagefacility/MessageLogger/MessageLogger.h"
34 namespace lar_cluster3d {
51 void configure(fhicl::ParameterSet
const &pset)
override;
83 reco::ClusterParametersList::iterator positionItr,
88 const Eigen::Vector3f&,
89 const Eigen::Vector3f&,
90 const Eigen::Vector3f&,
92 Eigen::Vector3f&)
const;
94 using Point = std::tuple<float,float,const reco::ClusterHit3D*>;
113 m_pcaAlg(pset.
get<fhicl::ParameterSet>(
"PrincipalComponentsAlg"))
130 m_clusterAlg = art::make_tool<lar_cluster3d::IClusterAlg>(pset.get<fhicl::ParameterSet>(
"ClusterAlg"));
152 cet::cpu_timer theClockBuildClusters;
157 int countClusters(0);
162 reco::ClusterParametersList::iterator clusterParametersListItr = clusterParametersList.begin();
164 while(clusterParametersListItr != clusterParametersList.end())
169 std::cout <<
"**> Looking at Cluster " << countClusters++ << std::endl;
185 std::cout <<
">>>>>>>>>>> Reclustered to " << reclusteredParameters.size() <<
" Clusters <<<<<<<<<<<<<<<" << std::endl;
188 if (!reclusteredParameters.empty())
191 for (
auto&
cluster : reclusteredParameters)
193 std::cout <<
"****> Calling breakIntoTinyBits" << std::endl;
198 std::cout <<
"****> Broke Cluster with " <<
cluster.getHitPairListPtr().size() <<
" into " <<
cluster.daughterList().size() <<
" sub clusters";
199 for(
auto& clus :
cluster.daughterList())
std::cout <<
", " << clus.getHitPairListPtr().size();
209 clusterParametersListItr++;
214 theClockBuildClusters.stop();
219 mf::LogDebug(
"Cluster3D") <<
">>>>> Cluster Path finding done" << std::endl;
225 reco::ClusterParametersList::iterator positionItr,
237 std::string pluses(level/2,
'+');
238 std::string indent(level/2,
' ');
242 reco::ClusterParametersList::iterator inputPositionItr = positionItr;
248 std::cout << indent <<
">>> breakIntoTinyBits with " << clusterToBreak.
getHitPairListPtr().size() <<
" input hits " << std::endl;
255 bool storeCurrentCluster(
true);
270 clusHitPairVector.sort([](
const auto&
left,
const auto&
right){
return left->getArclenToPoca() <
right->getArclenToPoca();});
274 std::vector<const reco::ClusterHit3D*> vertexHitVec;
276 std::cout << indent <<
"+> Breaking cluster, convex hull has " << bestEdgeList.size() <<
" edges to work with" << std::endl;
278 for(
const auto& edge : bestEdgeList)
280 vertexHitVec.push_back(std::get<0>(edge));
281 vertexHitVec.push_back(std::get<1>(edge));
285 std::sort(vertexHitVec.begin(),vertexHitVec.end(),[](
const auto&
left,
const auto&
right){
return left->getArclenToPoca() <
right->getArclenToPoca();});
288 using Hit3DItrPair = std::pair<reco::HitPairListPtr::iterator,reco::HitPairListPtr::iterator>;
289 using VertexPairList = std::list<Hit3DItrPair>;
291 VertexPairList vertexPairList;
292 reco::HitPairListPtr::iterator firstHitItr = clusHitPairVector.begin();
294 for(
const auto& hit3D : vertexHitVec)
296 reco::HitPairListPtr::iterator vertexItr = std::find(firstHitItr,clusHitPairVector.end(),hit3D);
298 if (vertexItr == clusHitPairVector.end())
300 std::cout << indent <<
">>>>>>>>>>>>>>>>> Hit not found in input list, cannot happen? <<<<<<<<<<<<<<<<<<<" << std::endl;
304 std::cout << indent <<
"+> -- Distance from first to current vertex point: " <<
std::distance(firstHitItr,vertexItr) <<
" first: " << *firstHitItr <<
", vertex: " << *vertexItr;
307 if (
std::distance(firstHitItr,vertexItr) > minimumClusterSize)
309 vertexPairList.emplace_back(Hit3DItrPair(firstHitItr,vertexItr));
310 firstHitItr = vertexItr;
321 std::cout << indent <<
"+> loop over vertices done, remant distance: " <<
std::distance(firstHitItr,clusHitPairVector.end()) << std::endl;
324 if (!vertexPairList.empty() &&
std::distance(firstHitItr,clusHitPairVector.end()) < minimumClusterSize)
325 vertexPairList.back().second = clusHitPairVector.end();
327 vertexPairList.emplace_back(Hit3DItrPair(firstHitItr,clusHitPairVector.end()));
330 std::cout << indent <<
"+> ---> breaking cluster into " << vertexPairList.size() <<
" subclusters" << std::endl;
332 if (vertexPairList.size() > 1)
334 storeCurrentCluster =
false;
337 for(
auto& hit3DItrPair : vertexPairList)
342 std::cout << indent <<
"+> -- building new cluster, size: " <<
std::distance(hit3DItrPair.first,hit3DItrPair.second) << std::endl;
345 hitPairListPtr.resize(
std::distance(hit3DItrPair.first,hit3DItrPair.second));
348 std::copy(hit3DItrPair.first,hit3DItrPair.second,hitPairListPtr.begin());
359 std::cout << indent <<
"+> -- >> cluster has a valid Full PCA" << std::endl;
366 if (fullPrimaryVec.dot(newPrimaryVec) < 0.)
368 for(
size_t vecIdx = 0; vecIdx < 3; vecIdx++) newFullPCA.
flipAxis(vecIdx);
374 positionItr =
breakIntoTinyBits(clusterParams, positionItr, outputClusterList, level+4);
381 if (storeCurrentCluster)
386 std::set<const reco::ClusterHit2D*> hitSet;
391 for(
const auto& hit2D : hit3D->getHits())
393 if (hit2D) hitSet.insert(hit2D);
398 for(
const auto& hit2D : hitSet)
404 std::cout << indent <<
"*********>>> storing new subcluster of size " << clusterToBreak.
getHitPairListPtr().size() << std::endl;
406 positionItr = outputClusterList.insert(positionItr,clusterToBreak);
411 else if (inputPositionItr == positionItr)
std::cout << indent <<
"***** DID NOT STORE A CLUSTER *****" << std::endl;
419 std::string minuses(level/2,
'-');
420 std::string indent(level/2,
' ');
432 using Point = std::tuple<float,float,const reco::ClusterHit3D*>;
439 Eigen::Vector3f pcaToHitVec(hit3D->getPosition()[0] - pcaCenter(0),
440 hit3D->getPosition()[1] - pcaCenter(1),
441 hit3D->getPosition()[2] - pcaCenter(2));
444 pointList.emplace_back(
dcel2d::Point(pcaToHit(1),pcaToHit(2),hit3D));
448 pointList.sort([](
const auto&
left,
const auto&
right){
return (
std::abs(std::get<0>(
left) - std::get<0>(
right)) > std::numeric_limits<float>::epsilon()) ? std::get<0>(
left) < std::get<0>(
right) : std::get<1>(
left) < std::get<1>(
right);});
451 std::vector<ConvexHull> convexHullVec;
452 std::vector<PointList> rejectedListVec;
453 bool increaseDepth(pointList.size() > 5);
454 float lastArea(std::numeric_limits<float>::max());
459 convexHullVec.push_back(
ConvexHull(pointList));
462 const ConvexHull& convexHull = convexHullVec.back();
463 PointList& rejectedList = rejectedListVec.back();
466 increaseDepth =
false;
470 std::cout << indent <<
"-> built convex hull, 3D hits: " << pointList.size() <<
" with " << convexHullPoints.size() <<
" vertices" <<
", area: " << convexHull.
getConvexHullArea() << std::endl;
472 for(
const auto& point : convexHullPoints)
473 std::cout <<
" (" << std::get<0>(point) <<
"," << std::get<1>(point) <<
")";
478 for(
auto& point : convexHullPoints)
480 pointList.remove(point);
481 rejectedList.emplace_back(point);
490 while(!convexHullVec.empty() && convexHullVec.back().getConvexHullArea() < 0.5)
492 convexHullVec.pop_back();
493 rejectedListVec.pop_back();
497 if (!convexHullVec.empty())
499 size_t nRejectedTotal(0);
502 for(
const auto& rejectedList : rejectedListVec)
504 nRejectedTotal += rejectedList.size();
506 for(
const auto& rejectedPoint : rejectedList)
508 std::cout << indent <<
"-> -- Point is " << convexHullVec.back().findNearestDistance(rejectedPoint) <<
" from nearest edge" << std::endl;
510 if (convexHullVec.back().findNearestDistance(rejectedPoint) > 0.5)
511 hitPairListPtr.remove(std::get<2>(rejectedPoint));
515 std::cout << indent <<
"-> Removed " << nRejectedTotal <<
" leaving " << pointList.size() <<
"/" << hitPairListPtr.size() <<
" points" << std::endl;
521 Point lastPoint = convexHullVec.back().getConvexHull().front();
523 for(
auto& curPoint : convexHullVec.back().getConvexHull())
525 if (curPoint == lastPoint)
continue;
530 float distBetweenPoints = (curPoint3D->getPosition()[0] - lastPoint3D->
getPosition()[0]) * (curPoint3D->getPosition()[0] - lastPoint3D->
getPosition()[0])
531 + (curPoint3D->getPosition()[1] - lastPoint3D->
getPosition()[1]) * (curPoint3D->getPosition()[1] - lastPoint3D->
getPosition()[1])
532 + (curPoint3D->getPosition()[2] - lastPoint3D->
getPosition()[2]) * (curPoint3D->getPosition()[2] - lastPoint3D->
getPosition()[2]);
534 distBetweenPoints = std::sqrt(distBetweenPoints);
538 edgeMap[lastPoint3D].push_back(edge);
539 edgeMap[curPoint3D].push_back(edge);
540 bestEdgeList.emplace_back(edge);
542 lastPoint = curPoint;
563 Eigen::Vector3f pcaToHitVec(hit3D->getPosition()[0] - pcaCenter(0),
564 hit3D->getPosition()[1] - pcaCenter(1),
565 hit3D->getPosition()[2] - pcaCenter(2));
568 pointList.emplace_back(
dcel2d::Point(pcaToHit(1),pcaToHit(2),hit3D));
572 pointList.sort([](
const auto&
left,
const auto&
right){
return (
std::abs(std::get<0>(
left) - std::get<0>(
right)) > std::numeric_limits<float>::epsilon()) ? std::get<0>(
left) < std::get<0>(
right) : std::get<1>(
left) < std::get<1>(
right);});
588 for(
auto&
vertex : vertexList)
590 Eigen::Vector3f coords = rotationMatrixInv *
vertex.getCoords();
614 dcel2d::Point lastPoint = voronoiDiagram.getConvexHull().front();
619 for(
auto& curPoint : voronoiDiagram.getConvexHull())
621 if (curPoint == lastPoint)
continue;
626 float distBetweenPoints = (curPoint3D->getPosition()[0] - lastPoint3D->
getPosition()[0]) * (curPoint3D->getPosition()[0] - lastPoint3D->
getPosition()[0])
627 + (curPoint3D->getPosition()[1] - lastPoint3D->
getPosition()[1]) * (curPoint3D->getPosition()[1] - lastPoint3D->
getPosition()[1])
628 + (curPoint3D->getPosition()[2] - lastPoint3D->
getPosition()[2]) * (curPoint3D->getPosition()[2] - lastPoint3D->
getPosition()[2]);
630 distBetweenPoints = std::sqrt(distBetweenPoints);
634 edgeMap[lastPoint3D].push_back(edge);
635 edgeMap[curPoint3D].push_back(edge);
636 bestEdgeList.emplace_back(edge);
638 lastPoint = curPoint;
641 std::cout <<
"****> vertexList containted " << vertexList.size() <<
" vertices" << std::endl;
648 const Eigen::Vector3f& u0,
649 const Eigen::Vector3f& P1,
650 const Eigen::Vector3f& u1,
651 Eigen::Vector3f& poca0,
652 Eigen::Vector3f& poca1)
const
655 Eigen::Vector3f w0 = P0 - P1;
661 float den(a * c - b * b);
663 float arcLen0 = (b *
e - c * d) / den;
664 float arcLen1 = (a *
e - b * d) / den;
666 poca0 = P0 + arcLen0 * u0;
667 poca1 = P1 + arcLen1 * u1;
669 return (poca0 - poca1).norm();
void PCAAnalysis_calc3DDocas(const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca) const
void PCAAnalysis_3D(const reco::HitPairListPtr &hitPairList, reco::PrincipalComponents &pca, bool skeletonOnly=false) const
void flipAxis(size_t axis)
const Eigen::Vector3f getPosition() const
void configure(fhicl::ParameterSet const &pset) override
Interface for configuring the particular algorithm tool.
reco::PrincipalComponents & getSkeletonPCA()
~ClusterPathFinder()
Destructor.
std::unique_ptr< lar_cluster3d::IClusterAlg > m_clusterAlg
Algorithm to do 3D space point clustering.
reco::EdgeList & getBestEdgeList()
float closestApproach(const Eigen::Vector3f &, const Eigen::Vector3f &, const Eigen::Vector3f &, const Eigen::Vector3f &, Eigen::Vector3f &, Eigen::Vector3f &) const
reco::Hit3DToEdgeMap & getHit3DToEdgeMap()
reco::HitPairListPtr & getHitPairListPtr()
std::tuple< float, float, const reco::ClusterHit3D * > Point
std::pair< MinMaxPoints, MinMaxPoints > MinMaxPointPair
IClusterModAlg interface class definiton.
ClusterParametersList & daughterList()
Implements a ConvexHull for use in clustering.
std::list< EdgeTuple > EdgeList
void buildConvexHull(reco::ClusterParameters &clusterParameters, int level=0) const
std::list< Point > PointList
void initializeHistograms(art::TFileDirectory &) override
Interface for initializing histograms if they are desired Note that the idea is to put hisgtograms in...
VoronoiDiagram class definiton.
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
reco::PrincipalComponents & getFullPCA()
float getTimeToExecute() const override
If monitoring, recover the time to execute a particular function.
std::tuple< double, double, const reco::ClusterHit3D * > Point
Definitions used by the VoronoiDiagram algorithm.
void buildVoronoiDiagram(const dcel2d::PointList &)
Given an input set of 2D points construct a 2D voronoi diagram.
const Eigen::Vector3f & getAvePosition() const
size_t m_minTinyClusterSize
Minimum size for a "tiny" cluster.
std::tuple< const reco::ClusterHit3D *, const reco::ClusterHit3D *, double > EdgeTuple
std::list< const reco::ClusterHit3D * > HitPairListPtr
void UpdateParameters(const reco::ClusterHit2D *hit)
float getConvexHullArea() const
recover the area of the convex hull
This provides an art tool interface definition for 3D Cluster algorithms.
const PointList & getConvexHull() const
recover the list of convex hull vertices
reco::ClusterParametersList::iterator breakIntoTinyBits(reco::ClusterParameters &cluster, reco::ClusterParametersList::iterator positionItr, reco::ClusterParametersList &outputClusterList, int level=0) const
Use PCA to try to find path in cluster.
std::pair< Point, Point > MinMaxPoints
ConvexHull class definiton.
This header file defines the interface to a principal components analysis designed to be used within ...
bool m_enableMonitoring
Data members to follow.
void buildVoronoiDiagram(reco::ClusterParameters &clusterParameters, int level=0) const
PrincipalComponentsAlg m_pcaAlg
dcel2d::FaceList & getFaceList()
This provides an art tool interface definition for 3D Cluster algorithms.
void ModifyClusters(reco::ClusterParametersList &) const override
Scan an input collection of clusters and modify those according to the specific implementing algorith...
std::list< Point > PointList
std::unordered_map< const reco::ClusterHit3D *, reco::EdgeList > Hit3DToEdgeMap
process_name physics producers generator hPHist_pi physics producers generator P0
ClusterPathFinder(const fhicl::ParameterSet &)
Constructor.
dcel2d::HalfEdgeList & getHalfEdgeList()
std::list< Vertex > VertexList
BEGIN_PROLOG don t mess with this pandoraTrackGausCryoW true
BEGIN_PROLOG could also be cout
std::list< ClusterParameters > ClusterParametersList
const EigenVectors & getEigenVectors() const
dcel2d::VertexList & getVertexList()