numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
bicgstab_rosenbrock_equation_solver.h
Go to the documentation of this file.
1/*
2 * Copyright 2021 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 <Eigen/Core>
23// IWYU pragma: no_include <Eigen/SparseCore>
24
25#include <cmath> // IWYU pragma: keep
26#include <limits>
27#include <optional>
28
37
39
45template <concepts::multi_variate_problem Problem>
47public:
49 using problem_type = Problem;
50
52 using variable_type = typename problem_type::variable_type;
53
55 using scalar_type = typename problem_type::scalar_type;
56
58 static constexpr bool use_time_derivative =
60
63
71 const scalar_type& inverted_jacobian_coeff)
72 : inverted_jacobian_coeff_(inverted_jacobian_coeff) {}
73
83 const scalar_type& time, const scalar_type& step_size,
84 const variable_type& variable) {
85 problem_ = &problem;
86 time_ = time;
87 step_size_ = step_size;
88 variable_ = variable;
89
90 problem.evaluate_on(time, variable,
91 evaluation_type{.diff_coeff = true,
92 .time_derivative = use_time_derivative,
93 .mass = use_mass});
94 if constexpr (use_time_derivative) {
95 time_derivative_ = problem.time_derivative();
96 }
97 }
98
107 void apply_jacobian(const Target& target, Result& result) {
109 problem_ != nullptr, "evaluate_and_update_jacobian is not called.");
110
111 const scalar_type target_norm = target.norm();
112 if (target_norm < std::numeric_limits<scalar_type>::min()) {
113 result = variable_type::Zero(target.size());
114 return;
115 }
116 const scalar_type diff_width =
117 std::sqrt(std::numeric_limits<scalar_type>::epsilon()) /
118 target_norm;
119
120 problem_->evaluate_on(time_, variable_ + diff_width * target,
121 evaluation_type{.diff_coeff = true});
122 result = problem_->diff_coeff();
123
124 problem_->evaluate_on(time_, variable_ - diff_width * target,
125 evaluation_type{.diff_coeff = true});
126 result -= problem_->diff_coeff();
127 result /= static_cast<scalar_type>(2) * diff_width;
128 }
129
138 const scalar_type& coeff, variable_type& target) {
139 if constexpr (use_time_derivative) {
140 if (time_derivative_) {
141 target += step_size * coeff * (*time_derivative_);
142 }
143 }
144 }
145
152 void solve(const variable_type& rhs, variable_type& result) {
153 const auto coeff_function = [this](const auto& target, auto& result) {
154 this->apply_jacobian(target, result);
156 if constexpr (use_mass) {
157 result += problem_->mass() * target;
158 } else {
159 result += target;
160 }
161 };
162 result = variable_type::Zero(rhs.size());
163 bicgstab_.solve(coeff_function, rhs, result);
164 }
165
175 return *this;
176 }
177
178private:
181
184
187
190
192 std::optional<variable_type> time_derivative_{};
193
196
199};
200
201} // namespace num_collect::ode::rosenbrock
Definition of bicgstab class.
Class of error tolerances hairer1993.
Class to solve linear equations using BiCGstab golub2013.
Definition bicgstab.h:48
auto tolerances(const error_tolerances< vector_type > &val) -> bicgstab &
Set the error tolerances.
Definition bicgstab.h:141
void solve(CoeffFunction &&coeff_function, const vector_type &rhs, vector_type &solution)
Solve.
Definition bicgstab.h:80
Class to solve equations in Rosenbrock methods using BiCGstab.
std::optional< variable_type > time_derivative_
Partial derivative with respect to time.
scalar_type inverted_jacobian_coeff_
Coefficient multiplied to Jacobian matrices in inverted matrices.
void solve(const variable_type &rhs, variable_type &result)
Solve a linear equation.
auto tolerances(const error_tolerances< variable_type > &val) -> bicgstab_rosenbrock_equation_solver &
Set the error tolerances.
static constexpr bool use_time_derivative
Whether to use partial derivative with respect to time.
void add_time_derivative_term(const scalar_type &step_size, const scalar_type &coeff, variable_type &target)
Add a term of partial derivative with respect to time.
void evaluate_and_update_jacobian(problem_type &problem, const scalar_type &time, const scalar_type &step_size, const variable_type &variable)
Update Jacobian matrix and internal parameters.
bicgstab_rosenbrock_equation_solver(const scalar_type &inverted_jacobian_coeff)
Constructor.
void apply_jacobian(const Target &target, Result &result)
Multiply Jacobian matrix to a vector.
Concept of Eigen's dense vectors with real scalars.
Concept of problems of ordinary differential equations with mass.
Concept of problems of ordinary differential equations differentiable by time variable.
Definition of error_tolerances class.
Definition of evaluation_type enumeration.
Definition of mass_problem concept.
Definition of multi_variate_problem concept.
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_dense_vector concept.
Struct to specify types of evaluations.
Definition of time_differentiable_problem concept.