Reference documentation for deal.II version GIT relicensing-487-ge9eb5ab491 2024-04-25 07:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
nox.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2022 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_trilinos_nox
16#define dealii_trilinos_nox
17
18#include <deal.II/base/config.h>
19
20#ifdef DEAL_II_TRILINOS_WITH_NOX
21
23
25
26# include <Teuchos_ParameterList.hpp>
27
28# include <exception>
29# include <functional>
30
32
33namespace TrilinosWrappers
34{
35 // Indicate that NOXSolver has not converged.
37
38
86 template <typename VectorType>
88 {
89 public:
96 {
97 public:
101 AdditionalData(const unsigned int max_iter = 10,
102 const double abs_tol = 1.e-20,
103 const double rel_tol = 1.e-5,
104 const unsigned int threshold_nonlinear_iterations = 1,
105 const unsigned int threshold_n_linear_iterations = 0,
106 const bool reuse_solver = false);
107
111 unsigned int max_iter;
112
119 double abs_tol;
120
127 double rel_tol;
128
134
150
158 };
159
171 const Teuchos::RCP<Teuchos::ParameterList> &parameters =
172 Teuchos::rcp(new Teuchos::ParameterList));
173
178
182 void
184
189 unsigned int
190 solve(VectorType &solution);
191
204 std::function<void(const VectorType &u, VectorType &F)> residual;
205
218 std::function<void(const VectorType &current_u)> setup_jacobian;
219
237 std::function<void(const VectorType &current_u)> setup_preconditioner;
238
259 std::function<void(const VectorType &x, VectorType &y)> apply_jacobian;
260
284 std::function<
285 void(const VectorType &y, VectorType &x, const double tolerance)>
287
320 std::function<
321 int(const VectorType &y, VectorType &x, const double tolerance)>
323
343 std::function<SolverControl::State(const unsigned int i,
344 const double norm_f,
345 const VectorType &current_u,
346 const VectorType &f)>
348
368
369 private:
374
380 const Teuchos::RCP<Teuchos::ParameterList> parameters;
381
386
391
396
401
407 mutable std::exception_ptr pending_exception;
408 };
409} // namespace TrilinosWrappers
410
412
413#endif
414
415#endif
std::function< void(const VectorType &x, VectorType &y)> apply_jacobian
Definition nox.h:259
std::function< SolverControl::State(const unsigned int i, const double norm_f, const VectorType &current_u, const VectorType &f)> check_iteration_status
Definition nox.h:347
unsigned int n_residual_evaluations
Definition nox.h:385
unsigned int n_nonlinear_iterations
Definition nox.h:395
std::function< void(const VectorType &current_u)> setup_preconditioner
Definition nox.h:237
unsigned int n_last_linear_iterations
Definition nox.h:400
std::function< void(const VectorType &y, VectorType &x, const double tolerance)> solve_with_jacobian
Definition nox.h:286
std::function< bool()> update_preconditioner_predicate
Definition nox.h:367
std::function< void(const VectorType &u, VectorType &F)> residual
Definition nox.h:204
AdditionalData additional_data
Definition nox.h:373
unsigned int solve(VectorType &solution)
NOXSolver(AdditionalData &additional_data, const Teuchos::RCP< Teuchos::ParameterList > &parameters=Teuchos::rcp(new Teuchos::ParameterList))
std::function< void(const VectorType &current_u)> setup_jacobian
Definition nox.h:218
const Teuchos::RCP< Teuchos::ParameterList > parameters
Definition nox.h:380
std::exception_ptr pending_exception
Definition nox.h:407
unsigned int n_jacobian_applications
Definition nox.h:390
std::function< int(const VectorType &y, VectorType &x, const double tolerance)> solve_with_jacobian_and_track_n_linear_iterations
Definition nox.h:322
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define DeclException0(Exception0)
Definition exceptions.h:471
static ::ExceptionBase & ExcNOXNoConvergence()
AdditionalData(const unsigned int max_iter=10, const double abs_tol=1.e-20, const double rel_tol=1.e-5, const unsigned int threshold_nonlinear_iterations=1, const unsigned int threshold_n_linear_iterations=0, const bool reuse_solver=false)