9 #include "Pandora/AlgorithmHeaders.h" 17 #include <Eigen/Dense> 24 VertexRefinementAlgorithm::VertexRefinementAlgorithm() :
28 m_twoDDistanceCut(10.
f)
37 PANDORA_RETURN_RESULT_IF(STATUS_CODE_SUCCESS, !=, PandoraContentApi::GetCurrentList(*
this, pInputVertexList));
39 if (!pInputVertexList || pInputVertexList->empty())
41 if (PandoraContentApi::GetSettings(*this)->ShouldDisplayAlgorithmInfo())
42 std::cout <<
"VertexRefinementAlgorithm: unable to find current vertex list " << std::endl;
44 return STATUS_CODE_SUCCESS;
48 std::string temporaryListName;
49 PANDORA_RETURN_RESULT_IF(STATUS_CODE_SUCCESS, !=, PandoraContentApi::CreateTemporaryListAndSetCurrent(*
this, pOutputVertexList, temporaryListName));
51 ClusterList clusterListU, clusterListV, clusterListW;
54 this->
RefineVertices(pInputVertexList, clusterListU, clusterListV, clusterListW);
56 if (!pOutputVertexList->empty())
58 PANDORA_RETURN_RESULT_IF(STATUS_CODE_SUCCESS, !=, PandoraContentApi::SaveList<Vertex>(*
this,
m_outputVertexListName));
59 PANDORA_RETURN_RESULT_IF(STATUS_CODE_SUCCESS, !=, PandoraContentApi::ReplaceCurrentList<Vertex>(*
this,
m_outputVertexListName));
62 return STATUS_CODE_SUCCESS;
68 const StringVector &inputClusterListNames, ClusterList &clusterListU, ClusterList &clusterListV, ClusterList &clusterListW)
const 70 for (
const std::string &clusterListName : inputClusterListNames)
72 const ClusterList *pClusterList(
nullptr);
73 PANDORA_THROW_RESULT_IF_AND_IF(
74 STATUS_CODE_SUCCESS, STATUS_CODE_NOT_INITIALIZED, !=, PandoraContentApi::GetList(*
this, clusterListName, pClusterList));
76 if (!pClusterList || pClusterList->empty())
78 if (PandoraContentApi::GetSettings(*this)->ShouldDisplayAlgorithmInfo())
79 std::cout <<
"VertexRefinementAlgorithm: unable to find cluster list " << clusterListName << std::endl;
86 if ((TPC_VIEW_U != hitType) && (TPC_VIEW_V != hitType) && (TPC_VIEW_W != hitType))
87 throw StatusCodeException(STATUS_CODE_INVALID_PARAMETER);
89 ClusterList &outputClusterList((TPC_VIEW_U == hitType) ? clusterListU : (TPC_VIEW_V == hitType) ? clusterListV : clusterListW);
90 outputClusterList.insert(outputClusterList.end(), pClusterList->begin(), pClusterList->end());
97 const ClusterList &clusterListV,
const ClusterList &clusterListW)
const 99 for (
const Vertex *
const pVertex : *pVertexList)
101 const CartesianVector originalPosition(pVertex->GetPosition());
103 const CartesianVector vtxU(
105 const CartesianVector vtxV(
107 const CartesianVector vtxW(
110 CartesianVector vtxUV(0.
f, 0.
f, 0.
f), vtxUW(0.
f, 0.
f, 0.
f), vtxVW(0.
f, 0.
f, 0.
f), vtx3D(0.
f, 0.
f, 0.
f), position3D(0.
f, 0.
f, 0.
f);
111 float chi2UV(0.
f), chi2UW(0.
f), chi2VW(0.
f), chi23D(0.
f), chi2(0.
f);
118 if (chi2UV < chi2UW && chi2UV < chi2VW && chi2UV < chi23D)
123 else if (chi2UW < chi2VW && chi2UW < chi23D)
128 else if (chi2VW < chi23D)
140 position3D = originalPosition;
142 if ((position3D - originalPosition).GetMagnitude() >
m_distanceCut)
143 position3D = originalPosition;
145 PandoraContentApi::Vertex::Parameters parameters;
146 parameters.m_position = position3D;
147 parameters.m_vertexLabel = VERTEX_INTERACTION;
148 parameters.m_vertexType = VERTEX_3D;
150 const Vertex *pNewVertex(NULL);
151 PANDORA_THROW_RESULT_IF(STATUS_CODE_SUCCESS, !=, PandoraContentApi::Vertex::Create(*
this, parameters, pNewVertex));
159 CartesianPointVector intercepts, directions;
162 for (
const Cluster *
const pCluster : clusterList)
170 CartesianVector centroid(0.
f, 0.
f, 0.
f);
173 CartesianPointVector pointVector;
179 directions.push_back(eigenVectors.at(0).GetUnitVector());
183 CartesianVector newVtxPos(originalVtxPos);
184 if (intercepts.size() > 1)
188 return originalVtxPos;
196 const FloatVector &weights, CartesianVector &bestFitPoint)
const 198 const int n(intercepts.size());
200 Eigen::VectorXd
d = Eigen::VectorXd::Zero(3 *
n);
201 Eigen::MatrixXd G = Eigen::MatrixXd::Zero(3 *
n,
n + 3);
203 for (
int i = 0; i <
n; ++i)
205 d(3 * i) = intercepts[i].GetX() * weights[i];
206 d(3 * i + 1) = intercepts[i].GetY() * weights[i];
207 d(3 * i + 2) = intercepts[i].GetZ() * weights[i];
209 G(3 * i, 0) = weights[i];
210 G(3 * i + 1, 1) = weights[i];
211 G(3 * i + 2, 2) = weights[i];
213 G(3 * i, i + 3) = -directions[i].GetX() * weights[i];
214 G(3 * i + 1, i + 3) = -directions[i].GetY() * weights[i];
215 G(3 * i + 2, i + 3) = -directions[i].GetZ() * weights[i];
218 if ((G.transpose() * G).determinant() < std::numeric_limits<float>::epsilon())
220 bestFitPoint = CartesianVector(std::numeric_limits<float>::max(), std::numeric_limits<float>::max(), std::numeric_limits<float>::max());
224 Eigen::VectorXd m = (G.transpose() * G).inverse() * G.transpose() *
d;
226 bestFitPoint = CartesianVector(m[0], m[1], m[2]);
233 PANDORA_RETURN_RESULT_IF(STATUS_CODE_SUCCESS, !=, XmlHelper::ReadVectorOfValues(xmlHandle,
"InputClusterListNames",
m_inputClusterListNames));
235 PANDORA_RETURN_RESULT_IF(STATUS_CODE_SUCCESS, !=, XmlHelper::ReadValue(xmlHandle,
"InputVertexListName",
m_inputVertexListName));
237 PANDORA_RETURN_RESULT_IF(STATUS_CODE_SUCCESS, !=, XmlHelper::ReadValue(xmlHandle,
"OutputVertexListName",
m_outputVertexListName));
239 PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
"ChiSquaredCut",
m_chiSquaredCut));
241 PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
"DistanceCut",
m_distanceCut));
243 PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
"MinimumHitsCut",
m_minimumHitsCut));
245 PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
"TwoDDistanceCut",
m_twoDDistanceCut));
247 return STATUS_CODE_SUCCESS;
float m_chiSquaredCut
The maximum chi2 value a refined vertex can have to be kept.
static void MergeTwoPositions3D(const pandora::Pandora &pandora, const pandora::HitType view1, const pandora::HitType view2, const pandora::CartesianVector &position1, const pandora::CartesianVector &position2, pandora::CartesianVector &position3D, float &chiSquared)
Merge 2D positions from two views to give unified 3D position.
pandora::CartesianVector EigenValues
static pandora::CartesianVector ProjectPosition(const pandora::Pandora &pandora, const pandora::CartesianVector &position3D, const pandora::HitType view)
Project 3D position into a given 2D view.
unsigned int m_minimumHitsCut
The minimum size of a cluster to be used in refinement.
Header file for the principal curve analysis helper class.
pandora::CartesianVector RefineVertexTwoD(const pandora::ClusterList &clusterList, const pandora::CartesianVector &originalVtxPos) const
Refine the position of a two dimensional projection of a vertex using the clusters in that view...
std::string m_outputVertexListName
The refined vertex list to be outputted.
static pandora::HitType GetClusterHitType(const pandora::Cluster *const pCluster)
Get the hit type associated with a two dimensional cluster.
float m_twoDDistanceCut
The maximum distance a cluster can be from the original position to be used in refinement.
void GetClusterLists(const pandora::StringVector &inputClusterListNames, pandora::ClusterList &clusterListU, pandora::ClusterList &clusterListV, pandora::ClusterList &clusterListW) const
Get the input cluster lists.
void GetBestFitPoint(const pandora::CartesianPointVector &intercepts, const pandora::CartesianPointVector &directions, const pandora::FloatVector &weights, pandora::CartesianVector &bestFitPoint) const
Calculate the best fit point of a set of lines using a matrix equation.
Header file for the geometry helper class.
Header file for the cluster helper class.
static void MergeThreePositions3D(const pandora::Pandora &pandora, const pandora::HitType view1, const pandora::HitType view2, const pandora::HitType view3, const pandora::CartesianVector &position1, const pandora::CartesianVector &position2, const pandora::CartesianVector &position3, pandora::CartesianVector &position3D, float &chiSquared)
Merge 2D positions from three views to give unified 3D position.
float m_distanceCut
The maximum distance a refined vertex can be from the original position to be kept.
void RefineVertices(const pandora::VertexList *const pVertexList, const pandora::ClusterList &clusterListU, const pandora::ClusterList &clusterListV, const pandora::ClusterList &clusterListW) const
Perform the refinement proceduce on a list of vertices.
std::string m_inputVertexListName
The initial vertex list.
pandora::StringVector m_inputClusterListNames
The list of input cluster list names.
static void RunPca(const T &t, pandora::CartesianVector ¢roid, EigenValues &outputEigenValues, EigenVectors &outputEigenVectors)
Run principal component analysis using input calo hits (TPC_VIEW_U,V,W or TPC_3D; all treated as 3D p...
std::vector< pandora::CartesianVector > EigenVectors
pandora::StatusCode ReadSettings(const pandora::TiXmlHandle xmlHandle)
static void GetCoordinateVector(const pandora::Cluster *const pCluster, pandora::CartesianPointVector &coordinateVector)
Get vector of hit coordinates from an input cluster.
boost::graph_traits< ModuleGraph >::vertex_descriptor Vertex
static pandora::CartesianVector GetClosestPosition(const pandora::CartesianVector &position, const pandora::ClusterList &clusterList)
Get closest position in a list of clusters to a specified input position vector.
Header file for the vertex refinement algorithm class.
std::list< Vertex > VertexList
pandora::StatusCode Run()
static float GetClosestDistance(const pandora::ClusterList &clusterList1, const pandora::ClusterList &clusterList2)
Get closest distance between clusters in a pair of cluster lists.