numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
avf4_formula.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
26#include "num_collect/constants/half.h" // IWYU pragma: keep
27#include "num_collect/constants/one.h" // IWYU pragma: keep
28#include "num_collect/constants/zero.h" // IWYU pragma: keep
37
38namespace num_collect::ode::avf {
39
46template <concepts::differentiable_problem Problem>
48public:
50 using problem_type = Problem;
51
53 using variable_type = typename problem_type::variable_type;
54
56 using scalar_type = typename problem_type::scalar_type;
57
59 using jacobian_type = typename problem_type::jacobian_type;
60
62 static constexpr index_type stages = 1;
63
65 static constexpr index_type order = 4;
66
68 static constexpr auto log_tag =
69 logging::log_tag_view("num_collect::ode::avf::avf4_formula");
70
78
87 void step(scalar_type time, scalar_type step_size,
88 const variable_type& current, variable_type& estimate) {
89 const auto dim = current.size();
90
91 problem().evaluate_on(time, current,
92 evaluation_type{.diff_coeff = true, .jacobian = true});
93 static constexpr scalar_type coeff_jacobi =
94 -static_cast<scalar_type>(1) / static_cast<scalar_type>(12);
95 jacobian_type coeff = step_size *
96 (jacobian_type::Identity(dim, dim) +
97 coeff_jacobi * step_size * step_size * problem().jacobian() *
98 problem().jacobian());
99
100 estimate = current + coeff * problem().diff_coeff();
101
102 integrand_.time(time);
103 integrand_.prev_var(current);
104 variable_type prev_estimate;
105 constexpr index_type max_loops = 10000;
106 for (index_type i = 0; i < max_loops; ++i) {
107 problem().evaluate_on(time,
108 constants::half<scalar_type> * (current + estimate),
109 evaluation_type{.diff_coeff = true, .jacobian = true});
110 coeff = step_size *
111 (jacobian_type::Identity(dim, dim) +
112 coeff_jacobi * step_size * step_size *
113 problem().jacobian() * problem().jacobian());
114
115 integrand_.next_var(estimate);
116 prev_estimate = estimate;
117 estimate = current +
118 coeff *
121 if (norm(estimate - prev_estimate) < tol_residual_norm_) {
122 return;
123 }
124 }
125 }
126
132 [[nodiscard]] auto problem() -> problem_type& {
133 return integrand_.problem();
134 }
135
141 [[nodiscard]] auto problem() const -> const problem_type& {
142 return integrand_.problem();
143 }
144
151 NUM_COLLECT_PRECONDITION(val > static_cast<scalar_type>(0),
152 "Tolerance of residual norm must be a positive value.");
153 tol_residual_norm_ = val;
154 }
155
156private:
159
161 static constexpr index_type integrator_degree = 5;
162
166
168 static constexpr auto default_tol_residual_norm =
169 static_cast<scalar_type>(1e-8);
170
173};
174
181template <concepts::differentiable_problem Problem>
183
190template <concepts::differentiable_problem Problem>
192
193} // namespace num_collect::ode::avf
Definition of avf_integrand class.
Class to perform numerical integration with Gauss-Legendre formula.
Class of tags of logs without memory management.
Class of 4th order average vector field (AVF) method quispel2008.
avf4_formula(const problem_type &problem=problem_type())
Constructor.
auto problem() const -> const problem_type &
Get the problem.
Problem problem_type
Type of problem.
static constexpr auto default_tol_residual_norm
Default tolerance of residual norm.
void step(scalar_type time, scalar_type step_size, const variable_type &current, variable_type &estimate)
Compute the next variable.
typename problem_type::variable_type variable_type
Type of variables.
auto problem() -> problem_type &
Get the problem.
impl::avf_integrand< problem_type > integrand_
Integrand.
static constexpr index_type order
Order of this formula.
void tol_residual_norm(scalar_type val)
Set tolerance of residual norm.
typename problem_type::scalar_type scalar_type
Type of scalars.
typename problem_type::jacobian_type jacobian_type
Type of Jacobian.
static constexpr auto log_tag
Log tag.
static constexpr index_type stages
Number of stages of this formula.
static constexpr index_type integrator_degree
Degree of integrator.
scalar_type tol_residual_norm_
Tolerance of residual norm.
integration::gauss_legendre_integrator< variable_type(scalar_type)> integrator_
Integrator.
Class of integrand for average vector field (AVF) method quispel2008.
void time(scalar_type val)
Set time.
void prev_var(const variable_type &var)
Set the previous variable.
void next_var(const variable_type &var)
Set the next variable.
auto problem() -> problem_type &
Get the problem.
Class of solvers of ODEs using embedded formulas.
Class of simple solver of ODEs.
Definition of differentiable_problem concept.
Definition of evaluation_type enumeration.
Definition of exceptions.
Definition of legendre_roots function.
Definition of half.
Definition of index_type type.
Definition of log_tag_view class.
Definition of macros for logging.
std::ptrdiff_t index_type
Type of indices in this library.
Definition index_type.h:33
auto norm(const Matrix &matrix)
Calculate norm of a matrix.
Definition norm.h:39
constexpr T zero
Value 0.
Definition zero.h:30
constexpr T half
Value 0.5.
Definition half.h:30
constexpr T one
Value 1.
Definition one.h:30
Namespace of average vector field (AVF) method.
Definition of non_embedded_formula_wrapper class.
Definition of norm class.
Definition of one.
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 simple_solver class.
Struct to specify types of evaluations.
Definition of zero.