OGS
NumLib::PETScNonlinearSolver Class Referencefinal

Detailed Description

Definition at line 17 of file PETScNonlinearSolver.h.

#include <PETScNonlinearSolver.h>

Inheritance diagram for NumLib::PETScNonlinearSolver:
[legend]
Collaboration diagram for NumLib::PETScNonlinearSolver:
[legend]

Public Types

using System = NonlinearSystem<NonlinearSolverTag::Newton>
 Type of the nonlinear equation system to be solved.

Public Member Functions

 PETScNonlinearSolver (GlobalLinearSolver &linear_solver, int maxiter, std::string prefix)
void setEquationSystem (System &eq, ConvergenceCriterion &conv_crit)
 Set the nonlinear equation system that will be solved.
void compensateNonEquilibriumInitialResiduum (bool const value)
void calculateNonEquilibriumInitialResiduum (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) override
NonlinearSolverStatus solve (std::vector< GlobalVector * > &x, std::vector< GlobalVector * > const &x_prev, std::function< void(int, bool, std::vector< GlobalVector * > const &)> const &postIterationCallback, int const process_id) override
Public Member Functions inherited from NumLib::NonlinearSolverBase
virtual ~NonlinearSolverBase ()=default

Private Attributes

SNES _snes_solver
System_equation_system = nullptr
ConvergenceCriterion_convergence_criterion = nullptr
std::size_t _residual_id = 0u
 ID of the residual vector.
std::size_t _jacobian_id = 0u
 ID of the Jacobian matrix.
std::size_t _petsc_x_id = 0u
std::size_t _petsc_jacobian_id = 0u
std::size_t _petsc_residual_id = 0u
bool _compensate_non_equilibrium_initial_residuum = false

Member Typedef Documentation

◆ System

Type of the nonlinear equation system to be solved.

Definition at line 21 of file PETScNonlinearSolver.h.

Constructor & Destructor Documentation

◆ PETScNonlinearSolver()

NumLib::PETScNonlinearSolver::PETScNonlinearSolver ( GlobalLinearSolver & linear_solver,
int maxiter,
std::string prefix )

Constructs a new instance.

Parameters
linear_solverthe linear solver used by this nonlinear solver.
maxiterthe maximum number of iterations used to solve the equation.
prefixused to set the SNES options.

Definition at line 88 of file PETScNonlinearSolver.cpp.

91{
92 SNESCreate(PETSC_COMM_WORLD, &_snes_solver);
93 if (!prefix.empty())
94 {
95 prefix = prefix + "_";
96 SNESSetOptionsPrefix(_snes_solver, prefix.c_str());
97 }
98 // force SNESSolve() to take at least one iteration regardless of the
99 // initial residual norm
100 SNESSetForceIteration(_snes_solver, PETSC_TRUE);
101
102 // Set the maximum iterations.
103 PetscReal atol, rtol, stol;
104 PetscInt maxf;
105 SNESGetTolerances(_snes_solver, &atol, &rtol, &stol, nullptr, &maxf);
106 SNESSetTolerances(_snes_solver, atol, rtol, stol, maxiter, maxf);
107
108 SNESSetFromOptions(_snes_solver);
109#ifndef NDEBUG
110 PetscOptionsView(nullptr, PETSC_VIEWER_STDOUT_WORLD);
111#endif // NDEBUG
112}

References _snes_solver.

Member Function Documentation

◆ calculateNonEquilibriumInitialResiduum()

void NumLib::PETScNonlinearSolver::calculateNonEquilibriumInitialResiduum ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id )
overridevirtual

Implements NumLib::NonlinearSolverBase.

Definition at line 127 of file PETScNonlinearSolver.cpp.

130{
132 {
133 return;
134 }
135
136 OGS_FATAL(
137 "Non-equilibrium initial residuum is not implemented for the "
138 "PETScNonlinearSolver.");
139}
#define OGS_FATAL(...)
Definition Error.h:19

References _compensate_non_equilibrium_initial_residuum, and OGS_FATAL.

◆ compensateNonEquilibriumInitialResiduum()

void NumLib::PETScNonlinearSolver::compensateNonEquilibriumInitialResiduum ( bool const value)

◆ setEquationSystem()

