numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
num_collect::rbf Namespace Reference

Namespace of RBF interpolation. More...

Namespaces

namespace  concepts
 Namespace of C++ concepts.
 
namespace  distance_functions
 Namespace of distance functions for RBF interpolation.
 
namespace  impl
 Namespace of internal implementations.
 
namespace  length_parameter_calculators
 Namespace of calculators of length parameters.
 
namespace  rbfs
 Namespace of RBFs.
 

Classes

class  gaussian_process_interpolator
 Class to interpolate using Gaussian process. More...
 
class  polynomial_calculator
 Class to calculate polynomials used with RBF interpolation. More...
 
class  polynomial_calculator< Variable, PolynomialDegree >
 Class to calculate polynomials used with RBF interpolation. More...
 
class  rbf_interpolator
 Class to interpolate using RBF. More...
 
class  rbf_interpolator< FunctionValue(Variable), RBF, KernelMatrixType, DistanceFunction, LengthParameterCalculator >
 Class to interpolate using RBF. More...
 
class  rbf_polynomial_interpolator
 Class to interpolate using RBF and polynomials. More...
 
class  rbf_polynomial_interpolator< FunctionValue(Variable), RBF, PolynomialDegree, KernelMatrixType, DistanceFunction, LengthParameterCalculator >
 Class to interpolate using RBF and polynomials. More...
 

Typedefs

template<typename FunctionSignature , concepts::rbf RBF = rbfs::gaussian_rbf<impl::get_default_scalar_type<FunctionSignature>>, kernel_matrix_type KernelMatrixType = kernel_matrix_type::dense, concepts::distance_function DistanceFunction = distance_functions::euclidean_distance_function< impl::get_variable_type_t<FunctionSignature>>>
using global_rbf_interpolator
 Class to interpolate using RBF with globally fixed length parameters.
 
template<typename FunctionSignature , concepts::rbf RBF = rbfs::gaussian_rbf<impl::get_default_scalar_type<FunctionSignature>>, index_type PolynomialDegree = 1, kernel_matrix_type KernelMatrixType = kernel_matrix_type::dense, concepts::distance_function DistanceFunction = distance_functions::euclidean_distance_function< impl::get_variable_type_t<FunctionSignature>>>
using global_rbf_polynomial_interpolator
 Class to interpolate using RBF and polynomials with globally fixed length parameters.
 
template<typename FunctionSignature , concepts::rbf RBF = rbfs::wendland_csrbf< impl::get_default_scalar_type<FunctionSignature>, 3, 1>, concepts::distance_function DistanceFunction = distance_functions::euclidean_distance_function< impl::get_variable_type_t<FunctionSignature>>>
using local_csrbf_interpolator
 Class to interpolate using compactly supported RBF with length parameters localized for sample points.
 
template<typename FunctionSignature , concepts::rbf RBF = rbfs::gaussian_rbf<impl::get_default_scalar_type<FunctionSignature>>, kernel_matrix_type KernelMatrixType = kernel_matrix_type::dense, concepts::distance_function DistanceFunction = distance_functions::euclidean_distance_function< impl::get_variable_type_t<FunctionSignature>>>
using local_rbf_interpolator
 Class to interpolate using RBF with length parameters localized for sample points.
 

Enumerations

enum class  kernel_matrix_type : std::uint8_t { dense , sparse }
 Enumeration of types of kernel matrices. More...
 

Functions

template<concepts::distance_function DistanceFunction, concepts::rbf RBF, concepts::length_parameter_calculator LengthParameterCalculator, base::concepts::dense_matrix KernelMatrix>
requires LengthParameterCalculator::uses_global_length_parameter
void compute_kernel_matrix (const DistanceFunction &distance_function, const RBF &rbf, LengthParameterCalculator &length_parameter_calculator, const std::vector< typename DistanceFunction::variable_type > &variables, KernelMatrix &kernel_matrix)
 Compute a kernel matrix.
 
