numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
ros34pw3_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<ros34pw3_formula<Problem, EquationSolver>,
45 Problem, EquationSolver> {
46public:
48 using base_type =
50 Problem, 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 = 4;
67
69 static constexpr index_type order = 4;
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::ros34pw3_formula");
77
91 static constexpr scalar_type a21 = coeff(2.5155456020628817e+00);
92 static constexpr scalar_type a31 = coeff(5.0777280103144085e-01);
93 static constexpr scalar_type a32 = coeff(7.5000000000000000e-01);
94 static constexpr scalar_type a41 = coeff(1.3959081404277204e-01);
95 static constexpr scalar_type a42 = coeff(-3.3111001065419338e-01);
96 static constexpr scalar_type a43 = coeff(8.2040559712714178e-01);
97
98 static constexpr scalar_type b1 = coeff(0);
99 static constexpr scalar_type b2 = a21;
100 static constexpr scalar_type b3 = a31 + a32;
101 static constexpr scalar_type b4 = a41 + a42 + a43;
102
103 static constexpr scalar_type g21 = coeff(-2.5155456020628817e+00);
104 static constexpr scalar_type g31 = coeff(-8.7991339217106512e-01);
105 static constexpr scalar_type g32 = coeff(-9.6014187766190695e-01);
106 static constexpr scalar_type g41 = coeff(-4.1731389379448741e-01);
107 static constexpr scalar_type g42 = coeff(4.1091047035857703e-01);
108 static constexpr scalar_type g43 = coeff(-1.3558873204765276e+00);
109 static constexpr scalar_type g = coeff(1.0685790213016289e+00);
110
111 static constexpr scalar_type g1 = g;
112 static constexpr scalar_type g2 = g21 + g;
113 static constexpr scalar_type g3 = g31 + g32 + g;
114 static constexpr scalar_type g4 = g41 + g42 + g43 + g;
115
116 static constexpr scalar_type c1 = coeff(2.2047681286931747e-01);
117 static constexpr scalar_type c2 = coeff(2.7828278331185935e-03);
118 static constexpr scalar_type c3 = coeff(7.1844787635140066e-03);
119 static constexpr scalar_type c4 = coeff(7.6955588053404989e-01);
120
121 static constexpr scalar_type cw1 = coeff(3.1300297285209688e-01);
122 static constexpr scalar_type cw2 = coeff(-2.8946895245112692e-01);
123 static constexpr scalar_type cw3 = coeff(9.7646597959903003e-01);
124
125 static constexpr scalar_type ce1 = c1 - cw1;
126 static constexpr scalar_type ce2 = c2 - cw2;
127 static constexpr scalar_type ce3 = c3 - cw3;
128 static constexpr scalar_type ce4 = c4;
130
138
140 void step(scalar_type time, scalar_type step_size,
141 const variable_type& current, variable_type& estimate) {
142 variable_type unused;
143 step_embedded(time, step_size, current, estimate, unused);
144 }
145
157 const variable_type& current, variable_type& estimate,
158 variable_type& error) {
159 equation_solver().evaluate_and_update_jacobian(
160 problem(), time, step_size, current);
161
162 // 1st stage
163 temp_rhs_ = problem().diff_coeff();
164 equation_solver().add_time_derivative_term(step_size, g1, temp_rhs_);
165 equation_solver().solve(temp_rhs_, k1_);
166
167 // 2nd stage
168 temp_var_ = g21 * k1_;
169 equation_solver().apply_jacobian(temp_var_, temp_rhs_);
170 temp_rhs_ *= step_size;
171 temp_var_ = current + step_size * (a21 * k1_);
172 problem().evaluate_on(time + b2 * step_size, temp_var_,
173 evaluation_type{.diff_coeff = true});
174 temp_rhs_ += problem().diff_coeff();
175 equation_solver().add_time_derivative_term(step_size, g2, temp_rhs_);
176 equation_solver().solve(temp_rhs_, k2_);
177
178 // 3rd stage
179 temp_var_ = g31 * k1_ + g32 * k2_;
180 equation_solver().apply_jacobian(temp_var_, temp_rhs_);
181 temp_rhs_ *= step_size;
182 temp_var_ = current + step_size * (a31 * k1_ + a32 * k2_);
183 problem().evaluate_on(time + b3 * step_size, temp_var_,
184 evaluation_type{.diff_coeff = true});
185 temp_rhs_ += problem().diff_coeff();
186 equation_solver().add_time_derivative_term(step_size, g3, temp_rhs_);
187 equation_solver().solve(temp_rhs_, k3_);
188
189 // 4th stage
190 temp_var_ = g41 * k1_ + g42 * k2_ + g43 * k3_;
191 equation_solver().apply_jacobian(temp_var_, temp_rhs_);
192 temp_rhs_ *= step_size;
193 temp_var_ = current + step_size * (a41 * k1_ + a42 * k2_ + a43 * k3_);
194 problem().evaluate_on(time + b4 * step_size, temp_var_,
195 evaluation_type{.diff_coeff = true});
196 temp_rhs_ += problem().diff_coeff();
197 equation_solver().add_time_derivative_term(step_size, g4, temp_rhs_);
198 equation_solver().solve(temp_rhs_, k4_);
199
200 estimate =
201 current + step_size * (c1 * k1_ + c2 * k2_ + c3 * k3_ + c4 * k4_);
202 error = step_size * (ce1 * k1_ + ce2 * k2_ + ce3 * k3_ + ce4 * k4_);
203 }
204
205private:
216
219
222};
223
229template <concepts::problem Problem>
231
232} // 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 ROS34PW3 formula in rang2005 for Rosenbrock method.
static constexpr index_type order
Order of this formula.
static constexpr index_type stages
Number of stages of this formula.
variable_type temp_rhs_
Temporary right-hand-side vector.
void step(scalar_type time, scalar_type step_size, const variable_type &current, variable_type &estimate)
Compute the next variable.
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.
static constexpr index_type lesser_order
Order of lesser coefficients of this formula.
variable_type temp_var_
Temporary variable.
auto equation_solver() const -> const equation_solver_type &
Access the solver of equations in Rosenbrock methods.
ros34pw3_formula(const problem_type &problem)
Constructor.
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.