opatIO-cpp 0.3.0a
Open Parametrized Array Table
|
Represents a lattice structure for interpolating data from an OPAT object. More...
#include <tableLattice.h>
Public Member Functions | |
TableLattice (const opat::OPAT &opat) | |
Constructs a TableLattice object from an OPAT object. | |
TableLattice (const opat::OPAT &opat, const InterpolationType &interpolationType) | |
Constructs a TableLattice object from an OPAT object with a specified interpolation type. | |
TableLattice (const TableLattice &)=default | |
TableLattice & | operator= (const TableLattice &)=delete |
TableLattice (TableLattice &&) noexcept=default | |
TableLattice & | operator= (TableLattice &&) noexcept=delete |
~TableLattice ()=default | |
DataCard | get (const FloatIndexVector &indexVector) const |
Retrieves interpolated data for a given index vector. | |
InterpolationType | getInterpolationType () const |
Gets the current interpolation type. | |
void | setInterpolationType (InterpolationType interpolationType) |
Sets the interpolation type. | |
void | dumpTriangulationToAscii (const std::string &points_file, const std::string &simplices_file) const |
Dumps the Delaunay triangulation to ASCII files. | |
Private Member Functions | |
void | initialize () |
Initializes the TableLattice internal structures. | |
void | buildDelaunay () |
Builds the Delaunay triangulation of the index vectors. | |
Simplex | findContainingSimplex (const FloatIndexVector &indexVector) const |
Finds the simplex containing the given index vector using a walk algorithm. | |
void | validateIndexVector (const FloatIndexVector &indexVector) const |
Validates if the given index vector is within the global bounds of the table and has the correct dimension. | |
std::vector< double > | calculateBarycentricWeights (const FloatIndexVector &queryPoint, const std::vector< FloatIndexVector > &simplexActualVertices) const |
Calculates the barycentric weights of a query point with respect to the vertices of a given simplex. | |
Private Attributes | |
const opat::OPAT & | m_opat |
Reference to the OPAT object. | |
std::size_t | m_indexVectorSize {} |
The dimensionality of the index vectors. | |
InterpolationType | m_interpolationType {InterpolationType::Linear} |
The type of interpolation to use. | |
std::vector< FloatIndexVector > | m_indexVectors |
Stores all unique index vectors from the OPAT file, serving as the vertices of the triangulation. | |
std::vector< std::vector< double > > | m_axisValues |
Stores the unique values for each axis/dimension (Not currently used by Delaunay approach). | |
std::size_t | m_numCorners {} |
The number of corners in a hypercube (2^m_indexVectorSize), relevant for hypercube-based approaches (not Delaunay). | |
std::vector< std::vector< std::size_t > > | m_simplices |
Stores the simplices of the Delaunay triangulation. Each inner vector is a list of global vertex indices (indices into m_indexVectors ). | |
std::vector< std::vector< std::size_t > > | m_simplexAdjacency |
Adjacency list for simplices. m_simplexAdjacency[i][j] stores the ID of the simplex adjacent to simplex i across the face opposite to its j -th local vertex. A value of static_cast<std::size_t>(-1) indicates no neighbor (boundary). | |
Simplex | m_lastFoundSimplex |
Stores the last found simplex (ID and barycentric weights) as a starting point for the findContainingSimplex walk algorithm, optimizing searches for spatially coherent query points. Initialized with an invalid ID. | |
Represents a lattice structure for interpolating data from an OPAT object.
The TableLattice class builds a Delaunay triangulation of the OPAT index points and provides methods to interpolate data for a given index vector.
Example:
|
explicit |
Constructs a TableLattice object from an OPAT object.
Initializes the lattice using the index vectors from the OPAT object and builds a Delaunay triangulation using Qhull. Defaults to Linear interpolation.
opat | The OPAT object containing the data. |
std::runtime_error | if Delaunay triangulation construction fails (e.g., due to Qhull errors, which could be caused by insufficient or degenerate input points from the OPAT file). Resolution: Ensure the OPAT file contains valid and sufficient index points for triangulation. |
Example:
opat::lattice::TableLattice::TableLattice | ( | const opat::OPAT & | opat, |
const InterpolationType & | interpolationType ) |
Constructs a TableLattice object from an OPAT object with a specified interpolation type.
Initializes the lattice using the index vectors from the OPAT object and builds a Delaunay triangulation using Qhull.
opat | The OPAT object containing the data. |
interpolationType | The type of interpolation to use. |
std::runtime_error | if the specified interpolationType is not InterpolationType::Linear , as other types are not yet implemented. |
std::runtime_error | if Delaunay triangulation construction fails (e.g., due to Qhull errors, which could be caused by insufficient or degenerate input points from the OPAT file). Resolution: Ensure the OPAT file contains valid and sufficient index points for triangulation. Use InterpolationType::Linear . |
Example:
|
default |
|
defaultnoexcept |
|
default |
|
private |
Builds the Delaunay triangulation of the index vectors.
Uses the Qhull library to perform Delaunay triangulation on the points stored in m_indexVectors
. Populates m_simplices
with the vertex indices of each simplex and m_simplexAdjacency
with neighbor information. This method is called by the constructors after initialize()
.
std::runtime_error | (wrapping orgQhull::QhullError ) if Qhull encounters an error during triangulation. This can happen due to issues like insufficient points for the dimension, all points being collinear/coplanar in a way that prevents triangulation, or other Qhull internal errors. Resolution: Check the validity and distribution of index vectors in the input OPAT file. Ensure there are enough non-degenerate points for the given number of dimensions. |
|
private |
Calculates the barycentric weights of a query point with respect to the vertices of a given simplex.
Solves a linear system of equations to find the weights. The sum of barycentric weights is 1. If a point is inside a simplex, all its weights are non-negative.
queryPoint | The point for which to calculate barycentric weights. |
simplexActualVertices | A vector of FloatIndexVector representing the coordinates of the simplex's vertices. The number of vertices must be m_indexVectorSize + 1 . |
double
containing the barycentric weights. The order of weights corresponds to the order of vertices in simplexActualVertices
. std::invalid_argument | if simplexActualVertices does not contain m_indexVectorSize + 1 vertices, or if queryPoint or any vertex in simplexActualVertices does not have m_indexVectorSize dimensions. Resolution: Ensure input vectors have correct sizes and dimensions matching the lattice. |
std::runtime_error | if solveLinearSystem fails, typically because the matrix formed by simplex vertices is singular (e.g., degenerate simplex - vertices are collinear/coplanar). Resolution: This indicates an issue with the geometry of the provided simplex vertices. |
std::logic_error | if solveLinearSystem returns an unexpected number of solution components. Resolution: This indicates an internal logic error in the linear solver or its usage. |
void opat::lattice::TableLattice::dumpTriangulationToAscii | ( | const std::string & | points_file, |
const std::string & | simplices_file ) const |
Dumps the Delaunay triangulation to ASCII files.
Creates two files: one for the points (vertices) and one for the simplices. This is useful for visualizing or debugging the triangulation.
points_file | The path to the file where points will be saved. |
simplices_file | The path to the file where simplices will be saved. |
std::ios_base::failure | if an error occurs during file I/O operations (e.g., cannot open file for writing). Resolution: Ensure the program has write permissions to the specified file paths and that the paths are valid. |
Example:
|
nodiscardprivate |
Finds the simplex containing the given index vector using a walk algorithm.
Starts from a known simplex (potentially the last found one, m_lastFoundSimplex
) and "walks" through adjacent simplices until the one containing the indexVector
is found.
indexVector | The index vector to locate. |
Simplex
struct containing the ID of the simplex and the barycentric weights of the indexVector
within it. std::out_of_range | if validateIndexVector fails (point outside overall bounds), or if the walk terminates at the convex hull and the point is determined to be outside. Also thrown for internal out-of-bounds access to m_simplices , m_indexVectors , or m_simplexAdjacency . Resolution: Ensure indexVector is within the data bounds and convex hull. Internal errors suggest issues with triangulation data. |
std::invalid_argument | if validateIndexVector fails (dimension mismatch) or if indexVector dimension doesn't match m_indexVectorSize during the walk. Resolution: Ensure indexVector has the correct number of dimensions. |
std::runtime_error | if the triangulation is empty (m_simplices.empty() ), if a cycle is detected during the walk, if the walk cannot determine an exit face from a simplex, if the maximum number of walk steps is exceeded, or if calculateBarycentricWeights fails (e.g. degenerate simplex). Resolution: These often indicate issues with the triangulation's integrity, adjacency information, numerical precision challenges, or a query point far outside the domain leading to extended walks. Verify the input OPAT data. |
std::logic_error | from calculateBarycentricWeights if an unexpected number of weights is returned. |
|
nodiscard |
Retrieves interpolated data for a given index vector.
Finds the containing simplex for the index vector and performs barycentric interpolation of the data from the simplex vertices.
indexVector | The index vector for which to interpolate data. |
std::out_of_range | if the indexVector is outside the bounds of the table data (as determined by validateIndexVector ), or if findContainingSimplex determines the point is outside the convex hull of the data points, or if an internal ID (vertex/simplex) is out of bounds. Resolution: Ensure indexVector values are within the min/max range for each dimension (check opat::OPAT::getBounds() ). Ensure the point lies within the convex hull of the OPAT index points. |
std::invalid_argument | if the indexVector has a dimension different from the OPAT data's index dimension (checked by validateIndexVector or findContainingSimplex ). Resolution: Ensure indexVector.size() matches opat_object.header.numIndex . |
std::runtime_error | if findContainingSimplex fails due to internal errors (e.g., no simplices, cycle in walk, cannot determine exit face, max walk steps reached), or if calculateBarycentricWeights fails due to a degenerate simplex (singular matrix). Resolution: These errors often indicate issues with the underlying triangulation or extreme numerical conditions. Verify the integrity of the input OPAT file's index data. |
std::logic_error | if calculateBarycentricWeights returns an unexpected number of weights (internal logic error). |
Example:
|
nodiscard |
Gets the current interpolation type.
|
private |
Initializes the TableLattice internal structures.
Extracts unique index vectors from the OPAT object and stores them in m_indexVectors
. Sets m_indexVectorSize
based on OPAT header. This method is called by the constructors.
|
delete |
|
deletenoexcept |
void opat::lattice::TableLattice::setInterpolationType | ( | InterpolationType | interpolationType | ) |
Sets the interpolation type.
interpolationType | The new InterpolationType to set. |
std::runtime_error | if interpolationType is not InterpolationType::Linear , as other types are not currently implemented. Resolution: Only use InterpolationType::Linear . |
Example:
|
private |
Validates if the given index vector is within the global bounds of the table and has the correct dimension.
indexVector | The index vector to validate. |
std::invalid_argument | if indexVector.size() does not match m_opat.header.numIndex . Resolution: Ensure the query indexVector has the same number of dimensions as the OPAT data. |
std::out_of_range | if any component of indexVector is outside the min/max bounds defined by m_opat.getBounds() for that dimension. Resolution: Ensure all components of the indexVector are within the established bounds of the OPAT data. |
Example:
|
private |
Stores the unique values for each axis/dimension (Not currently used by Delaunay approach).
|
private |
Stores all unique index vectors from the OPAT file, serving as the vertices of the triangulation.
|
private |
The dimensionality of the index vectors.
|
private |
The type of interpolation to use.
|
mutableprivate |
Stores the last found simplex (ID and barycentric weights) as a starting point for the findContainingSimplex
walk algorithm, optimizing searches for spatially coherent query points. Initialized with an invalid ID.
|
private |
The number of corners in a hypercube (2^m_indexVectorSize), relevant for hypercube-based approaches (not Delaunay).
|
private |
Reference to the OPAT object.
|
private |
Adjacency list for simplices. m_simplexAdjacency[i][j]
stores the ID of the simplex adjacent to simplex i
across the face opposite to its j
-th local vertex. A value of static_cast<std::size_t>(-1)
indicates no neighbor (boundary).
|
private |
Stores the simplices of the Delaunay triangulation. Each inner vector is a list of global vertex indices (indices into m_indexVectors
).