numerical-collection-cpp 0.10.0
A collection of algorithms in numerical analysis implemented in C++
Loading...
Searching...
No Matches
inexact_newton_update_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
46 "num_collect::ode::inexact_newton_update_equation_solver");
47
64template <concepts::differentiable_problem Problem>
66
83template <concepts::single_variate_differentiable_problem Problem>
85 : public iterative_solver_base<
86 inexact_newton_update_equation_solver<Problem>> {
87public:
90
92 using problem_type = Problem;
93
95 using variable_type = typename problem_type::variable_type;
96
98 using scalar_type = typename problem_type::scalar_type;
99
101 using jacobian_type = typename problem_type::jacobian_type;
102
103 static_assert(!problem_type::allowed_evaluations.mass,
104 "Mass matrix is not supported.");
105 // TODO: Actually, support is not difficult, but I want to test after
106 // implementation, so I postpone the implementation.
107
112
124 scalar_type step_size, const variable_type& variable,
125 scalar_type slope_coeff) {
126 problem_ = &problem;
127 time_ = time;
128 step_size_ = step_size;
129 variable_ = variable;
130 slope_coeff_ = slope_coeff;
131
132 problem_->evaluate_on(time_, variable_,
133 evaluation_type{.diff_coeff = true, .jacobian = true});
134
135 coeff_inverse_ = static_cast<jacobian_type>(1) /
136 (static_cast<jacobian_type>(1) -
137 step_size_ * slope_coeff_ * problem_->jacobian());
138 using std::isfinite;
139 if (!isfinite(coeff_inverse_)) {
141 algorithm_failure, "Failed to calculate inverse.");
142 }
143 }
144
154 solution_ = &solution;
155 update_norm_.reset();
157 constexpr auto exponent = static_cast<scalar_type>(0.8);
158 constexpr auto min_rate = static_cast<scalar_type>(0.5);
159 using std::pow;
161 if (*update_reduction_rate_ < min_rate) {
162 *update_reduction_rate_ = min_rate;
163 }
164 }
165 iterations_ = 0;
166 }
167
177 variable_type& solution) {
178 time_ = time;
179 init(solution_offset, solution);
180 }
181
188 void iterate() {
189 NUM_COLLECT_PRECONDITION(problem_ != nullptr && solution_ != nullptr,
190 this->logger(), "Initialization must be done before iterations.");
191
192 temp_variable_ = variable_ + (*solution_);
193
194 problem_->evaluate_on(
195 time_, temp_variable_, evaluation_type{.diff_coeff = true});
196 residual_ = (*solution_) -
197 step_size_ * slope_coeff_ * problem_->diff_coeff() -
200 *solution_ += update_;
201
203 tolerances().calc_norm(variable_, update_);
204 if (update_norm_) {
205 update_reduction_rate_ = update_norm / (*update_norm_);
206 }
208
209 ++iterations_;
210 }
211
217 [[nodiscard]] auto is_stop_criteria_satisfied() const -> bool {
218 bool converged = false;
220 *update_reduction_rate_ < static_cast<scalar_type>(1)) {
221 converged =
223 (static_cast<scalar_type>(1) - *update_reduction_rate_)) *
224 (*update_norm_) <=
226 }
227
228 constexpr index_type max_iterations = 1000; // safe guard
229 return converged || (iterations_ > max_iterations);
230 }
231
239 const {
240 iteration_logger.template append<index_type>(
241 "Iter.", &this_type::iterations);
242 iteration_logger.template append<scalar_type>(
243 "Update", &this_type::update_norm);
244 }
245
251 [[nodiscard]] auto solution_offset() const -> const variable_type& {
252 return solution_offset_;
253 }
254
260 [[nodiscard]] auto update_norm() const -> scalar_type {
261 if (!update_norm_) {
262 return static_cast<scalar_type>(0);
263 }
264 return *update_norm_;
265 }
266
272 [[nodiscard]] auto iterations() const -> index_type { return iterations_; }
273
282 tolerances_ = val;
283 return *this;
284 }
285
291 [[nodiscard]] auto tolerances() const
293 return tolerances_;
294 }
295
296private:
299
302
305
308
311
314
317
320
323
326
329
331 std::optional<scalar_type> update_norm_{};
332
334 std::optional<scalar_type> update_reduction_rate_{};
335
337 static constexpr auto default_tolerance_rate =
338 static_cast<scalar_type>(1e-2);
339
342
345
348};
349
366template <concepts::multi_variate_differentiable_problem Problem>
367class inexact_newton_update_equation_solver<Problem>
368 : public iterative_solver_base<
369 inexact_newton_update_equation_solver<Problem>> {
370public:
373
375 using problem_type = Problem;
376
378 using variable_type = typename problem_type::variable_type;
379
381 using scalar_type = typename problem_type::scalar_type;
382
384 using jacobian_type = typename problem_type::jacobian_type;
385
386 static_assert(!problem_type::allowed_evaluations.mass,
387 "Mass matrix is not supported.");
388 // TODO: Actually, support is not difficult, but I want to test after
389 // implementation, so I postpone the implementation.
390
395
407 template <typename VariableExpression>
409 scalar_type step_size,
410 const Eigen::MatrixBase<VariableExpression>& variable,
411 scalar_type slope_coeff) {
412 problem_ = &problem;
413 time_ = time;
414 step_size_ = step_size;
415 variable_ = variable;
416 slope_coeff_ = slope_coeff;
417
418 problem_->evaluate_on(time_, variable_,
419 evaluation_type{.diff_coeff = true, .jacobian = true});
420
421 const index_type dim = variable_.size();
422 lu_.compute(jacobian_type::Identity(dim, dim) -
423 step_size_ * slope_coeff_ * problem_->jacobian());
424 }
425
435 template <typename OffsetExpression>
436 void init(const Eigen::MatrixBase<OffsetExpression>& solution_offset,
437 variable_type& solution) {
439 solution_ = &solution;
440 update_norm_.reset();
442 constexpr auto exponent = static_cast<scalar_type>(0.8);
443 constexpr auto min_rate = static_cast<scalar_type>(0.5);
444 using std::pow;
446 if (*update_reduction_rate_ < min_rate) {
447 *update_reduction_rate_ = min_rate;
448 }
449 }
450 iterations_ = 0;
451 }
452
463 template <typename OffsetExpression>
464 void init(const scalar_type& time,
465 const Eigen::MatrixBase<OffsetExpression>& solution_offset,
466 variable_type& solution) {
467 time_ = time;
468 init(solution_offset, solution);
469 }
470
477 void iterate() {
478 NUM_COLLECT_PRECONDITION(problem_ != nullptr && solution_ != nullptr,
479 this->logger(), "Initialization must be done before iterations.");
480
481 temp_variable_ = variable_ + (*solution_);
482
483 problem_->evaluate_on(
484 time_, temp_variable_, evaluation_type{.diff_coeff = true});
485 residual_ = (*solution_) -
486 step_size_ * slope_coeff_ * problem_->diff_coeff() -
488 update_ = -lu_.solve(residual_);
489 if (!update_.array().isFinite().all()) {
491 "Failed to solve an equation. step_size={}, cond={}.",
492 step_size_, lu_.rcond());
493 }
494 *solution_ += update_;
495
497 tolerances().calc_norm(variable_, update_);
498 if (update_norm_) {
499 update_reduction_rate_ = update_norm / (*update_norm_);
500 }
502
503 ++iterations_;
504 }
505
511 [[nodiscard]] auto is_stop_criteria_satisfied() const -> bool {
512 bool converged = false;
514 *update_reduction_rate_ < static_cast<scalar_type>(1)) {
515 converged =
517 (static_cast<scalar_type>(1) - *update_reduction_rate_)) *
518 (*update_norm_) <=
520 }
521
522 constexpr index_type max_iterations = 100; // safe guard
523 return converged || (iterations_ > max_iterations);
524 }
525
533 const {
534 iteration_logger.template append<index_type>(
535 "Iter.", &this_type::iterations);
536 iteration_logger.template append<scalar_type>(
537 "Update", &this_type::update_norm);
538 }
539
545 [[nodiscard]] auto solution_offset() const -> const variable_type& {
546 return solution_offset_;
547 }
548
554 [[nodiscard]] auto update_norm() const -> scalar_type {
555 if (!update_norm_) {
556 return static_cast<scalar_type>(0);
557 }
558 return *update_norm_;
559 }
560
566 [[nodiscard]] auto iterations() const -> index_type { return iterations_; }
567
576 tolerances_ = val;
577 return *this;
578 }
579
585 [[nodiscard]] auto tolerances() const
587 return tolerances_;
588 }
589
590private:
592 problem_type* problem_{nullptr};
593
595 scalar_type time_{};
596
598 scalar_type step_size_{};
599
601 scalar_type slope_coeff_{};
602
604 variable_type variable_{};
605
607 variable_type solution_offset_{};
608
610 variable_type* solution_{nullptr};
611
613 Eigen::PartialPivLU<jacobian_type> lu_{};
614
616 variable_type temp_variable_{};
617
619 variable_type residual_{};
620
622 variable_type update_{};
623
625 std::optional<scalar_type> update_norm_{};
626
628 std::optional<scalar_type> update_reduction_rate_{};
629
631 static constexpr auto default_tolerance_rate =
632 static_cast<scalar_type>(1e-2);
633
635 scalar_type tolerance_rate_{default_tolerance_rate};
636
638 index_type iterations_{0};
639
641 error_tolerances<variable_type> tolerances_{};
642};
643
644} // 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.
void update_jacobian(problem_type &problem, scalar_type time, scalar_type step_size, const Eigen::MatrixBase< VariableExpression > &variable, scalar_type slope_coeff)
Update Jacobian and internal parameters.
static constexpr auto default_tolerance_rate
Default rate of tolerance in this solver.
std::optional< scalar_type > update_reduction_rate_
Rate in which update is reduced from the previous step.
variable_type solution_offset_
Offset of the solution added to the term of slopes.
void update_jacobian(problem_type &problem, scalar_type time, scalar_type step_size, const variable_type &variable, scalar_type slope_coeff)
Update Jacobian and internal parameters.
Eigen::PartialPivLU< jacobian_type > lu_
LU decomposition of the current coefficient matrix.
auto solution_offset() const -> const variable_type &
Get the offset of the solution added to the term of slopes.
void init(const Eigen::MatrixBase< OffsetExpression > &solution_offset, variable_type &solution)
Initialize for solving an equation.
void configure_iteration_logger(logging::iterations::iteration_logger< this_type > &iteration_logger) const
Configure an iteration logger.
void configure_iteration_logger(logging::iterations::iteration_logger< inexact_newton_update_equation_solver< Problem > > &iteration_logger) const
Configure an iteration logger.
scalar_type slope_coeff_
Coefficient to multiply to slope in the equation.
auto tolerances(const error_tolerances< variable_type > &val) -> inexact_newton_update_equation_solver &
Set the error tolerances.
auto tolerances() const -> const error_tolerances< variable_type > &
Get the error tolerances.
void init(const scalar_type &time, const Eigen::MatrixBase< OffsetExpression > &solution_offset, variable_type &solution)
Initialize for solving an equation.
inexact_newton_update_equation_solver< Problem > this_type
This class.
void init(const variable_type &solution_offset, variable_type &solution)
Initialize for solving an equation.
void init(const scalar_type &time, const variable_type &solution_offset, variable_type &solution)
Initialize for solving an equation.
auto is_stop_criteria_satisfied() const -> bool
Determine if stopping criteria of the algorithm are satisfied.
Class to solve equations of implicit updates 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_update_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.