numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
ros3w_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
30
32
40template <concepts::problem Problem,
41 concepts::rosenbrock_equation_solver EquationSolver =
44 : public rosenbrock_formula_base<ros3w_formula<Problem, EquationSolver>,
45 Problem, EquationSolver> {
46public:
48 using base_type =
50 EquationSolver>;
51
53 using typename base_type::problem_type;
54 using typename base_type::scalar_type;
55 using typename base_type::variable_type;
56
60
61protected:
62 using base_type::coeff;
63
64public:
66 static constexpr index_type stages = 3;
67
69 static constexpr index_type order = 3;
70
72 static constexpr index_type lesser_order = 2;
73
75 static constexpr auto log_tag =
76 logging::log_tag_view("num_collect::ode::rosenbrock::ros3w_formula");
77
91 static constexpr scalar_type a21 = coeff(6.666666666666666e-01);
92 static constexpr scalar_type a31 = coeff(6.666666666666666e-01);
93 static constexpr scalar_type a32 = coeff(0);
94
95 static constexpr scalar_type b1 = coeff(0);
96 static constexpr scalar_type b2 = a21;
97 static constexpr scalar_type b3 = a31 + a32;
98
99 static constexpr scalar_type g21 = coeff(3.635068368900681e-01);
100 static constexpr scalar_type g31 = coeff(-8.996866791992636e-01);
101 static constexpr scalar_type g32 = coeff(-1.537997822626885e-01);
102 static constexpr scalar_type g = coeff(4.358665215084590e-01);
103
104 static constexpr scalar_type g1 = g;
105 static constexpr scalar_type g2 = g21 + g;
106 static constexpr scalar_type g3 = g31 + g32 + g;
107
108 static constexpr scalar_type c1 = coeff(2.500000000000000e-01);
109 static constexpr scalar_type c2 = coeff(2.500000000000000e-01);
110 static constexpr scalar_type c3 = coeff(5.000000000000000e-01);
111
112 static constexpr scalar_type cw1 = coeff(7.467047032740110e-01);
113 static constexpr scalar_type cw2 = coeff(1.144064078371002e-01);
114 static constexpr scalar_type cw3 = coeff(1.388888888888889e-01);
115
116 static constexpr scalar_type ce1 = c1 - cw1;
117 static constexpr scalar_type ce2 = c2 - cw2;
118 static constexpr scalar_type ce3 = c3 - cw3;
120
128
130 void step(scalar_type time, scalar_type step_size,
131 const variable_type& current, variable_type& estimate) {
132 variable_type unused;
133 step_embedded(time, step_size, current, estimate, unused);
134 }
135
147 const variable_type& current, variable_type& estimate,
148 variable_type& error) {
149 equation_solver().evaluate_and_update_jacobian(
150 problem(), time, step_size, current);
151
152 // 1st stage
153 temp_rhs_ = problem().diff_coeff();
154 equation_solver().add_time_derivative_term(step_size, g1, temp_rhs_);
155 equation_solver().solve(temp_rhs_, k1_);
156
157 // 2nd stage
158 temp_var_ = g21 * k1_;
159 equation_solver().apply_jacobian(temp_var_, temp_rhs_);
160 temp_rhs_ *= step_size;
161 temp_var_ = current + step_size * (a21 * k1_);
162 problem().evaluate_on(time + b2 * step_size, temp_var_,
163 evaluation_type{.diff_coeff = true});
164 temp_rhs_ += problem().diff_coeff();
165 equation_solver().add_time_derivative_term(step_size, g2, temp_rhs_);
166 equation_solver().solve(temp_rhs_, k2_);
167
168 // 3rd stage
169 temp_var_ = g31 * k1_ + g32 * k2_;
170 equation_solver().apply_jacobian(temp_var_, temp_rhs_);
171 temp_rhs_ *= step_size;
172 temp_var_ = current + step_size * (a31 * k1_ + a32 * k2_);
173 problem().evaluate_on(time + b3 * step_size, temp_var_,
174 evaluation_type{.diff_coeff = true});
175 temp_rhs_ += problem().diff_coeff();
176 equation_solver().add_time_derivative_term(step_size, g3, temp_rhs_);
177 equation_solver().solve(temp_rhs_, k3_);
178
179 estimate = current + step_size * (c1 * k1_ + c2 * k2_ + c3 * k3_);
180 error = step_size * (ce1 * k1_ + ce2 * k2_ + ce3 * k3_);
181 }
182
183private:
193
196
199};
200
206template <concepts::problem Problem>
208
209} // namespace num_collect::ode::rosenbrock
Class of tags of logs without memory management.
Class of solvers of ODEs using embedded formulas.
Base class of formulas in ODE solvers.
auto problem() -> problem_type &
Get the problem.
Class of ROS3w formula in rang2005 for Rosenbrock method.
variable_type temp_rhs_
Temporary right-hand-side vector.
static constexpr index_type lesser_order
Order of lesser coefficients of this formula.
static constexpr auto log_tag
Log tag.
void step(scalar_type time, scalar_type step_size, const variable_type &current, variable_type &estimate)
Compute the next variable.
static constexpr index_type order
Order of this formula.
variable_type temp_var_
Temporary variable.
ros3w_formula(const problem_type &problem)
Constructor.
static constexpr index_type stages
Number of stages of this formula.
auto equation_solver() const -> const equation_solver_type &
Access the solver of equations in Rosenbrock methods.
void step_embedded(scalar_type time, scalar_type step_size, const variable_type &current, variable_type &estimate, variable_type &error)
Compute the next variable and weak estimate of it with embedded formula.
Base class of formulas in Rosenbrock method.
typename problem_type::variable_type variable_type
Type of variables.
typename problem_type::scalar_type scalar_type
Type of scalars.
EquationSolver equation_solver_type
Type of class to solve equations in Rosenbrock methods.
auto equation_solver() const -> const equation_solver_type &
Access the solver of equations in Rosenbrock methods.
static constexpr auto coeff(T val) -> scalar_type
Convert coefficients.
Definition of default_rosenbrock_equation_solver class.
Definition of embedded_solver class.
Definition of evaluation_type enumeration.
Definition of index_type type.
Definition of log_tag_view class.
std::ptrdiff_t index_type
Type of indices in this library.
Definition index_type.h:33
typename default_rosenbrock_equation_solver< Problem >::type default_rosenbrock_equation_solver_t
Default class to solve equations in Rosenbrock methods.
Definition of problem concept.
Definition of rosenbrock_equation_solver concept.
Definition of rosenbrock_formula_base class.
Struct to specify types of evaluations.