template<concepts::distance_function DistanceFunction, concepts::rbf RBF, concepts::length_parameter_calculator LengthParameterCalculator, base::concepts::dense_matrix KernelMatrix>
requires std::is_same_v< typename LengthParameterCalculator::distance_function_type, DistanceFunction> && std::is_same_v<typename DistanceFunction::value_type, typename RBF::scalar_type> && std::is_same_v<typename DistanceFunction::value_type, typename KernelMatrix::Scalar> && (!LengthParameterCalculator::uses_global_length_parameter)
void compute_kernel_matrix (const DistanceFunction &distance_function, const RBF &rbf, LengthParameterCalculator &length_parameter_calculator, const std::vector< typename DistanceFunction::variable_type > &variables, KernelMatrix &kernel_matrix)
 Compute a kernel matrix.
 
template<concepts::distance_function DistanceFunction, concepts::csrbf RBF, concepts::length_parameter_calculator LengthParameterCalculator, base::concepts::sparse_matrix_of< typename DistanceFunction::value_type > KernelMatrix>
requires std::is_same_v< typename LengthParameterCalculator::distance_function_type, DistanceFunction> && std::is_same_v<typename DistanceFunction::value_type, typename RBF::scalar_type> && std::is_same_v<typename DistanceFunction::value_type, typename KernelMatrix::Scalar>
void compute_kernel_matrix (const DistanceFunction &distance_function, const RBF &rbf, LengthParameterCalculator &length_parameter_calculator, const std::vector< typename DistanceFunction::variable_type > &variables, KernelMatrix &kernel_matrix)
 Compute a kernel matrix.
 
template<index_type PolynomialDegree, base::concepts::real_scalar Variable, base::concepts::dense_matrix_of< Variable > Matrix>
requires (PolynomialDegree >= 0)
void compute_polynomial_term_matrix (const std::vector< Variable > &variables, Matrix &matrix)
 Compute a matrix of polynomial terms in RBF interpolation.
 
template<base::concepts::real_scalar Scalar>
auto generate_1d_halton_nodes (index_type num_nodes) -> std::vector< Scalar >
 Generate Halton nodes in 1 dimension [6].
 
template<base::concepts::real_scalar Scalar, index_type Dimensions>
auto generate_halton_nodes (index_type num_nodes) -> std::vector< Eigen::Vector< Scalar, Dimensions > >
 Generate Halton nodes [6].
 

Variables

constexpr auto rbf_interpolator_tag
 Tag of rbf_interpolator.
 
constexpr auto rbf_polynomial_interpolator_tag
 Tag of rbf_polynomial_interpolator.
 

Detailed Description

Namespace of RBF interpolation.

Typedef Documentation

◆ global_rbf_interpolator

template<typename FunctionSignature , concepts::rbf RBF = rbfs::gaussian_rbf<impl::get_default_scalar_type<FunctionSignature>>, kernel_matrix_type KernelMatrixType = kernel_matrix_type::dense, concepts::distance_function DistanceFunction = distance_functions::euclidean_distance_function< impl::get_variable_type_t<FunctionSignature>>>
using num_collect::rbf::global_rbf_interpolator
Initial value:
rbf_interpolator<FunctionSignature, RBF, KernelMatrixType, DistanceFunction,
length_parameter_calculators::global_length_parameter_calculator<
DistanceFunction>>

Class to interpolate using RBF with globally fixed length parameters.

Template Parameters
FunctionSignatureSignature of the function to interpolate. (Example: double(double), double(Eigen::Vector3d), ...)
RBFType of the RBF.
KernelMatrixTypeType of kernel matrices.
DistanceFunctionType of the distance function.

Definition at line 388 of file rbf_interpolator.h.

◆ global_rbf_polynomial_interpolator

template<typename FunctionSignature , concepts::rbf RBF = rbfs::gaussian_rbf<impl::get_default_scalar_type<FunctionSignature>>, index_type PolynomialDegree = 1, kernel_matrix_type KernelMatrixType = kernel_matrix_type::dense, concepts::distance_function DistanceFunction = distance_functions::euclidean_distance_function< impl::get_variable_type_t<FunctionSignature>>>
using num_collect::rbf::global_rbf_polynomial_interpolator
Initial value:
rbf_polynomial_interpolator<FunctionSignature, RBF, PolynomialDegree,
KernelMatrixType, DistanceFunction,
length_parameter_calculators::global_length_parameter_calculator<
DistanceFunction>>

Class to interpolate using RBF and polynomials with globally fixed length parameters.

Template Parameters
FunctionSignatureSignature of the function to interpolate. (Example: double(double), double(Eigen::Vector3d), ...)
RBFType of the RBF.
PolynomialDegreeDegree of the polynomial.
KernelMatrixTypeType of kernel matrices.
DistanceFunctionType of the distance function.

