File snes.hxx#

Functions

BOUT_ENUM_CLASS(BoutSnesEquationForm, pseudo_transient, rearranged_backward_euler, backward_euler, direct_newton)#
BOUT_ENUM_CLASS(BoutPTCStrategy, inverse_residual, history_based, hybrid)#

Combine inverse_residual and history_based strategies.

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 *opts = nullptr)#
~SNESSolver() override = default#
virtual int init() override#

Initialise the solver.

virtual int run() override#

Run the solver, calling monitors nout times, at intervals of tstep. This function is called by solve(), and is specific to each solver type

This should probably be protected, since it shouldn’t be called by users.

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
The form depends on equation_form

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

PetscErrorCode scaleJacobian(Mat B)#

Scale an approximate Jacobian, and update the internal RHS scaling factors This is called by SNESComputeJacobianScaledColor with the finite difference approximated Jacobian.

virtual void outputVars(Options &output_options, bool save_repeat = true) override#

Save diagnostics to output.

Private Functions

PetscErrorCode FDJinitialise()#

Finite Difference Jacobian initialise.

PetscErrorCode FDJpruneJacobian()#

Remove small elements from the Jacobian.

PetscErrorCode FDJrestoreFromPruning()#

Restore Jacobian to original pattern.

PetscErrorCode rhs_function(Vec x, Vec f, bool linear)#

Call the physics model RHS function

Parameters:
  • x[in] The state vector. Will be scaled if scale_vars=true

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

  • linear[in] Specifies that the SNES solver is in a linear (KSP) inner loop

PetscErrorCode initPseudoTimestepping()#

Initialize the Pseudo-Transient Continuation method.

PetscErrorCode updatePseudoTimestepping(Vec x)#

Update dt_vec based on new solution x.

BoutReal updatePseudoTimestep(BoutReal previous_timestep, BoutReal previous_residual, BoutReal current_residual)#

Decide the next pseudo-timestep. Called by updatePseudoTimestepping.

Switched Evolution Relaxation (SER)

BoutReal updatePseudoTimestep_inverse_residual(BoutReal previous_timestep, BoutReal current_residual)#

Timestep inversely proportional to the residual, clipped to avoid rapid changes in timestep.

BoutReal updatePseudoTimestep_history_based(BoutReal previous_timestep, BoutReal previous_residual, BoutReal current_residual)#
BoutReal pid(BoutReal timestep, int nl_its, BoutReal max_dt)#

Updates the timestep.

void updateColoring()#

Updates the coloring using Jfd.

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.

BoutSnesEquationForm equation_form#

Form of the equation to solve.

std::string snes_type#
BoutReal atol#

Absolute tolerance.

BoutReal rtol#

Relative tolerance.

BoutReal stol#

Convergence tolerance.

int maxf#

Maximum number of function evaluations allowed in the solver (default: 10000)

int maxits#

Maximum nonlinear iterations.

int lower_its#
int upper_its#

Limits on iterations for timestep adjustment.

int max_snes_failures#

Maximum number of consecutive SNES failures before abort.

BoutReal timestep_factor_on_failure#
BoutReal timestep_factor_on_upper_its#
BoutReal timestep_factor_on_lower_its#
BoutPTCStrategy pseudo_strategy#

Strategy to use when setting timesteps.

BoutReal pseudo_alpha#

dt = alpha / residual

BoutReal pseudo_growth_factor#

Timestep increase 1.1 - 1.2.

BoutReal pseudo_reduction_factor#

Timestep decrease 0.5.

BoutReal pseudo_max_ratio#

Maximum timestep ratio between neighboring cells.

Vec dt_vec#

Each quantity can have its own timestep.

Field3D pseudo_residual#

Diagnostic output.

Field2D pseudo_residual_2d#
Field3D pseudo_timestep#

PID controller parameters.

bool pid_controller#

Use PID controller?

int target_its#

Use with caution! Not tested values.

Target number of nonlinear iterations for the PID controller.

BoutReal kP#

(0.6 - 0.8) Proportional parameter (main response to current step)

BoutReal kI#

(0.2 - 0.4) Integral parameter (smooths history of changes)

BoutReal kD#

(0.1 - 0.3) Derivative (dampens oscillation - optional)

bool pid_consider_failures#

Reduce timestep increases if recent solves have failed.

BoutReal recent_failure_rate#

Rolling average of recent failure rate.

const BoutReal inv_failure_window = 0.1#

1 / number of recent solves

int nl_its_prev#
int nl_its_prev2#
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.

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.

Vec output_x#

Solution to output. Used if interpolating.

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.

Mat Jfd#

Finite Difference Jacobian.

MatFDColoring fdcoloring = {nullptr}#

Matrix coloring context Jacobian evaluation

bool use_precon#

Use preconditioner.

std::string ksp_type#

Linear solver type.

bool kspsetinitialguessnonzero#

Set initial guess to non-zero.

int maxl#

Maximum linear iterations.

std::string pc_type#

Preconditioner type.

std::string pc_hypre_type#

Hypre preconditioner type.

std::string line_search_type#

Line search type.

bool matrix_free#

Use matrix free Jacobian.

bool matrix_free_operator#

Use matrix free Jacobian in the operator?

int lag_jacobian#

Re-use Jacobian.

bool use_coloring#

Use matrix coloring.

bool jacobian_recalculated#

Flag set when Jacobian is recalculated.

bool prune_jacobian#

Remove small elements in the Jacobian?

BoutReal prune_abstol#

Prune values with absolute values smaller than this.

BoutReal prune_fraction#

Prune if fraction of small elements is larger than this.

bool jacobian_pruned = {false}#

Has the Jacobian been pruned?

Mat Jfd_original#

Used to reset the Jacobian if over-pruned.

bool scale_rhs#

Scale time derivatives?

Vec rhs_scaling_factors#

Factors to multiply RHS function.

Vec jac_row_inv_norms#

1 / Norm of the rows of the Jacobian

bool scale_vars#

Scale individual variables?

Vec var_scaling_factors#

Factors to multiply variables when passing to user.

Vec scaled_x#

The values passed to the user RHS.

bool asinh_vars#

Evolve asinh(vars) to compress magnitudes while preserving signs.

const BoutReal asinh_scale = 1e-5#