26#include <Eigen/Eigenvalues>
49template <
typename KernelValue,
typename FunctionValue,
60template <base::concepts::real_scalar KernelValue,
typename FunctionValue>
89 num_variables_ = kernel_matrix.rows();
91 "Kernel matrix must be a square matrix.");
93 "Matrix of additional terms must have the same number of rows as "
94 "the kernel matrix.");
95 num_additional_terms_ = additional_matrix.cols();
97 "The number of variables must be larger than the number of "
99 kernel_subspace_dimensions_ = num_variables_ - num_additional_terms_;
101 qr_decomposition_.compute(additional_matrix);
102 if (qr_decomposition_.rank() != additional_matrix.cols()) {
104 "The matrix of additional terms must have full "
105 "column rank. (columns={}, rand={})",
106 additional_matrix.cols(), qr_decomposition_.rank());
108 q_matrix_ = qr_decomposition_.householderQ();
110 transformed_kernel_matrix_ =
111 q_matrix_.rightCols(kernel_subspace_dimensions_).transpose() *
112 kernel_matrix * q_matrix_.rightCols(kernel_subspace_dimensions_);
114 eigen_value_decomposition_.compute(
115 transformed_kernel_matrix_, Eigen::ComputeEigenvectors);
116 data_transformation_matrix_ =
117 eigen_value_decomposition_.eigenvectors().transpose() *
118 q_matrix_.rightCols(kernel_subspace_dimensions_).transpose();
119 spectre_ = data_transformation_matrix_ * data;
121 kernel_matrix_ = &kernel_matrix;
137 "compute() must be called before solve().");
139 reg_param = correct_reg_param_if_needed(reg_param);
141 kernel_coefficients = data_transformation_matrix_.transpose() *
142 (eigen_value_decomposition_.eigenvalues().array() + reg_param)
148 additional_coefficients = qr_decomposition_.solve(
149 (*data_) - (*kernel_matrix_) * kernel_coefficients);
163 reg_param = correct_reg_param_if_needed(reg_param);
165 constexpr scalar_type limit = std::numeric_limits<scalar_type>::max() *
167 if (eigen_value_decomposition_.eigenvalues()(0) + reg_param <=
174 static_cast<scalar_type>(kernel_subspace_dimensions_) *
175 log(calc_reg_term(reg_param)) +
176 calc_log_determinant(reg_param);
192 return (spectre_.array().abs2().rowwise().sum() /
193 (eigen_value_decomposition_.eigenvalues().array() + reg_param))
206 return (eigen_value_decomposition_.eigenvalues().array() + reg_param)
220 eigen_value_decomposition_.eigenvalues()(0);
222 eigen_value_decomposition_.eigenvalues()(
223 eigen_value_decomposition_.eigenvalues().size() - 1);
225 largest_eigenvalue * std::numeric_limits<scalar_type>::epsilon();
227 eigenvalue_safe_limit - smallest_eigenvalue;
229 if (reg_param < reg_param_safe_limit) {
230 return reg_param_safe_limit;
236 Eigen::SelfAdjointEigenSolver<kernel_matrix_type>
237 eigen_value_decomposition_{};
240 Eigen::ColPivHouseholderQR<additional_matrix_type> qr_decomposition_{};
Class of exception on failure in algorithm.
KernelValue scalar_type
Type of scalars.
void compute(const kernel_matrix_type &kernel_matrix, const additional_matrix_type &additional_matrix, const vector_type &data)
Compute internal matrices.
auto calc_mle_objective(scalar_type reg_param) const -> scalar_type
Calculate maximum likelihood estimation (MLE) objective function scheuerer2011.
void solve(vector_type &kernel_coefficients, vector_type &additional_coefficients, scalar_type reg_param) const
Solve the linear equation with a regularization parameter.
auto calc_log_determinant(const scalar_type ®_param) const -> scalar_type
Calculate the logarithm of the determinant of kernel matrix plus regularization parameter.
Eigen::MatrixX< KernelValue > additional_matrix_type
Type of matrices of additional terms.
Eigen::MatrixX< KernelValue > kernel_matrix_type
Type of matrices of kernels.
auto correct_reg_param_if_needed(const scalar_type ®_param) const noexcept -> scalar_type
Correct regularization parameter if needed.
general_spline_equation_solver()=default
Constructor.
Eigen::VectorX< FunctionValue > vector_type
Type of vectors.
auto calc_reg_term(const scalar_type ®_param) const -> scalar_type
Calculate the regularization term.
Class to solve linear equations of kernel matrices and matrices of additional terms in RBF interpolat...
Definition of exceptions.
Definition of index_type type.
Definition of kernel_matrix_type enumeration.
Definition of macros for logging.
#define NUM_COLLECT_LOG_AND_THROW(EXCEPTION_TYPE,...)
Write an error log and throw an exception for an error.
std::ptrdiff_t index_type
Type of indices in this library.
Namespace of internal implementations.
kernel_matrix_type
Enumeration of types of kernel matrices.
Definition of NUM_COLLECT_PRECONDITION macro.
#define NUM_COLLECT_PRECONDITION(CONDITION,...)
Check whether a precondition is satisfied and throw an exception if not.
Definition of real_scalar concept.