numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
gauss_seidel_iterative_solver.h
Go to the documentation of this file.
1/*
2 * Copyright 2023 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// IWYU pragma: no_include <string_view>
23
24#include <cmath>
25
26#include <Eigen/Core>
27
36
37namespace num_collect::linear {
38
39template <base::concepts::sparse_matrix Matrix>
40class gauss_seidel_iterative_solver;
41
42namespace impl {
43
49template <base::concepts::sparse_matrix Matrix>
54
55} // namespace impl
56
63template <base::concepts::sparse_matrix Matrix>
65 : public iterative_solver_base<gauss_seidel_iterative_solver<Matrix>> {
66 static_assert(Matrix::IsRowMajor == 1, "Row major matrix is required.");
68 "Complex matrices are not supported.");
69
70public:
72 using base_type =
74
75 using typename base_type::matrix_type;
76 using typename base_type::real_scalar_type;
77 using typename base_type::scalar_type;
78 using typename base_type::storage_index_type;
79
80protected:
81 using base_type::coeff;
82
83public:
85 using vector_type = Eigen::VectorX<scalar_type>;
86
91
97 void compute(const matrix_type& coeff) {
99 diag_ = coeff.diagonal();
100 inv_diag_ = diag_.cwiseInverse();
101 NUM_COLLECT_PRECONDITION(inv_diag_.array().isFinite().all(),
102 "All diagonal elements of the coefficient matrix must not be "
103 "zero.");
104 }
105
114 template <base::concepts::dense_vector_of<scalar_type> Right,
115 base::concepts::dense_vector_of<scalar_type> Solution>
116 void solve_vector_in_place(const Right& right, Solution& solution) const {
117 const auto& coeff_ref = coeff();
118
119 NUM_COLLECT_PRECONDITION(coeff_ref.rows() == coeff_ref.cols(),
120 "Coefficient matrix must be a square matrix.");
121 NUM_COLLECT_PRECONDITION(right.rows() == coeff_ref.cols(),
122 "Right-hand-side vector must have the number of elements same as "
123 "the size of the coefficient matrix.");
124 NUM_COLLECT_PRECONDITION(solution.rows() == coeff_ref.cols(),
125 "Solution vector must have the number of elements same as the size "
126 "of the coefficient matrix.");
127
128 iterations_ = 0;
129 const scalar_type right_norm = right.squaredNorm();
131 while (iterations_ < max_iterations) {
132 iterate(coeff_ref, right, solution);
133 ++iterations_;
134 using std::sqrt;
135 residual_rate_ = sqrt(residual_ / right_norm);
137 break;
138 }
139 }
140 }
141
149 [[nodiscard]] auto iterations() const noexcept -> index_type {
150 return iterations_;
151 }
152
160 [[nodiscard]] auto residual_rate() const noexcept -> scalar_type {
161 return residual_rate_;
162 }
163
164private:
174 template <base::concepts::dense_vector_of<scalar_type> Right,
175 base::concepts::dense_vector_of<scalar_type> Solution>
176 void iterate(const matrix_type& coeff_ref, const Right& right,
177 Solution& solution) const {
178 const index_type size = coeff_ref.rows();
179 residual_ = static_cast<scalar_type>(0);
180 for (index_type i = 0; i < size; ++i) {
181 scalar_type numerator = right(i);
182 for (typename matrix_type::InnerIterator iter(coeff_ref, i); iter;
183 ++iter) {
184 if (iter.index() != i) {
185 numerator -= iter.value() * solution(iter.index());
186 }
187 }
188 const scalar_type row_residual = numerator - diag_(i) * solution(i);
189 solution(i) = numerator * inv_diag_(i);
190 residual_ += row_residual * row_residual;
191 }
192 }
193
196
199
202
205
208};
209
210} // namespace num_collect::linear
Class to solve linear equations using Gauss-Seidel iteration golub2013.
auto iterations() const noexcept -> index_type
Get the number of iterations.
void compute(const matrix_type &coeff)
Prepare to solve.
void iterate(const matrix_type &coeff_ref, const Right &right, Solution &solution) const
Iterate once.
Eigen::VectorX< scalar_type > vector_type
Type of vectors.
vector_type inv_diag_
Inverse of diagonal coefficients.
void solve_vector_in_place(const Right &right, Solution &solution) const
Iterate repeatedly until stop criterion is satisfied for a vector.
auto residual_rate() const noexcept -> scalar_type
Get the rate of the last residual.
auto compute(const matrix_type &coeff) -> gauss_seidel_iterative_solver< Matrix > &
Definition of dense_vector_of concept.
Definition of exceptions.
Definition of index_type type.
Definition of iterative_solver_base class.
Definition of macros for logging.
std::ptrdiff_t index_type
Type of indices in this library.
Definition index_type.h:33
Namespace of solvers of linear equations.
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.
Definition of sparse_matrix concept.