numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
rodaspr_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<rodaspr_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 = 6;
67
69 static constexpr index_type order = 4;
70
76 static constexpr index_type lesser_order = 3;
77
79 static constexpr auto log_tag =
80 logging::log_tag_view("num_collect::ode::rosenbrock::rodaspr_formula");
81
95 static constexpr scalar_type a21 = coeff(0.75);
96 static constexpr scalar_type a31 = coeff(7.5162877593868457e-2);
97 static constexpr scalar_type a32 = coeff(2.4837122406131545e-2);
98 static constexpr scalar_type a41 = coeff(1.6532708886396510);
99 static constexpr scalar_type a42 = coeff(0.21545706385445562);
100 static constexpr scalar_type a43 = coeff(-1.3157488872766792);
101 static constexpr scalar_type a51 = coeff(19.385003738039885);
102 static constexpr scalar_type a52 = coeff(1.2007117225835324);
103 static constexpr scalar_type a53 = coeff(-19.337924059522791);
104 static constexpr scalar_type a54 = coeff(-0.24779140110062559);
105 static constexpr scalar_type a61 = coeff(-7.3844531665375115);
106 static constexpr scalar_type a62 = coeff(-0.30593419030174646);
107 static constexpr scalar_type a63 = coeff(7.8622074209377981);
108 static constexpr scalar_type a64 = coeff(0.57817993590145966);
109 static constexpr scalar_type a65 = coeff(0.25);
110
111 static constexpr scalar_type b1 = coeff(0);
112 static constexpr scalar_type b2 = a21;
113 static constexpr scalar_type b3 = a31 + a32;
114 static constexpr scalar_type b4 = a41 + a42 + a43;
115 static constexpr scalar_type b5 = a51 + a52 + a53 + a54;
116 static constexpr scalar_type b6 = a61 + a62 + a63 + a64 + a65;
117
118 static constexpr scalar_type g21 = coeff(-0.75);
119 static constexpr scalar_type g31 = coeff(-8.8644359075349941e-2);
120 static constexpr scalar_type g32 = coeff(-2.8688974257983398e-2);
121 static constexpr scalar_type g41 = coeff(-4.8470034585330284);
122 static constexpr scalar_type g42 = coeff(-0.31583244269672095);
123 static constexpr scalar_type g43 = coeff(4.9536568360123221);
124 static constexpr scalar_type g51 = coeff(-26.769456904577400);
125 static constexpr scalar_type g52 = coeff(-1.5066459128852787);
126 static constexpr scalar_type g53 = coeff(27.200131480460591);
127 static constexpr scalar_type g54 = coeff(0.82597133700208525);
128 static constexpr scalar_type g61 = coeff(6.5876206496361416);
129 static constexpr scalar_type g62 = coeff(0.36807059172993878);
130 static constexpr scalar_type g63 = coeff(-6.7423520694658121);
131 static constexpr scalar_type g64 = coeff(-0.10619631475741095);
132 static constexpr scalar_type g65 = coeff(-0.35714285714285715);
133 static constexpr scalar_type g = coeff(0.25);
134
135 static constexpr scalar_type g1 = g;
136 static constexpr scalar_type g2 = g21 + g;
137 static constexpr scalar_type g3 = g31 + g32 + g;
138 static constexpr scalar_type g4 = g41 + g42 + g43 + g;
139 static constexpr scalar_type g5 = g51 + g52 + g53 + g54 + g;
140 static constexpr scalar_type g6 = g61 + g62 + g63 + g64 + g65 + g;
141
142 static constexpr scalar_type c1 = coeff(-0.79683251690137014);
143 static constexpr scalar_type c2 = coeff(6.2136401428192344e-02);
144 static constexpr scalar_type c3 = coeff(1.1198553514719862);
145 static constexpr scalar_type c4 = coeff(0.47198362114404874);
146 static constexpr scalar_type c5 = coeff(-0.10714285714285714);
147 static constexpr scalar_type c6 = coeff(0.25);
148
149 static constexpr scalar_type cw1 = coeff(-7.3844531665375115);
150 static constexpr scalar_type cw2 = coeff(-0.30593419030174646);
151 static constexpr scalar_type cw3 = coeff(7.8622074209377981);
152 static constexpr scalar_type cw4 = coeff(0.57817993590145966);
153 static constexpr scalar_type cw5 = coeff(0.25);
154
155 static constexpr scalar_type ce1 = c1 - cw1;
156 static constexpr scalar_type ce2 = c2 - cw2;
157 static constexpr scalar_type ce3 = c3 - cw3;
158 static constexpr scalar_type ce4 = c4 - cw4;
159 static constexpr scalar_type ce5 = c5 - cw5;
160 static constexpr scalar_type ce6 = c6;
162
170
172 void step(scalar_type time, scalar_type step_size,
173 const variable_type& current, variable_type& estimate) {
174 variable_type unused;
175 step_embedded(time, step_size, current, estimate, unused);
176 }
177
189 const variable_type& current, variable_type& estimate,
190 variable_type& error) {
191 equation_solver().evaluate_and_update_jacobian(
192 problem(), time, step_size, current);
193
194 // 1st stage
195 temp_rhs_ = problem().diff_coeff();
196 equation_solver().add_time_derivative_term(step_size, g1, temp_rhs_);
197 equation_solver().solve(temp_rhs_, k1_);
198
199 // 2nd stage
200 temp_var_ = g21 * k1_;
201 equation_solver().apply_jacobian(temp_var_, temp_rhs_);
202 temp_rhs_ *= step_size;
203 temp_var_ = current + step_size * (a21 * k1_);
204 problem().evaluate_on(time + b2 * step_size, temp_var_,
205 evaluation_type{.diff_coeff = true});
206 temp_rhs_ += problem().diff_coeff();
207 equation_solver().add_time_derivative_term(step_size, g2, temp_rhs_);
208 equation_solver().solve(temp_rhs_, k2_);
209
210 // 3rd stage
211 temp_var_ = g31 * k1_ + g32 * k2_;
212 equation_solver().apply_jacobian(temp_var_, temp_rhs_);
213 temp_rhs_ *= step_size;
214 temp_var_ = current + step_size * (a31 * k1_ + a32 * k2_);
215 problem().evaluate_on(time + b3 * step_size, temp_var_,
216 evaluation_type{.diff_coeff = true});
217 temp_rhs_ += problem().diff_coeff();
218 equation_solver().add_time_derivative_term(step_size, g3, temp_rhs_);
219 equation_solver().solve(temp_rhs_, k3_);
220
221 // 4th stage
222 temp_var_ = g41 * k1_ + g42 * k2_ + g43 * k3_;
223 equation_solver().apply_jacobian(temp_var_, temp_rhs_);
224 temp_rhs_ *= step_size;
225 temp_var_ = current + step_size * (a41 * k1_ + a42 * k2_ + a43 * k3_);
226 problem().evaluate_on(time + b4 * step_size, temp_var_,
227 evaluation_type{.diff_coeff = true});
228 temp_rhs_ += problem().diff_coeff();
229 equation_solver().add_time_derivative_term(step_size, g4, temp_rhs_);
230 equation_solver().solve(temp_rhs_, k4_);
231
232 // 5th stage
233 temp_var_ = g51 * k1_ + g52 * k2_ + g53 * k3_ + g54 * k4_;
234 equation_solver().apply_jacobian(temp_var_, temp_rhs_);
235 temp_rhs_ *= step_size;
236 temp_var_ = current +
237 step_size * (a51 * k1_ + a52 * k2_ + a53 * k3_ + a54 * k4_);
238 problem().evaluate_on(time + b5 * step_size, temp_var_,
239 evaluation_type{.diff_coeff = true});
240 temp_rhs_ += problem().diff_coeff();
241 equation_solver().add_time_derivative_term(step_size, g5, temp_rhs_);
242 equation_solver().solve(temp_rhs_, k5_);
243
244 // 6th stage
245 temp_var_ = g61 * k1_ + g62 * k2_ + g63 * k3_ + g64 * k4_ + g65 * k5_;
246 equation_solver().apply_jacobian(temp_var_, temp_rhs_);
247 temp_rhs_ *= step_size;
248 temp_var_ = current +
249 step_size *
250 (a61 * k1_ + a62 * k2_ + a63 * k3_ + a64 * k4_ + a65 * k5_);
251 problem().evaluate_on(time + b6 * step_size, temp_var_,
252 evaluation_type{.diff_coeff = true});
253 temp_rhs_ += problem().diff_coeff();
254 equation_solver().add_time_derivative_term(step_size, g6, temp_rhs_);
255 equation_solver().solve(temp_rhs_, k6_);
256
257 estimate = current +
258 step_size *
259 (c1 * k1_ + c2 * k2_ + c3 * k3_ + c4 * k4_ + c5 * k5_ +
260 c6 * k6_);
261 error = step_size *
262 (ce1 * k1_ + ce2 * k2_ + ce3 * k3_ + ce4 * k4_ + ce5 * k5_ +
263 ce6 * k6_);
264 }
265
266private:
279
282
285};
286
292template <concepts::problem Problem>
294
295} // 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 RODASPR formula rang2015 for Rosenbrock method.
static constexpr index_type lesser_order
Order of lesser coefficients of this formula.
variable_type temp_var_
Temporary 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 stages
Number of stages of this formula.
variable_type temp_rhs_
Temporary right-hand-side vector.
rodaspr_formula(const problem_type &problem)
Constructor.
void step(scalar_type time, scalar_type step_size, const variable_type &current, variable_type &estimate)
Compute the next variable.
auto equation_solver() const -> const equation_solver_type &
Access the solver of equations in Rosenbrock methods.
static constexpr index_type order
Order of this 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.