Definition at line 391 of file rbf_polynomial_interpolator.h.

◆ local_csrbf_interpolator

template<typename FunctionSignature , concepts::rbf RBF = rbfs::wendland_csrbf< impl::get_default_scalar_type<FunctionSignature>, 3, 1>, concepts::distance_function DistanceFunction = distance_functions::euclidean_distance_function< impl::get_variable_type_t<FunctionSignature>>>
using num_collect::rbf::local_csrbf_interpolator
Initial value:
rbf_interpolator<FunctionSignature, RBF,
kernel_matrix_type::sparse, DistanceFunction,
length_parameter_calculators::local_length_parameter_calculator<
DistanceFunction>>

Class to interpolate using compactly supported RBF with length parameters localized for sample points.

Template Parameters
FunctionSignatureSignature of the function to interpolate. (Example: double(double), double(Eigen::Vector3d), ...)
RBFType of the RBF.
DistanceFunctionType of the distance function.

Definition at line 408 of file rbf_interpolator.h.

◆ local_rbf_interpolator

template<typename FunctionSignature , concepts::rbf RBF = rbfs::gaussian_rbf<impl::get_default_scalar_type<FunctionSignature>>, kernel_matrix_type KernelMatrixType = kernel_matrix_type::dense, concepts::distance_function DistanceFunction = distance_functions::euclidean_distance_function< impl::get_variable_type_t<FunctionSignature>>>
using num_collect::rbf::local_rbf_interpolator
Initial value:
rbf_interpolator<FunctionSignature, RBF, KernelMatrixType, DistanceFunction,
length_parameter_calculators::local_length_parameter_calculator<
DistanceFunction>>

Class to interpolate using RBF with length parameters localized for sample points.

Template Parameters
FunctionSignatureSignature of the function to interpolate. (Example: double(double), double(Eigen::Vector3d), ...)
RBFType of the RBF.
KernelMatrixTypeType of kernel matrices.
DistanceFunctionType of the distance function.

Definition at line 367 of file rbf_interpolator.h.

Enumeration Type Documentation

◆ kernel_matrix_type

enum class num_collect::rbf::kernel_matrix_type : std::uint8_t
strong

Enumeration of types of kernel matrices.

Enumerator
dense 

Dense matrix.

sparse 

Sparse matrix.

Definition at line 29 of file kernel_matrix_type.h.

Function Documentation

◆ compute_kernel_matrix() [1/3]

template<concepts::distance_function DistanceFunction, concepts::rbf RBF, concepts::length_parameter_calculator LengthParameterCalculator, base::concepts::dense_matrix KernelMatrix>
requires LengthParameterCalculator::uses_global_length_parameter
void num_collect::rbf::compute_kernel_matrix ( const DistanceFunction & distance_function,
const RBF & rbf,
LengthParameterCalculator & length_parameter_calculator,
const std::vector< typename DistanceFunction::variable_type > & variables,
KernelMatrix & kernel_matrix )
inline

Compute a kernel matrix.

Template Parameters
DistanceFunctionType of the distance function.
RBFType of the RBF.
LengthParameterCalculatorType of the calculator of length parameters.
KernelMatrixType of the kernel matrix.
Parameters
[in]distance_functionDistance function.
[in]rbfRBF.
[in,out]length_parameter_calculatorCalculator of length parameters.
[in]variablesVariables.
[out]kernel_matrixKernel matrix.

Definition at line 66 of file compute_kernel_matrix.h.

◆ compute_kernel_matrix() [2/3]

template<concepts::distance_function DistanceFunction, concepts::rbf RBF, concepts::length_parameter_calculator LengthParameterCalculator, base::concepts::dense_matrix KernelMatrix>
requires std::is_same_v< typename LengthParameterCalculator::distance_function_type, DistanceFunction> && std::is_same_v<typename DistanceFunction::value_type, typename RBF::scalar_type> && std::is_same_v<typename DistanceFunction::value_type, typename KernelMatrix::Scalar> && (!LengthParameterCalculator::uses_global_length_parameter)
void num_collect::rbf::compute_kernel_matrix ( const DistanceFunction & distance_function,
const RBF & rbf,
LengthParameterCalculator & length_parameter_calculator,
const std::vector< typename DistanceFunction::variable_type > & variables,
KernelMatrix & kernel_matrix )
inline

Compute a kernel matrix.

Template Parameters
DistanceFunctionType of the distance function.
RBFType of the RBF.
LengthParameterCalculatorType of the calculator of length parameters.
KernelMatrixType of the kernel matrix.
Parameters
[in]distance_functionDistance function.
[in]rbfRBF.
[in,out]length_parameter_calculatorCalculator of length parameters.
[in]variablesVariables.
[out]kernel_matrixKernel matrix.

Definition at line 125 of file compute_kernel_matrix.h.

◆ compute_kernel_matrix() [3/3]

template<concepts::distance_function DistanceFunction, concepts::csrbf RBF, concepts::length_parameter_calculator LengthParameterCalculator, base::concepts::sparse_matrix_of< typename DistanceFunction::value_type > KernelMatrix>
requires std::is_same_v< typename LengthParameterCalculator::distance_function_type, DistanceFunction> && std::is_same_v<typename DistanceFunction::value_type, typename RBF::scalar_type> && std::is_same_v<typename DistanceFunction::value_type, typename KernelMatrix::Scalar>
void num_collect::rbf::compute_kernel_matrix ( const DistanceFunction & distance_function,
const RBF & rbf,
LengthParameterCalculator & length_parameter_calculator,
const std::vector< typename DistanceFunction::variable_type > & variables,
KernelMatrix & kernel_matrix )
inline

Compute a kernel matrix.

Template Parameters
DistanceFunctionType of the distance function.
RBFType of the RBF.
LengthParameterCalculatorType of the calculator of length parameters.
KernelMatrixType of the kernel matrix.
Parameters
[in]distance_functionDistance function.
[in]rbfRBF.
[in,out]length_parameter_calculatorCalculator of length parameters.
[in]variablesVariables.
[out]kernel_matrixKernel matrix.

Definition at line 178 of file compute_kernel_matrix.h.

◆ compute_polynomial_term_matrix()

template<index_type PolynomialDegree, base::concepts::real_scalar Variable, base::concepts::dense_matrix_of< Variable > Matrix>
requires (PolynomialDegree >= 0)
void num_collect::rbf::compute_polynomial_term_matrix ( const std::vector< Variable > & variables,
Matrix & matrix )
inline

Compute a matrix of polynomial terms in RBF interpolation.

Template Parameters
PolynomialDegreeDegree of polynomials.
VariableType of variables.
MatrixType of the matrix.
Parameters
[in]variablesVariables.
[out]matrixMatrix.

Definition at line 47 of file compute_polynomial_term_matrix.h.

◆ generate_1d_halton_nodes()

template<base::concepts::real_scalar Scalar>
auto num_collect::rbf::generate_1d_halton_nodes ( index_type num_nodes) -> std::vector<Scalar>
nodiscard

Generate Halton nodes in 1 dimension [6].

Template Parameters
ScalarType of scalars.
Parameters
[in]num_nodesNumber of nodes.
Returns
Generated nodes.

Definition at line 122 of file generate_halton_nodes.h.

◆ generate_halton_nodes()

template<base::concepts::real_scalar Scalar, index_type Dimensions>
auto num_collect::rbf::generate_halton_nodes ( index_type num_nodes) -> std::vector<Eigen::Vector<Scalar, Dimensions>>
nodiscard

Generate Halton nodes [6].

Template Parameters
ScalarType of scalars.
DimensionsNumber of dimensions of the space in which nodes are generated.
Parameters
[in]num_nodesNumber of nodes.
Returns
Generated nodes.

Definition at line 90 of file generate_halton_nodes.h.

Variable Documentation

◆ rbf_interpolator_tag

auto num_collect::rbf::rbf_interpolator_tag
constexpr
Initial value:
=
logging::log_tag_view("num_collect::rbf::rbf_interpolator")

Tag of rbf_interpolator.

Definition at line 52 of file rbf_interpolator.h.

◆ rbf_polynomial_interpolator_tag

auto num_collect::rbf::rbf_polynomial_interpolator_tag
constexpr
Initial value:
=
logging::log_tag_view("num_collect::rbf::rbf_polynomial_interpolator")

Tag of rbf_polynomial_interpolator.

Definition at line 54 of file rbf_polynomial_interpolator.h.