void NumLib::PETScNonlinearSolver::setEquationSystem ( System & eq,
ConvergenceCriterion & conv_crit )

Set the nonlinear equation system that will be solved.

Definition at line 114 of file PETScNonlinearSolver.cpp.

116{
117 _equation_system = &eq;
118 _convergence_criterion = &conv_crit;
119}
ConvergenceCriterion * _convergence_criterion

References _convergence_criterion, and _equation_system.

◆ solve()

NonlinearSolverStatus NumLib::PETScNonlinearSolver::solve ( std::vector< GlobalVector * > & x,
std::vector< GlobalVector * > const & x_prev,
std::function< void(int, bool, std::vector< GlobalVector * > const &)> const & postIterationCallback,
int const process_id )
overridevirtual

Assemble and solve the equation system.

Parameters
xin: the initial guess, out: the solution.
x_prevprevious time step solution.
postIterationCallbackcalled after each iteration if set.
process_idusually used in staggered schemes.
Return values
trueif the equation system could be solved
falseotherwise

Implements NumLib::NonlinearSolverBase.

Definition at line 141 of file PETScNonlinearSolver.cpp.

148{
149 DBUG("PETScNonlinearSolver: solve()");
150 using TimeDiscretizedSystem =
151 TimeDiscretizedODESystem<ODESystemTag::FirstOrderImplicitQuasilinear,
153
154 auto* system = static_cast<TimeDiscretizedSystem*>(_equation_system);
155
156 DBUG("PETScNonlinearSolver: create vectors");
157 // r and J on which the ogs assembly operates.
159 system->getMatrixSpecifications(process_id), _jacobian_id);
161 system->getMatrixSpecifications(process_id), _residual_id);
162
163 // Vectors and matrices used by SNES for solutions. These will be locked
164 // during the SNES' solve call.
166 system->getMatrixSpecifications(process_id), _petsc_jacobian_id);
167 BaseLib::RunTime timer_dirichlet;
168 double time_dirichlet = 0.0;
169
170 timer_dirichlet.start();
171 system->computeKnownSolutions(*x[process_id], process_id);
172 system->applyKnownSolutions(*x[process_id]);
173 time_dirichlet += timer_dirichlet.elapsed();
174 INFO("[time] Applying Dirichlet BCs took {} s.", time_dirichlet);
175
177 system->getMatrixSpecifications(process_id), _petsc_residual_id);
179 system->getMatrixSpecifications(process_id), _petsc_x_id);
180 MathLib::LinAlg::copy(*x[process_id], x_snes); // Initial guess.
181
182 PetscContext petsc_context{_equation_system, x, x_prev, &r, &J, process_id};
183
184 DBUG("PETScNonlinearSolver: set function");
185 SNESSetFunction(_snes_solver, r_snes.getRawVector(), updateResidual,
186 &petsc_context);
187
188 DBUG("PETScNonlinearSolver: set jacobian");
189 // The jacobian and the preconditioner matrices point to the same location.
190 SNESSetJacobian(_snes_solver, J_snes.getRawMatrix(), J_snes.getRawMatrix(),
191 updateJacobian, &petsc_context);
192
193 std::unique_ptr<GlobalVector> xl = nullptr;
194 std::unique_ptr<GlobalVector> xu = nullptr;
195
196 SNESType snes_type;
197 SNESGetType(_snes_solver, &snes_type);
198 if ((std::strcmp(snes_type, SNESVINEWTONRSLS) == 0) ||
199 (std::strcmp(snes_type, SNESVINEWTONSSLS) == 0))
200 {
201 // Set optional constraints via callback.
202 DBUG("PETScNonlinearSolver: set constraints");
203 xl = MathLib::MatrixVectorTraits<GlobalVector>::newInstance(
204 system->getMatrixSpecifications(process_id));
205 xu = MathLib::MatrixVectorTraits<GlobalVector>::newInstance(
206 system->getMatrixSpecifications(process_id));
207
208 system->updateConstraints(*xl, *xu, process_id);
211
212 SNESVISetVariableBounds(_snes_solver, xl->getRawVector(),
213 xu->getRawVector());
214 }
215
216 DBUG("PETScNonlinearSolver: call SNESSolve");
217 SNESSolve(_snes_solver, nullptr, x_snes.getRawVector());
218
219 SNESConvergedReason reason;
220 SNESGetConvergedReason(_snes_solver, &reason);
222
223 PetscInt iterations;
224 SNESGetIterationNumber(_snes_solver, &iterations);
225 INFO("PETScSNES used {} iterations.", iterations);
226
227 // Copy back the solution.
228 MathLib::LinAlg::copy(x_snes, *x[process_id]);
229
235
236 return {reason >= 0, iterations};
237}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:31
void start()
Start the timer.
Definition RunTime.h:21
virtual void releaseMatrix(GlobalMatrix const &A)=0
virtual GlobalMatrix & getMatrix(std::size_t &id)=0
Get an uninitialized matrix with the given id.
std::size_t _jacobian_id
ID of the Jacobian matrix.
std::size_t _residual_id
ID of the residual vector.
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
virtual void releaseVector(GlobalVector const &x)=0
@ FirstOrderImplicitQuasilinear
Definition Types.h:27
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:30
void finalizeVectorAssembly(VEC_T &)
General function to finalize the vector assembly.
void printConvergenceReason(SNES const &snes_solver, SNESConvergedReason const &reason)
static NUMLIB_EXPORT MatrixProvider & provider
static NUMLIB_EXPORT VectorProvider & provider

