opatIO-cpp 0.3.0a
Open Parametrized Array Table
Loading...
Searching...
No Matches
opat::lattice::TableLattice Class Reference

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
 
TableLatticeoperator= (const TableLattice &)=delete
 
 TableLattice (TableLattice &&) noexcept=default
 
TableLatticeoperator= (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::OPATm_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< FloatIndexVectorm_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.
 

Detailed Description

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:

#include "opatIO.h" // For opat::readOPAT, opat::OPAT, opat::DataCard
#include "indexVector.h" // For FloatIndexVector
#include "tableLattice.h" // For opat::lattice::TableLattice
#include <iostream>
#include <string>
#include <vector>
int main() {
// Define the path to your OPAT file
std::string opat_filename = "gs98Hz.opat"; // Replace with your actual file path
// 1. Load the OPAT file
opat::OPAT opat_data = opat::readOPAT(opat_filename);
std::cout << "Successfully loaded OPAT file: " << opat_filename << std::endl;
// 2. Create a TableLattice object from the loaded OPAT data
// This will automatically build the Delaunay triangulation.
// Default interpolation is Linear. (Also the only one currently implemented)
std::cout << "TableLattice initialized." << std::endl;
// Optionally, set a different interpolation type (not yet implemented)
// lattice.setInterpolationType(opat::lattice::InterpolationType::Quadratic);
// 3. Define an index vector for which you want to interpolate data.
// The dimension of this vector must match opat_data.header.numIndex.
// For this example, let's assume a 2D OPAT file, and we want data for (0.54421, 0.077585).
FloatIndexVector query_point({0.54421, 0.077585});
std::cout << "Querying for index vector: " << query_point << std::endl;
// 4. Get the interpolated DataCard for the query_point.
opat::DataCard interpolated_card = lattice.get(query_point);
// 5. Access a specific table from the interpolated DataCard.
std::string target_table_tag = "data";
const opat::OPATTable& interpolated_table = interpolated_card[target_table_tag];
std::cout << "Successfully accessed interpolated table: " << target_table_tag << std::endl;
// You can now use the interpolated_table, e.g., print it or extract values
// interpolated_table.print();
// Example of dumping the triangulation (useful for debugging)
// lattice.dumpTriangulationToAscii("points_dump.txt", "simplices_dump.txt");
// std::cout << "Triangulation dumped to files." << std::endl;
return 0;
}
Definition indexVector.h:38
Represents a lattice structure for interpolating data from an OPAT object.
Definition tableLattice.h:110
Header file defining the FloatIndexVector class for handling floating-point index vectors.
Namespace for table lattice interpolation of OPAT files.
Definition tableLattice.cpp:28
OPAT readOPAT(const std::string &filename)
Reads an OPAT file and returns its contents as an OPAT structure.
Definition opatIO.cpp:76
Header file for the OPAT I/O library, providing structures and functions for reading and manipulating...
Structure to hold a DataCard, which contains multiple tables.
Definition opatIO.h:409
Structure to hold the entire OPAT file.
Definition opatIO.h:502
Structure to hold the data of an OPAT table.
Definition opatIO.h:259

Constructor & Destructor Documentation

◆ TableLattice() [1/4]

opat::lattice::TableLattice::TableLattice ( const opat::OPAT & opat)
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.

Parameters
opatThe OPAT object containing the data.
Exceptions
std::runtime_errorif 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::OPAT my_opat_data = opat::readOPAT("my_data.opat");
// Use the lattice...

◆ TableLattice() [2/4]

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.

Parameters
opatThe OPAT object containing the data.
interpolationTypeThe type of interpolation to use.
Exceptions
std::runtime_errorif the specified interpolationType is not InterpolationType::Linear, as other types are not yet implemented.
std::runtime_errorif 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:

opat::OPAT my_opat_data = opat::readOPAT("my_data.opat");
// Use the lattice
@ Linear
Linear interpolation.
Definition tableLattice.h:47

◆ TableLattice() [3/4]

opat::lattice::TableLattice::TableLattice ( const TableLattice & )
default

◆ TableLattice() [4/4]

opat::lattice::TableLattice::TableLattice ( TableLattice && )
defaultnoexcept

◆ ~TableLattice()

opat::lattice::TableLattice::~TableLattice ( )
default

Member Function Documentation

◆ buildDelaunay()

void opat::lattice::TableLattice::buildDelaunay ( )
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().

Exceptions
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.

◆ calculateBarycentricWeights()

std::vector< double > opat::lattice::TableLattice::calculateBarycentricWeights ( const FloatIndexVector & queryPoint,
const std::vector< FloatIndexVector > & simplexActualVertices ) const
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.

Parameters
queryPointThe point for which to calculate barycentric weights.
simplexActualVerticesA vector of FloatIndexVector representing the coordinates of the simplex's vertices. The number of vertices must be m_indexVectorSize + 1.
Returns
A vector of double containing the barycentric weights. The order of weights corresponds to the order of vertices in simplexActualVertices.
Exceptions
std::invalid_argumentif 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_errorif 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_errorif solveLinearSystem returns an unexpected number of solution components. Resolution: This indicates an internal logic error in the linear solver or its usage.

◆ dumpTriangulationToAscii()

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.

Parameters
points_fileThe path to the file where points will be saved.
simplices_fileThe path to the file where simplices will be saved.
Exceptions
std::ios_base::failureif 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:

// Assuming 'lattice' is an initialized TableLattice object
lattice.dumpTriangulationToAscii("my_points.txt", "my_simplices.txt");
// This will create two files:
// "my_points.txt" containing the coordinates of the triangulation vertices.
// "my_simplices.txt" containing the vertex indices for each simplex.

◆ findContainingSimplex()

Simplex opat::lattice::TableLattice::findContainingSimplex ( const FloatIndexVector & indexVector) const
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.

Parameters
indexVectorThe index vector to locate.
Returns
A Simplex struct containing the ID of the simplex and the barycentric weights of the indexVector within it.
Exceptions
std::out_of_rangeif 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_argumentif 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_errorif 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_errorfrom calculateBarycentricWeights if an unexpected number of weights is returned.

◆ get()

DataCard opat::lattice::TableLattice::get ( const FloatIndexVector & indexVector) const
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.

Parameters
indexVectorThe index vector for which to interpolate data.
Returns
A DataCard containing the interpolated data.
Exceptions
std::out_of_rangeif 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_argumentif 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_errorif 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_errorif calculateBarycentricWeights returns an unexpected number of weights (internal logic error).

Example:

// Assuming 'lattice' is an initialized TableLattice object
// and 'my_opat_object' is the OPAT object it was built from.
// For an OPAT with 2 index dimensions:
FloatIndexVector target_iv({0.25, 0.75});
opat::DataCard result = lattice.get(target_iv);
// Use 'result' which contains interpolated data

◆ getInterpolationType()

InterpolationType opat::lattice::TableLattice::getInterpolationType ( ) const
nodiscard

Gets the current interpolation type.

Returns
The current InterpolationType.

◆ initialize()

void opat::lattice::TableLattice::initialize ( )
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.

◆ operator=() [1/2]

TableLattice & opat::lattice::TableLattice::operator= ( const TableLattice & )
delete

◆ operator=() [2/2]

TableLattice & opat::lattice::TableLattice::operator= ( TableLattice && )
deletenoexcept

◆ setInterpolationType()

void opat::lattice::TableLattice::setInterpolationType ( InterpolationType interpolationType)

Sets the interpolation type.

Parameters
interpolationTypeThe new InterpolationType to set.
Exceptions
std::runtime_errorif interpolationType is not InterpolationType::Linear, as other types are not currently implemented. Resolution: Only use InterpolationType::Linear.

Example:

// Assuming 'lattice' is an initialized TableLattice object

◆ validateIndexVector()

void opat::lattice::TableLattice::validateIndexVector ( const FloatIndexVector & indexVector) const
private

Validates if the given index vector is within the global bounds of the table and has the correct dimension.

Parameters
indexVectorThe index vector to validate.
Exceptions
std::invalid_argumentif 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_rangeif 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:

// Assuming 'lattice' is an initialized TableLattice object
// and 'opat_data' is the OPAT object it was built from.
FloatIndexVector valid_point({0.5, 0.5}); // Assuming 2D and bounds contain this
lattice.validateIndexVector(valid_point); // Should pass
// lattice.validateIndexVector(FloatIndexVector({0.5})); // Would throw std::invalid_argument if numIndex != 1
// lattice.validateIndexVector(FloatIndexVector({100.0, 100.0})); // Would throw std::out_of_range if bounds are smaller

Member Data Documentation

◆ m_axisValues

std::vector<std::vector<double> > opat::lattice::TableLattice::m_axisValues
private

Stores the unique values for each axis/dimension (Not currently used by Delaunay approach).

◆ m_indexVectors

std::vector<FloatIndexVector> opat::lattice::TableLattice::m_indexVectors
private

Stores all unique index vectors from the OPAT file, serving as the vertices of the triangulation.

◆ m_indexVectorSize

std::size_t opat::lattice::TableLattice::m_indexVectorSize {}
private

The dimensionality of the index vectors.

◆ m_interpolationType

InterpolationType opat::lattice::TableLattice::m_interpolationType {InterpolationType::Linear}
private

The type of interpolation to use.

◆ m_lastFoundSimplex

Simplex opat::lattice::TableLattice::m_lastFoundSimplex
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.

◆ m_numCorners

std::size_t opat::lattice::TableLattice::m_numCorners {}
private

The number of corners in a hypercube (2^m_indexVectorSize), relevant for hypercube-based approaches (not Delaunay).

◆ m_opat

const opat::OPAT& opat::lattice::TableLattice::m_opat
private

Reference to the OPAT object.

◆ m_simplexAdjacency

std::vector<std::vector<std::size_t> > opat::lattice::TableLattice::m_simplexAdjacency
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).

◆ m_simplices

std::vector<std::vector<std::size_t> > opat::lattice::TableLattice::m_simplices
private

Stores the simplices of the Delaunay triangulation. Each inner vector is a list of global vertex indices (indices into m_indexVectors).


The documentation for this class was generated from the following files: