28 m_twoDDistanceCut(10.f)
39 if (!pInputVertexList || pInputVertexList->empty())
42 std::cout <<
"VertexRefinementAlgorithm: unable to find current vertex list " << std::endl;
44 return STATUS_CODE_SUCCESS;
48 std::string temporaryListName;
51 ClusterList clusterListU, clusterListV, clusterListW;
54 this->
RefineVertices(pInputVertexList, clusterListU, clusterListV, clusterListW);
56 if (!pOutputVertexList->empty())
62 return STATUS_CODE_SUCCESS;
70 for (
const std::string &clusterListName : inputClusterListNames)
76 if (!pClusterList || pClusterList->empty())
79 std::cout <<
"VertexRefinementAlgorithm: unable to find cluster list " << clusterListName << std::endl;
90 outputClusterList.insert(outputClusterList.end(), pClusterList->begin(), pClusterList->end());
99 for (
const Vertex *
const pVertex : *pVertexList)
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;
150 const Vertex *pNewVertex(NULL);
151 PANDORA_THROW_RESULT_IF(STATUS_CODE_SUCCESS, !=, PandoraContentApi::Vertex::Create(*
this, parameters, pNewVertex));
162 for (
const Cluster *
const pCluster : clusterList)
179 directions.push_back(eigenVectors.at(0).GetUnitVector());
184 if (intercepts.size() > 1)
188 return originalVtxPos;
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();
214 G(3 * i + 1, i + 3) = -directions[i].GetY();
215 G(3 * i + 2, i + 3) = -directions[i].GetZ();
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;
247 return STATUS_CODE_SUCCESS;
Header file for the cluster helper class.
Header file for the geometry helper class.
Header file for the principal curve analysis helper class.
#define PANDORA_THROW_RESULT_IF_AND_IF(StatusCode1, StatusCode2, Operator, Command)
#define PANDORA_THROW_RESULT_IF(StatusCode1, Operator, Command)
#define PANDORA_RETURN_RESULT_IF_AND_IF(StatusCode1, StatusCode2, Operator, Command)
#define PANDORA_RETURN_RESULT_IF(StatusCode1, Operator, Command)
Header file for the vertex refinement algorithm class.
static pandora::StatusCode ReplaceCurrentList(const pandora::Algorithm &algorithm, const std::string &newListName)
Replace the current list with a pre-saved list; use this new list as a permanent replacement for the ...
static pandora::StatusCode CreateTemporaryListAndSetCurrent(const pandora::Algorithm &algorithm, const T *&pT, std::string &temporaryListName)
Create a temporary list and set it to be the current list, enabling object creation.
static pandora::StatusCode GetCurrentList(const pandora::Algorithm &algorithm, const T *&pT)
Get the current list.
static pandora::StatusCode GetList(const pandora::Algorithm &algorithm, const std::string &listName, const T *&pT)
Get a named list.
static const pandora::PandoraSettings * GetSettings(const pandora::Algorithm &algorithm)
Get the pandora settings instance.
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.
static pandora::HitType GetClusterHitType(const pandora::Cluster *const pCluster)
Get the hit type associated with a two dimensional cluster.
static void GetCoordinateVector(const pandora::Cluster *const pCluster, pandora::CartesianPointVector &coordinateVector)
Get vector of hit coordinates from an input cluster.
static float GetClosestDistance(const pandora::ClusterList &clusterList1, const pandora::ClusterList &clusterList2)
Get closest distance between clusters in a pair of cluster lists.
static pandora::CartesianVector ProjectPosition(const pandora::Pandora &pandora, const pandora::CartesianVector &position3D, const pandora::HitType view)
Project 3D position into a given 2D view.
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.
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.
std::vector< pandora::CartesianVector > EigenVectors
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...
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.
pandora::StringVector m_inputClusterListNames
The list of input cluster list names.
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.
std::string m_outputVertexListName
The refined vertex list to be outputted.
pandora::StatusCode Run()
Run the algorithm.
pandora::StatusCode ReadSettings(const pandora::TiXmlHandle xmlHandle)
Read the algorithm settings.
float m_distanceCut
The maximum distance a refined vertex can be from the original position to be kept.
unsigned int m_minimumHitsCut
The minimum size of a cluster to be used in refinement.
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.
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.
float m_chiSquaredCut
The maximum chi2 value a refined vertex can have to be kept.
VertexRefinementAlgorithm()
Default constructor.
std::string m_inputVertexListName
The initial vertex list.
bool ShouldDisplayAlgorithmInfo() const
Whether to display algorithm information during processing.
const Pandora & GetPandora() const
Get the associated pandora instance.
StatusCodeException class.
static StatusCode ReadVectorOfValues(const TiXmlHandle &xmlHandle, const std::string &xmlElementName, std::vector< T > &vector)
Read a vector of values from a (space separated) list in an xml element.
static StatusCode ReadValue(const TiXmlHandle &xmlHandle, const std::string &xmlElementName, T &t)
Read a value from an xml element.
HitType
Calorimeter hit type enum.
MANAGED_CONTAINER< const Cluster * > ClusterList
std::vector< std::string > StringVector
std::vector< CartesianVector > CartesianPointVector
std::vector< float > FloatVector
MANAGED_CONTAINER< const Vertex * > VertexList
StatusCode
The StatusCode enum.