numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
implicit_euler_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
23#include "num_collect/constants/one.h" // IWYU pragma: keep
30
32
39template <concepts::problem Problem,
40 concepts::slope_equation_solver FormulaSolver =
41 inexact_newton_slope_equation_solver<Problem>>
43 : public implicit_formula_base<
44 implicit_euler_formula<Problem, FormulaSolver>, Problem,
45 FormulaSolver> {
46public:
49
51 using base_type =
53 Problem, FormulaSolver>;
54
55 using typename base_type::problem_type;
56 using typename base_type::scalar_type;
57 using typename base_type::variable_type;
58
59 static_assert(!problem_type::allowed_evaluations.mass,
60 "Mass matrix is not supported.");
61
64
65protected:
66 using base_type::coeff;
68
69public:
71 static constexpr index_type stages = 1;
72
74 static constexpr index_type order = 2;
75
77 static constexpr auto log_tag = logging::log_tag_view(
78 "num_collect::ode::runge_kutta::implicit_euler_formula");
79
81 void step(scalar_type time, scalar_type step_size,
82 const variable_type& current, variable_type& estimate) {
83 formula_solver().update_jacobian(problem(), time + step_size, step_size,
84 current, static_cast<scalar_type>(1));
85 slope_ = problem().diff_coeff();
86 formula_solver().init(slope_);
87 formula_solver().solve();
88 estimate = current + step_size * slope_;
89 }
90
91private:
94};
95
101template <concepts::problem Problem>
103
104} // namespace num_collect::ode::runge_kutta
Class of tags of logs without memory management.
Base class of formulas in ODE solvers.
auto problem() -> problem_type &
Get the problem.
static constexpr auto coeff(T val) -> scalar_type
Convert coefficients.
static constexpr index_type stages
Number of stages of this formula.
auto formula_solver() -> formula_solver_type &
Get solver of formula.
static constexpr index_type order
Order of this formula.
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.
typename problem_type::scalar_type scalar_type
Type of scalars.
auto formula_solver() -> formula_solver_type &
Get solver of formula.
Class of simple solver of ODEs.
Definition of implicit_formula_base class.
Definition of index_type type.
Definition of inexact_newton_slope_equation_solver class.
Definition of log_tag_view class.
std::ptrdiff_t index_type
Type of indices in this library.
Definition index_type.h:33
Namespace of Runge-Kutta method.
Definition of one.
Definition of problem concept.
Definition of simple_solver class.
Definition of slope_equation_solver concept.