File cvode.hxx#

class CvodeSolver : public Solver#

Public Functions

explicit CvodeSolver(Options *opts = nullptr)#
~CvodeSolver()#
inline virtual BoutReal getCurrentTimestep() override#

Return the current internal timestep.

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.

BoutReal run(BoutReal tout)#
virtual void resetInternalFields() override#

Should wipe out internal field vector and reset from current field object data.

void rhs(BoutReal t, BoutReal *udata, BoutReal *dudata)#
void pre(BoutReal t, BoutReal gamma, BoutReal delta, BoutReal *udata, BoutReal *rvec, BoutReal *zvec)#
void jac(BoutReal t, BoutReal *ydata, BoutReal *vdata, BoutReal *Jvdata)#

Private Functions

void set_vector_option_values(BoutReal *option_data, std::vector<BoutReal> &f2dtols, std::vector<BoutReal> &f3dtols)#
void loop_vector_option_values_op(Ind2D i2d, BoutReal *option_data, int &p, std::vector<BoutReal> &f2dtols, std::vector<BoutReal> &f3dtols, bool bndry)#
template<class FieldType>
std::vector<BoutReal> create_constraints(const std::vector<VarStr<FieldType>> &fields)#

Private Members

BoutReal hcur#
bool diagnose = {false}#
N_Vector uvec = {nullptr}#
void *cvode_mem = {nullptr}#
BoutReal pre_Wtime = {0.0}#
int pre_ncalls = {0}#
bool adams_moulton#

Use Adams Moulton implicit multistep. Otherwise BDF method.

bool func_iter#

Use functional iteration instead of Newton.

int max_order#

Maximum order of method to use. < 0 means no limit.

bool stablimdet#
BoutReal abstol#

Absolute tolerance.

BoutReal reltol#

Relative tolerance.

bool use_vector_abstol#

Use separate absolute tolerance for each field.

int mxsteps#

Maximum number of internal steps between outputs.

BoutReal max_timestep#

Maximum time step size.

BoutReal min_timestep#

Minimum time step size.

BoutReal start_timestep#

Starting time step. < 0 then chosen by CVODE.

int mxorder#

Maximum order.

int max_nonlinear_iterations#

Maximum number of nonlinear iterations allowed by CVODE before reducing timestep. CVODE default (used if this option is negative) is 3

bool apply_positivity_constraints#

Use CVODE function CVodeSetConstraints to constrain variables - the constraint to be applied is set by the positivity_constraint option in the subsection for each variable

int maxl#

Maximum number of linear iterations.

bool use_precon#

Use preconditioner?

bool rightprec#

Use right preconditioner? Otherwise use left.

bool use_jacobian#
BoutReal cvode_nonlinear_convergence_coef#
BoutReal cvode_linear_convergence_coef#
int nsteps = {0}#
int nfevals = {0}#
int nniters = {0}#
int npevals = {0}#
int nliters = {0}#
BoutReal last_step = {0.0}#
int last_order = {0}#
int num_fails = {0}#
int nonlin_fails = {0}#
int stab_lims = {0}#
bool cvode_initialised = false#
SUNLinearSolver sun_solver = {nullptr}#

SPGMR solver structure.

SUNNonlinearSolver nonlinear_solver = {nullptr}#

Solver for functional iterations for Adams-Moulton.

sundials::Context suncontext#

Context for SUNDIALS memory allocations.