numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
inexact_newton_slope_equation_solver.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
22#include <cmath>
23#include <optional>
24
25#include <Eigen/Core>
26#include <Eigen/LU>
27
40
41namespace num_collect::ode {
42
45 "num_collect::ode::inexact_newton_slope_equation_solver");
46
61template <concepts::differentiable_problem Problem>
63
78template <concepts::single_variate_differentiable_problem Problem>
80 : public iterative_solver_base<
81 inexact_newton_slope_equation_solver<Problem>> {
82public:
85
87 using problem_type = Problem;
88
90 using variable_type = typename problem_type::variable_type;
91
93 using scalar_type = typename problem_type::scalar_type;
94
96 using jacobian_type = typename problem_type::jacobian_type;
97
98 static_assert(!problem_type::allowed_evaluations.mass,
99 "Mass matrix is not supported.");
100 // TODO: Actually, support is not difficult, but I want to test after
101 // implementation, so I postpone the implementation.
102
107
119 scalar_type step_size, const variable_type& variable,
120 scalar_type solution_coeff) {
121 problem_ = &problem;
122 time_ = time;
123 step_size_ = step_size;
124 variable_ = variable;
125 solution_coeff_ = solution_coeff;
126
127 problem_->evaluate_on(time_, variable_,
128 evaluation_type{.diff_coeff = true, .jacobian = true});
129
130 coeff_inverse_ = static_cast<jacobian_type>(1) /
131 (static_cast<jacobian_type>(1) -
132 step_size_ * solution_coeff_ * problem_->jacobian());
133 using std::isfinite;
134 if (!isfinite(coeff_inverse_)) {
136 algorithm_failure, "Failed to calculate inverse.");
137 }
138 }
139
145 void init(variable_type& solution) {
146 solution_ = &solution;
147 update_norm_.reset();
148 if (update_reduction_rate_) {
149 constexpr auto exponent = static_cast<scalar_type>(0.8);
150 constexpr auto min_rate = static_cast<scalar_type>(0.5);
151 using std::pow;
152 *update_reduction_rate_ = pow(*update_reduction_rate_, exponent);
153 if (*update_reduction_rate_ < min_rate) {
154 *update_reduction_rate_ = min_rate;
155 }
156 }
157 iterations_ = 0;
158 }
159
166 void iterate() {
167 NUM_COLLECT_PRECONDITION(problem_ != nullptr && solution_ != nullptr,
168 this->logger(), "Initialization must be done before iterations.");
169
170 temp_variable_ =
171 variable_ + step_size_ * solution_coeff_ * (*solution_);
172
173 problem_->evaluate_on(
174 time_, temp_variable_, evaluation_type{.diff_coeff = true});
175 residual_ = (*solution_) - problem_->diff_coeff();
176 update_ = -coeff_inverse_ * residual_;
177 *solution_ += update_;
178
179 const scalar_type update_norm =
180 tolerances().calc_norm(variable_, update_);
181 if (update_norm_) {
182 update_reduction_rate_ = update_norm / (*update_norm_);
183 }
184 update_norm_ = update_norm;
185
186 ++iterations_;
187 }
188
194 [[nodiscard]] auto is_stop_criteria_satisfied() const -> bool {
195 bool converged = false;
196 if (update_norm_ && update_reduction_rate_ &&
197 *update_reduction_rate_ < static_cast<scalar_type>(1)) {
198 converged =
199 (*update_reduction_rate_ /
200 (static_cast<scalar_type>(1) - *update_reduction_rate_)) *
201 (*update_norm_) <=
202 tolerance_rate_;
203 }
204
205 constexpr index_type max_iterations = 1000; // safe guard
206 return converged || (iterations_ > max_iterations);
207 }
208
216 const {
217 iteration_logger.template append<index_type>(
218 "Iter.", &this_type::iterations);
219 iteration_logger.template append<scalar_type>(
220 "Update", &this_type::update_norm);
221 }
222
228 [[nodiscard]] auto update_norm() const -> scalar_type {
229 if (!update_norm_) {
230 return static_cast<scalar_type>(0);
231 }
232 return *update_norm_;
233 }
234
240 [[nodiscard]] auto iterations() const -> index_type { return iterations_; }
241
250 tolerances_ = val;
251 return *this;
252 }
253
259 [[nodiscard]] auto tolerances() const
261 return tolerances_;
262 }
263
264private:
266 problem_type* problem_{nullptr};
267
269 scalar_type time_{};
270
272 scalar_type step_size_{};
273
275 scalar_type solution_coeff_{};
276
278 variable_type variable_{};
279
281 variable_type* solution_{nullptr};
282
284 variable_type coeff_inverse_{};
285
287 variable_type temp_variable_{};
288
290 variable_type residual_{};
291
293 variable_type update_{};
294
296 std::optional<scalar_type> update_norm_{};
297
299 std::optional<scalar_type> update_reduction_rate_{};
300
302 static constexpr auto default_tolerance_rate =
303 static_cast<scalar_type>(1e-2);
304
306 scalar_type tolerance_rate_{default_tolerance_rate};
307
309 index_type iterations_{0};
310
313};
314
329template <concepts::multi_variate_differentiable_problem Problem>
331 : public iterative_solver_base<
332 inexact_newton_slope_equation_solver<Problem>> {
333public:
336
338 using problem_type = Problem;
339
341 using variable_type = typename problem_type::variable_type;
342
344 using scalar_type = typename problem_type::scalar_type;
345
347 using jacobian_type = typename problem_type::jacobian_type;
348
349 static_assert(!problem_type::allowed_evaluations.mass,
350 "Mass matrix is not supported.");
351 // TODO: Actually, support is not difficult, but I want to test after
352 // implementation, so I postpone the implementation.
353
358
370 template <typename VariableExpression>
372 scalar_type step_size,
373 const Eigen::MatrixBase<VariableExpression>& variable,
374 scalar_type solution_coeff) {
375 problem_ = &problem;
376 time_ = time;
377 step_size_ = step_size;
378 variable_ = variable;
379 solution_coeff_ = solution_coeff;
380
381 problem_->evaluate_on(time_, variable_,
382 evaluation_type{.diff_coeff = true, .jacobian = true});
383
384 const index_type dim = variable_.size();
385 lu_.compute(jacobian_type::Identity(dim, dim) -
386 step_size_ * solution_coeff_ * problem_->jacobian());
387 }
388
394 void init(variable_type& solution) {
395 solution_ = &solution;
396 update_norm_.reset();
397 if (update_reduction_rate_) {
398 constexpr auto exponent = static_cast<scalar_type>(0.8);
399 constexpr auto min_rate = static_cast<scalar_type>(0.5);
400 using std::pow;
401 *update_reduction_rate_ = pow(*update_reduction_rate_, exponent);
402 if (*update_reduction_rate_ < min_rate) {
403 *update_reduction_rate_ = min_rate;
404 }
405 }
406 iterations_ = 0;
407 }
408
415 void iterate() {
416 NUM_COLLECT_PRECONDITION(problem_ != nullptr && solution_ != nullptr,
417 this->logger(), "Initialization must be done before iterations.");
418
419 temp_variable_ =
420 variable_ + step_size_ * solution_coeff_ * (*solution_);
421
422 problem_->evaluate_on(
423 time_, temp_variable_, evaluation_type{.diff_coeff = true});
424 residual_ = (*solution_) - problem_->diff_coeff();
425 update_ = -lu_.solve(residual_);
426 if (!update_.array().isFinite().all()) {
428 "Failed to solve an equation. step_size={}, cond={}.",
429 step_size_, lu_.rcond());
430 }
431 *solution_ += update_;
432
433 const scalar_type update_norm =
434 tolerances().calc_norm(variable_, update_);
435 if (update_norm_) {
436 update_reduction_rate_ = update_norm / (*update_norm_);
437 }
438 update_norm_ = update_norm;
439
440 ++iterations_;
441 }
442
448 [[nodiscard]] auto is_stop_criteria_satisfied() const -> bool {
449 bool converged = false;
450 if (update_norm_ && update_reduction_rate_ &&
451 *update_reduction_rate_ < static_cast<scalar_type>(1)) {
452 converged =
453 (*update_reduction_rate_ /
454 (static_cast<scalar_type>(1) - *update_reduction_rate_)) *
455 (*update_norm_) <=
456 tolerance_rate_;
457 }
458
459 constexpr index_type max_iterations = 100; // safe guard
460 return converged || (iterations_ > max_iterations);
461 }
462
470 const {
471 iteration_logger.template append<index_type>(
472 "Iter.", &this_type::iterations);
473 iteration_logger.template append<scalar_type>(
474 "Update", &this_type::update_norm);
475 }
476
482 [[nodiscard]] auto update_norm() const -> scalar_type {
483 if (!update_norm_) {
484 return static_cast<scalar_type>(0);
485 }
486 return *update_norm_;
487 }
488
494 [[nodiscard]] auto iterations() const -> index_type { return iterations_; }
495
504 tolerances_ = val;
505 return *this;
506 }
507
513 [[nodiscard]] auto tolerances() const
515 return tolerances_;
516 }
517
518private:
520 problem_type* problem_{nullptr};
521
523 scalar_type time_{};
524
526 scalar_type step_size_{};
527
529 scalar_type solution_coeff_{};
530
532 variable_type variable_{};
533
535 variable_type* solution_{nullptr};
536
538 Eigen::PartialPivLU<jacobian_type> lu_{};
539
541 variable_type temp_variable_{};
542
544 variable_type residual_{};
545
547 variable_type update_{};
548
550 std::optional<scalar_type> update_norm_{};
551
553 std::optional<scalar_type> update_reduction_rate_{};
554
556 static constexpr auto default_tolerance_rate =
557 static_cast<scalar_type>(1e-2);
558
560 scalar_type tolerance_rate_{default_tolerance_rate};
561
563 index_type iterations_{0};
564
566 error_tolerances<variable_type> tolerances_{};
567};
568
569} // namespace num_collect::ode
Definition of iterative_solver_base class.
Class of exception on failure in algorithm.
Definition exception.h:93
Base class of iterative solvers.
Class of tags of logs without memory management.
Class of error tolerances hairer1993.
Class to solve equations of implicit slopes using inexact Newton method.
auto tolerances(const error_tolerances< variable_type > &val) -> inexact_newton_slope_equation_solver &
Set the error tolerances.
auto is_stop_criteria_satisfied() const -> bool
Determine if stopping criteria of the algorithm are satisfied.
void configure_iteration_logger(logging::iterations::iteration_logger< this_type > &iteration_logger) const
Configure an iteration logger.
void update_jacobian(problem_type &problem, scalar_type time, scalar_type step_size, const Eigen::MatrixBase< VariableExpression > &variable, scalar_type solution_coeff)
Update Jacobian and internal parameters.
void update_jacobian(problem_type &problem, scalar_type time, scalar_type step_size, const variable_type &variable, scalar_type solution_coeff)
Update Jacobian and internal parameters.
void init(variable_type &solution)
Initialize for solving an equation.
auto tolerances() const -> const error_tolerances< variable_type > &
Get the error tolerances.
Class to solve equations of implicit slopes using inexact Newton method.
Definition of differentiable_problem concept.
Definition of error_tolerances class.
Definition of evaluation_type enumeration.
Definition of exceptions.
Definition of index_type type.
Definition of iteration_logger class.
Definition of log_tag_view class.
Definition of macros for logging.
#define NUM_COLLECT_LOG_AND_THROW(EXCEPTION_TYPE,...)
Write an error log and throw an exception for an error.
Definition of multi_variate_differentiable_problem concept.
std::ptrdiff_t index_type
Type of indices in this library.
Definition index_type.h:33
auto isfinite(const T &val) -> bool
Check whether a number is finite.
Definition isfinite.h:39
Namespace of solvers of ordinary differential equations (ODE).
constexpr auto inexact_newton_slope_equation_solver_tag
Log tag.
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 single_variate_differentiable_problem concept.
Struct to specify types of evaluations.