References _equation_system, _jacobian_id, _petsc_jacobian_id, _petsc_residual_id, _petsc_x_id, _residual_id, _snes_solver, MathLib::LinAlg::copy(), DBUG(), BaseLib::RunTime::elapsed(), MathLib::finalizeVectorAssembly(), NumLib::FirstOrderImplicitQuasilinear, INFO(), NumLib::Newton, NumLib::GlobalMatrixProvider::provider, NumLib::GlobalVectorProvider::provider, and BaseLib::RunTime::start().

Member Data Documentation

◆ _compensate_non_equilibrium_initial_residuum

bool NumLib::PETScNonlinearSolver::_compensate_non_equilibrium_initial_residuum = false
private

Enables computation of the non-equilibrium initial residuum \( r_{\rm neq} \) before the first time step. The forces are zero if the external forces are in equilibrium with the initial state/initial conditions. During the simulation the new residuum reads \( \tilde r = r - r_{\rm neq} \).

Definition at line 66 of file PETScNonlinearSolver.h.

Referenced by calculateNonEquilibriumInitialResiduum(), and compensateNonEquilibriumInitialResiduum().

◆ _convergence_criterion

ConvergenceCriterion* NumLib::PETScNonlinearSolver::_convergence_criterion = nullptr
private

Definition at line 55 of file PETScNonlinearSolver.h.

Referenced by setEquationSystem().

◆ _equation_system

System* NumLib::PETScNonlinearSolver::_equation_system = nullptr
private

Definition at line 54 of file PETScNonlinearSolver.h.

Referenced by setEquationSystem(), and solve().

◆ _jacobian_id

std::size_t NumLib::PETScNonlinearSolver::_jacobian_id = 0u
private

ID of the Jacobian matrix.

Definition at line 58 of file PETScNonlinearSolver.h.

Referenced by solve().

◆ _petsc_jacobian_id

std::size_t NumLib::PETScNonlinearSolver::_petsc_jacobian_id = 0u
private

Definition at line 61 of file PETScNonlinearSolver.h.

Referenced by solve().

◆ _petsc_residual_id

std::size_t NumLib::PETScNonlinearSolver::_petsc_residual_id = 0u
private

Definition at line 62 of file PETScNonlinearSolver.h.

Referenced by solve().

◆ _petsc_x_id

std::size_t NumLib::PETScNonlinearSolver::_petsc_x_id = 0u
private

Definition at line 60 of file PETScNonlinearSolver.h.

Referenced by solve().

◆ _residual_id

std::size_t NumLib::PETScNonlinearSolver::_residual_id = 0u
private

ID of the residual vector.

Definition at line 57 of file PETScNonlinearSolver.h.

Referenced by solve().

◆ _snes_solver

SNES NumLib::PETScNonlinearSolver::_snes_solver
private

Definition at line 52 of file PETScNonlinearSolver.h.

Referenced by PETScNonlinearSolver(), and solve().


The documentation for this class was generated from the following files: