numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
kernel_matrix_solver.h
Go to the documentation of this file.
1/*
2 * Copyright 2024 MusicScience37 (Kenta Kabashima)
3 *
4 * Licensed under the Apache License, Version 2.0 (the "License");
5 * you may not use this file except in compliance with the License.
6 * You may obtain a copy of the License at
7 *
8 * http://www.apache.org/licenses/LICENSE-2.0
9 *
10 * Unless required by applicable law or agreed to in writing, software
11 * distributed under the License is distributed on an "AS IS" BASIS,
12 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13 * See the License for the specific language governing permissions and
14 * limitations under the License.
15 */
20#pragma once
21
22#include <cmath>
23#include <limits>
24
25#include <Eigen/Core>
26#include <Eigen/Eigenvalues>
27#include <Eigen/IterativeLinearSolvers>
28#include <Eigen/LU>
29#include <Eigen/SparseCore>
30#include <Eigen/src/Core/util/Constants.h>
31
37
38namespace num_collect::rbf::impl {
39
49template <typename KernelValue, typename FunctionValue,
50 kernel_matrix_type KernelMatrixType, bool UsesGlobalLengthParameter>
52
59template <typename KernelValue, typename FunctionValue>
60class kernel_matrix_solver<KernelValue, FunctionValue,
62public:
64 using kernel_matrix_type = Eigen::MatrixX<KernelValue>;
65
67 using vector_type = Eigen::VectorX<FunctionValue>;
68
70 using scalar_type = KernelValue;
71
74
81 void compute(
82 const kernel_matrix_type& kernel_matrix, const vector_type& data) {
83 decomposition_.compute(kernel_matrix, Eigen::ComputeEigenvectors);
84 spectre_ = decomposition_.eigenvectors().adjoint() * data;
85 }
86
94 void solve(vector_type& coefficients, scalar_type reg_param,
95 const vector_type& data) const {
96 (void)data; // not used in this version.
97
98 reg_param = correct_reg_param_if_needed(reg_param);
99
100 coefficients = decomposition_.eigenvectors() *
101 (decomposition_.eigenvalues().array() + reg_param)
102 .inverse()
103 .matrix()
104 .asDiagonal() *
105 spectre_;
106 }
107
117 [[nodiscard]] auto calc_mle_objective(scalar_type reg_param) const
118 -> scalar_type {
119 reg_param = correct_reg_param_if_needed(reg_param);
120
121 constexpr scalar_type limit = std::numeric_limits<scalar_type>::max() *
122 static_cast<scalar_type>(1e-20);
123 if (decomposition_.eigenvalues()(0) + reg_param <=
124 static_cast<scalar_type>(0)) {
125 return limit;
126 }
127
128 using std::log;
129 const scalar_type value = static_cast<scalar_type>(spectre_.rows()) *
130 log(calc_reg_term(reg_param)) +
131 calc_log_determinant(reg_param);
132 if (value < limit) {
133 return value;
134 }
135 return limit;
136 }
137
144 [[nodiscard]] auto calc_common_coeff(scalar_type reg_param) const
145 -> scalar_type {
146 reg_param = correct_reg_param_if_needed(reg_param);
147
148 return calc_reg_term(reg_param) /
149 static_cast<scalar_type>(spectre_.rows());
150 }
151
160 template <base::concepts::dense_vector InputData>
161 [[nodiscard]] auto calc_reg_term(
162 const InputData& data, scalar_type reg_param) const -> scalar_type {
163 reg_param = correct_reg_param_if_needed(reg_param);
164
165 return ((decomposition_.eigenvectors().adjoint() * data)
166 .array()
167 .abs2()
168 .rowwise()
169 .sum() /
170 (decomposition_.eigenvalues().array() + reg_param))
171 .sum();
172 }
173
179 [[nodiscard]] auto eigenvalues() const -> decltype(auto) {
180 return decomposition_.eigenvalues();
181 }
182
183private:
190 [[nodiscard]] auto calc_reg_term(const scalar_type& reg_param) const
191 -> scalar_type {
192 return (spectre_.array().abs2().rowwise().sum() /
193 (decomposition_.eigenvalues().array() + reg_param))
194 .sum();
195 }
196
204 [[nodiscard]] auto calc_log_determinant(const scalar_type& reg_param) const
205 -> scalar_type {
206 return (decomposition_.eigenvalues().array() + reg_param).log().sum();
207 }
208
215 [[nodiscard]] auto correct_reg_param_if_needed(
216 const scalar_type& reg_param) const noexcept -> scalar_type {
217 const scalar_type smallest_eigenvalue = decomposition_.eigenvalues()(0);
218 const scalar_type largest_eigenvalue = decomposition_.eigenvalues()(
219 decomposition_.eigenvalues().size() - 1);
220 const scalar_type eigenvalue_safe_limit =
221 largest_eigenvalue * std::numeric_limits<scalar_type>::epsilon();
222 const scalar_type reg_param_safe_limit =
223 eigenvalue_safe_limit - smallest_eigenvalue;
224
225 if (reg_param < reg_param_safe_limit) {
226 return reg_param_safe_limit;
227 }
228 return reg_param;
229 }
230
232 Eigen::SelfAdjointEigenSolver<kernel_matrix_type> decomposition_;
233
236};
237
244template <typename KernelValue, typename FunctionValue>
245class kernel_matrix_solver<KernelValue, FunctionValue,
246 kernel_matrix_type::dense, false> {
247public:
249 using kernel_matrix_type = Eigen::MatrixX<KernelValue>;
250
252 using vector_type = Eigen::VectorX<FunctionValue>;
253
255 using scalar_type = KernelValue;
256
259
267 const kernel_matrix_type& kernel_matrix, const vector_type& data) {
268 (void)data; // not used in this version.
269
270 solver_.compute(kernel_matrix);
271 }
272
280 void solve(vector_type& coefficients, scalar_type reg_param,
281 const vector_type& data) const {
282 NUM_COLLECT_PRECONDITION(reg_param == static_cast<scalar_type>(0),
283 "Non-zero regularization parameter cannot be used in this "
284 "configuration.");
285
286 coefficients = solver_.solve(data);
287 }
288
289private:
291 Eigen::PartialPivLU<kernel_matrix_type> solver_;
292};
293
300template <typename KernelValue, typename FunctionValue>
301class kernel_matrix_solver<KernelValue, FunctionValue,
302 kernel_matrix_type::sparse, false> {
303public:
306 Eigen::SparseMatrix<KernelValue, Eigen::RowMajor>;
307
309 using vector_type = Eigen::VectorX<FunctionValue>;
310
312 using scalar_type = KernelValue;
313
316
324 const kernel_matrix_type& kernel_matrix, const vector_type& data) {
325 (void)data; // not used in this version.
326
327 solver_.compute(kernel_matrix);
328 }
329
337 void solve(vector_type& coefficients, scalar_type reg_param,
338 const vector_type& data) const {
339 NUM_COLLECT_PRECONDITION(reg_param == static_cast<scalar_type>(0),
340 "Non-zero regularization parameter cannot be used in this "
341 "configuration.");
342
343 coefficients = solver_.solve(data);
344 }
345
346private:
348 Eigen::BiCGSTAB<kernel_matrix_type> solver_;
349};
350
351} // namespace num_collect::rbf::impl
void compute(const kernel_matrix_type &kernel_matrix, const vector_type &data)
Compute internal matrices.
void solve(vector_type &coefficients, scalar_type reg_param, const vector_type &data) const
Solve the linear equation with a regularization parameter.
void compute(const kernel_matrix_type &kernel_matrix, const vector_type &data)
Compute internal matrices.
void solve(vector_type &coefficients, scalar_type reg_param, const vector_type &data) const
Solve the linear equation with a regularization parameter.
auto calc_log_determinant(const scalar_type &reg_param) const -> scalar_type
Calculate the logarithm of the determinant of kernel matrix plus regularization parameter.
auto calc_reg_term(const scalar_type &reg_param) const -> scalar_type
Calculate the regularization term.
auto calc_common_coeff(scalar_type reg_param) const -> scalar_type
Calculate the coefficient of the kernel common in variables.
auto calc_reg_term(const InputData &data, scalar_type reg_param) const -> scalar_type
Calculate the regularization term for a vector.
auto correct_reg_param_if_needed(const scalar_type &reg_param) const noexcept -> scalar_type
Correct regularization parameter if needed.
void compute(const kernel_matrix_type &kernel_matrix, const vector_type &data)
Compute internal matrices.
void solve(vector_type &coefficients, scalar_type reg_param, const vector_type &data) const
Solve the linear equation with a regularization parameter.
Eigen::SelfAdjointEigenSolver< kernel_matrix_type > decomposition_
Eigen decomposition of the kernel matrix.
auto calc_mle_objective(scalar_type reg_param) const -> scalar_type
Calculate maximum likelihood estimation (MLE) objective function scheuerer2011.
Class to solve linear equations of kernel matrices.
Definition of dense_vector concept.
Definition of exceptions.
Definition of kernel_matrix_type enumeration.
Definition of macros for logging.
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.