File snes.hxx

Functions

BOUT_ENUM_CLASS(BoutSnesEquationForm, pseudo_transient, rearranged_backward_euler, backward_euler, direct_newton)
class SNESSolver : public Solver
#include <snes.hxx>

Uses PETSc’s SNES interface to find a steady state solution to a nonlinear ODE by integrating in time with Backward Euler

Public Functions

explicit SNESSolver(Options *opt = nullptr)
inline ~SNESSolver()
virtual int init(int nout, BoutReal tstep) override

Initialise solver. Must be called once and only once

Parameters
  • nout[in] Number of outputs

  • tstep[in] Time between outputs. NB: Not internal timestep

virtual int run() override

Run the simulation.

PetscErrorCode snes_function(Vec x, Vec f, bool linear)

Nonlinear function.

Nonlinear function. This is called by PETSc SNES object via a static C-style function. For implicit time integration this function calculates:

f = (x - gamma*G(x)) - rhs

Parameters
  • x[in] The state vector

  • f[out] The vector for the result f(x)

  • linear[in] Specifies that the SNES solver is in a linear (KSP) inner loop, so the operator should be linearised if possible

PetscErrorCode precon(Vec x, Vec f)

Preconditioner. Called by PCapply via a C-style static function.

Parameters
  • x[in] The vector to be operated on

  • f[out] The result of the operation

Private Members

BoutReal timestep

Internal timestep.

BoutReal dt

Current timestep used in snes_function.

BoutReal dt_min_reset

If dt falls below this, reset solve.

BoutReal max_timestep

Maximum timestep.

int lower_its
int upper_its

Limits on iterations for timestep adjustment.

BoutReal out_timestep

Output timestep.

int nsteps

Number of steps to take.

bool diagnose

Output additional diagnostics.

bool diagnose_failures

Print diagnostics on SNES failures.

int nlocal

Number of variables on local processor.

int neq

Number of variables in total.

BoutSnesEquationForm equation_form

Form of the equation to solve.

PetscLib lib

Handles initialising, finalising PETSc.

Vec snes_f

Used by SNES to store function.

Vec snes_x

Result of SNES.

Vec x0

Solution at start of current timestep.

Vec delta_x

Change in solution.

bool predictor

Use linear predictor?

Vec x1

Previous solution.

BoutReal time1 = {-1.0}

Time of previous solution.

SNES snes

SNES context.

Mat Jmf

Matrix-free Jacobian.

MatFDColoring fdcoloring

Matrix coloring context, used for finite difference Jacobian evaluation