Cantera  3.1.0a1
IdasIntegrator Class Reference

Wrapper for Sundials IDAS solver. More...

#include <IdasIntegrator.h>

Inheritance diagram for IdasIntegrator:
[legend]

Detailed Description

Wrapper for Sundials IDAS solver.

See also
FuncEval.h. Classes that use IdasIntegrator: FlowReactor

Definition at line 25 of file IdasIntegrator.h.

Public Member Functions

 IdasIntegrator ()
 Constructor. More...
 
void setTolerances (double reltol, size_t n, double *abstol) override
 Set error tolerances. More...
 
void setTolerances (double reltol, double abstol) override
 Set error tolerances. More...
 
void setSensitivityTolerances (double reltol, double abstol) override
 Set the sensitivity error tolerances. More...
 
void setLinearSolverType (const string &linearSolverType) override
 Set the linear solver type. More...
 
void initialize (double t0, FuncEval &func) override
 Initialize the integrator for a new problem. More...
 
void reinitialize (double t0, FuncEval &func) override
 
void integrate (double tout) override
 Integrate the system of equations. More...
 
double step (double tout) override
 Integrate the system of equations. More...
 
double & solution (size_t k) override
 The current value of the solution of equation k. More...
 
double * solution () override
 The current value of the solution of the system of equations. More...
 
int nEquations () const override
 The number of equations. More...
 
int maxOrder () const override
 
void setMaxOrder (int n) override
 Set the maximum integration order that will be used. More...
 
void setMaxStepSize (double hmax) override
 Set the maximum step size. More...
 
void setMaxSteps (int nmax) override
 Set the maximum number of time-steps the integrator can take before reaching the next output time. More...
 
int maxSteps () override
 Returns the maximum number of time-steps the integrator can take before reaching the next output time. More...
 
void setMaxErrTestFails (int n) override
 Set the maximum permissible number of error test failures. More...
 
AnyMap solverStats () const override
 Get solver stats from integrator. More...
 
int nSensParams () override
 
double sensitivity (size_t k, size_t p) override
 
string getErrorInfo (int N)
 Returns a string listing the weighted error estimates associated with each solution component. More...
 
int maxNonlinIterations () const override
 
void setMaxNonlinIterations (int n) override
 
int maxNonlinConvFailures () const override
 
void setMaxNonlinConvFailures (int n) override
 
bool algebraicInErrorTest () const override
 
void includeAlgebraicInErrorTest (bool yesno) override
 
void setMethod (MethodType t) override
 Set the solution method. More...
 
- Public Member Functions inherited from Integrator
 Integrator ()
 Default Constructor. More...
 
virtual ~Integrator ()
 Destructor. More...
 
virtual void setPreconditioner (shared_ptr< PreconditionerBase > preconditioner)
 Set preconditioner used by the linear solver. More...
 
virtual void preconditionerSolve (size_t stateSize, double *rhs, double *output)
 Solve a linear system Ax=b where A is the preconditioner. More...
 
virtual PreconditionerSide preconditionerSide ()
 Return the side of the system on which the preconditioner is applied. More...
 
virtual shared_ptr< PreconditionerBasepreconditioner ()
 Return preconditioner reference to object. More...
 
virtual string linearSolverType () const
 Return the integrator problem type. More...
 
virtual double * derivative (double tout, int n)
 n-th derivative of the output function at time tout. More...
 
virtual int lastOrder () const
 Order used during the last solution step. More...
 
virtual int nEvals () const
 The number of function evaluations. More...
 
virtual void setMinStepSize (double hmin)
 Set the minimum step size. More...
 
virtual void setBandwidth (int N_Upper, int N_Lower)
 

Public Attributes

string m_error_message
 Error message information provide by IDAS. More...
 

Protected Member Functions

void applyOptions ()
 Applies user-specified options to the underlying IDAS solver. More...
 

Private Member Functions

void sensInit (double t0, FuncEval &func)
 
void checkError (long flag, const string &ctMethod, const string &idaMethod) const
 Check whether an IDAS method indicated an error. More...
 

Private Attributes

size_t m_neq
 Number of equations/variables in the system. More...
 
void * m_ida_mem = nullptr
 Pointer to the IDA memory for the problem. More...
 
void * m_linsol = nullptr
 Sundials linear solver object. More...
 
void * m_linsol_matrix = nullptr
 matrix used by Sundials More...
 
SundialsContext m_sundials_ctx
 SUNContext object for Sundials>=6.0. More...
 
FuncEvalm_func = nullptr
 Object implementing the DAE residual function \( f(t, y, \dot{y}) = 0 \). More...
 
double m_t0 = 0.0
 The start time for the integrator. More...
 
double m_time
 The current integrator time, corresponding to m_y. More...
 
double m_tInteg
 The latest time reached by the integrator. May be greater than m_time. More...
 
N_Vector m_y = nullptr
 The current system state. More...
 
N_Vector m_ydot = nullptr
 The time derivatives of the system state. More...
 
N_Vector m_abstol = nullptr
 Absolute tolerances for each state variable. More...
 
string m_type = "DENSE"
 The linear solver type. More...
 
int m_itol
 Flag indicating whether scalar (IDA_SS) or vector (IDA_SV) absolute tolerances are being used. More...
 
int m_maxord = 0
 Maximum order allowed for the BDF method. More...
 
double m_reltol = 1.0e-9
 Relative tolerance for all state variables. More...
 
double m_abstols = 1.0e-15
 Scalar absolute tolerance. More...
 
double m_reltolsens
 Scalar relative tolerance for sensitivities. More...
 
double m_abstolsens
 Scalar absolute tolerance for sensitivities. More...
 
size_t m_nabs = 0
 ! Number of variables for which absolute tolerances were provided More...
 
double m_hmax = 0.0
 Maximum time step size. Zero means infinity. More...
 
int m_maxsteps = 20000
 Maximum number of internal steps taken in a call to integrate() More...
 
int m_maxErrTestFails = -1
 Maximum number of error test failures in attempting one step. More...
 
size_t m_np
 Number of sensitivity parameters. More...
 
N_Vector * m_yS = nullptr
 Sensitivities of y, size m_np by m_neq. More...
 
N_Vector * m_ySdot = nullptr
 Sensitivities of ydot, size m_np by m_neq. More...
 
N_Vector m_constraints = nullptr
 
bool m_sens_ok
 Indicates whether the sensitivities stored in m_yS and m_ySdot have been updated for the current integrator time. More...
 
int m_maxNonlinIters = 0
 Maximum number of nonlinear solver iterations at one solution. More...
 
int m_maxNonlinConvFails = -1
 Maximum number of nonlinear convergence failures. More...
 
bool m_setSuppressAlg = false
 If true, the algebraic variables don't contribute to error tolerances. More...
 
double m_init_step = 1e-14
 Initial IDA step size. More...
 

Additional Inherited Members

- Protected Attributes inherited from Integrator
shared_ptr< PreconditionerBasem_preconditioner
 Pointer to preconditioner object used in integration which is set by setPreconditioner and initialized inside of ReactorNet::initialize() More...
 
PreconditionerSide m_prec_side = PreconditionerSide::NO_PRECONDITION
 Type of preconditioning used in applyOptions. More...
 

Constructor & Destructor Documentation

◆ IdasIntegrator()

Constructor.

Default settings: dense Jacobian, no user-supplied Jacobian function, Newton iteration.

Definition at line 59 of file IdasIntegrator.cpp.

Member Function Documentation

◆ setTolerances() [1/2]

void setTolerances ( double  reltol,
size_t  n,
double *  abstol 
)
overridevirtual

Set error tolerances.

Parameters
reltolscalar relative tolerance
nNumber of equations
abstolarray of N absolute tolerance values

Reimplemented from Integrator.

Definition at line 93 of file IdasIntegrator.cpp.

◆ setTolerances() [2/2]

void setTolerances ( double  reltol,
double  abstol 
)
overridevirtual

Set error tolerances.

Parameters
reltolscalar relative tolerance
abstolscalar absolute tolerance

Reimplemented from Integrator.

Definition at line 109 of file IdasIntegrator.cpp.

◆ setSensitivityTolerances()

void setSensitivityTolerances ( double  reltol,
double  abstol 
)
overridevirtual

Set the sensitivity error tolerances.

Parameters
reltolscalar relative tolerance
abstolscalar absolute tolerance

Reimplemented from Integrator.

Definition at line 116 of file IdasIntegrator.cpp.

◆ setLinearSolverType()

void setLinearSolverType ( const string &  linSolverType)
overridevirtual

Set the linear solver type.

Parameters
linSolverTypeType of the linear solver

Reimplemented from Integrator.

Definition at line 123 of file IdasIntegrator.cpp.

◆ initialize()

void initialize ( double  t0,
FuncEval func 
)
overridevirtual

Initialize the integrator for a new problem.

Call after all options have been set.

Parameters
t0initial time
funcRHS evaluator object for system of equations.

Reimplemented from Integrator.

Definition at line 218 of file IdasIntegrator.cpp.

◆ integrate()

void integrate ( double  tout)
overridevirtual

Integrate the system of equations.

Parameters
toutIntegrate to this time. Note that this is the absolute time value, not a time interval.

Reimplemented from Integrator.

Definition at line 445 of file IdasIntegrator.cpp.

◆ step()

double step ( double  tout)
overridevirtual

Integrate the system of equations.

Parameters
toutintegrate to this time. Note that this is the absolute time value, not a time interval.

Reimplemented from Integrator.

Definition at line 482 of file IdasIntegrator.cpp.

◆ solution() [1/2]

double & solution ( size_t  k)
overridevirtual

The current value of the solution of equation k.

Reimplemented from Integrator.

Definition at line 83 of file IdasIntegrator.cpp.

◆ solution() [2/2]

double * solution ( )
overridevirtual

The current value of the solution of the system of equations.

Reimplemented from Integrator.

Definition at line 88 of file IdasIntegrator.cpp.

◆ nEquations()

int nEquations ( ) const
inlineoverridevirtual

The number of equations.

Reimplemented from Integrator.

Definition at line 44 of file IdasIntegrator.h.

◆ setMaxOrder()

void setMaxOrder ( int  n)
overridevirtual

Set the maximum integration order that will be used.

Reimplemented from Integrator.

Definition at line 128 of file IdasIntegrator.cpp.

◆ setMaxStepSize()

void setMaxStepSize ( double  hmax)
overridevirtual

Set the maximum step size.

Reimplemented from Integrator.

Definition at line 137 of file IdasIntegrator.cpp.

◆ setMaxSteps()

void setMaxSteps ( int  nmax)
overridevirtual

Set the maximum number of time-steps the integrator can take before reaching the next output time.

Parameters
nmaxThe maximum number of steps, setting this value to zero disables this option.

Reimplemented from Integrator.

Definition at line 146 of file IdasIntegrator.cpp.

◆ maxSteps()

int maxSteps ( )
overridevirtual

Returns the maximum number of time-steps the integrator can take before reaching the next output time.

Reimplemented from Integrator.

Definition at line 154 of file IdasIntegrator.cpp.

◆ setMaxErrTestFails()

void setMaxErrTestFails ( int  n)
overridevirtual

Set the maximum permissible number of error test failures.

Reimplemented from Integrator.

Definition at line 159 of file IdasIntegrator.cpp.

◆ solverStats()

AnyMap solverStats ( ) const
overridevirtual

Get solver stats from integrator.

Reimplemented from Integrator.

Definition at line 167 of file IdasIntegrator.cpp.

◆ getErrorInfo()

string getErrorInfo ( int  N)

Returns a string listing the weighted error estimates associated with each solution component.

This information can be used to identify which variables are responsible for integrator failures or unexpected small timesteps.

Definition at line 524 of file IdasIntegrator.cpp.

◆ setMethod()

void setMethod ( MethodType  t)
overridevirtual

Set the solution method.

Reimplemented from Integrator.

Definition at line 564 of file IdasIntegrator.cpp.

◆ applyOptions()

void applyOptions ( )
protected

Applies user-specified options to the underlying IDAS solver.

Called during integrator initialization or reinitialization.

Definition at line 324 of file IdasIntegrator.cpp.

◆ checkError()

void checkError ( long  flag,
const string &  ctMethod,
const string &  idaMethod 
) const
private

Check whether an IDAS method indicated an error.

If so, throw an exception containing the method name and the error code stashed by the ida_err() function.

Definition at line 548 of file IdasIntegrator.cpp.

Member Data Documentation

◆ m_error_message

string m_error_message

Error message information provide by IDAS.

Definition at line 68 of file IdasIntegrator.h.

◆ m_neq

size_t m_neq
private

Number of equations/variables in the system.

Definition at line 100 of file IdasIntegrator.h.

◆ m_ida_mem

void* m_ida_mem = nullptr
private

Pointer to the IDA memory for the problem.

Definition at line 101 of file IdasIntegrator.h.

◆ m_linsol

void* m_linsol = nullptr
private

Sundials linear solver object.

Definition at line 102 of file IdasIntegrator.h.

◆ m_linsol_matrix

void* m_linsol_matrix = nullptr
private

matrix used by Sundials

Definition at line 103 of file IdasIntegrator.h.

◆ m_sundials_ctx

SundialsContext m_sundials_ctx
private

SUNContext object for Sundials>=6.0.

Definition at line 104 of file IdasIntegrator.h.

◆ m_func

FuncEval* m_func = nullptr
private

Object implementing the DAE residual function \( f(t, y, \dot{y}) = 0 \).

Definition at line 107 of file IdasIntegrator.h.

◆ m_t0

double m_t0 = 0.0
private

The start time for the integrator.

Definition at line 109 of file IdasIntegrator.h.

◆ m_time

double m_time
private

The current integrator time, corresponding to m_y.

Definition at line 110 of file IdasIntegrator.h.

◆ m_tInteg

double m_tInteg
private

The latest time reached by the integrator. May be greater than m_time.

Definition at line 113 of file IdasIntegrator.h.

◆ m_y

N_Vector m_y = nullptr
private

The current system state.

Definition at line 115 of file IdasIntegrator.h.

◆ m_ydot

N_Vector m_ydot = nullptr
private

The time derivatives of the system state.

Definition at line 116 of file IdasIntegrator.h.

◆ m_abstol

N_Vector m_abstol = nullptr
private

Absolute tolerances for each state variable.

Definition at line 117 of file IdasIntegrator.h.

◆ m_type

string m_type = "DENSE"
private

The linear solver type.

See also
setLinearSolverType()

Definition at line 118 of file IdasIntegrator.h.

◆ m_itol

int m_itol
private

Flag indicating whether scalar (IDA_SS) or vector (IDA_SV) absolute tolerances are being used.

Definition at line 122 of file IdasIntegrator.h.

◆ m_maxord

int m_maxord = 0
private

Maximum order allowed for the BDF method.

Definition at line 124 of file IdasIntegrator.h.

◆ m_reltol

double m_reltol = 1.0e-9
private

Relative tolerance for all state variables.

Definition at line 125 of file IdasIntegrator.h.

◆ m_abstols

double m_abstols = 1.0e-15
private

Scalar absolute tolerance.

Definition at line 126 of file IdasIntegrator.h.

◆ m_reltolsens

double m_reltolsens
private

Scalar relative tolerance for sensitivities.

Definition at line 127 of file IdasIntegrator.h.

◆ m_abstolsens

double m_abstolsens
private

Scalar absolute tolerance for sensitivities.

Definition at line 128 of file IdasIntegrator.h.

◆ m_nabs

size_t m_nabs = 0
private

! Number of variables for which absolute tolerances were provided

Definition at line 131 of file IdasIntegrator.h.

◆ m_hmax

double m_hmax = 0.0
private

Maximum time step size. Zero means infinity.

Definition at line 133 of file IdasIntegrator.h.

◆ m_maxsteps

int m_maxsteps = 20000
private

Maximum number of internal steps taken in a call to integrate()

Definition at line 136 of file IdasIntegrator.h.

◆ m_maxErrTestFails

int m_maxErrTestFails = -1
private

Maximum number of error test failures in attempting one step.

Definition at line 139 of file IdasIntegrator.h.

◆ m_np

size_t m_np
private

Number of sensitivity parameters.

Definition at line 141 of file IdasIntegrator.h.

◆ m_yS

N_Vector* m_yS = nullptr
private

Sensitivities of y, size m_np by m_neq.

Definition at line 142 of file IdasIntegrator.h.

◆ m_ySdot

N_Vector* m_ySdot = nullptr
private

Sensitivities of ydot, size m_np by m_neq.

Definition at line 143 of file IdasIntegrator.h.

◆ m_sens_ok

bool m_sens_ok
private

Indicates whether the sensitivities stored in m_yS and m_ySdot have been updated for the current integrator time.

Definition at line 148 of file IdasIntegrator.h.

◆ m_maxNonlinIters

int m_maxNonlinIters = 0
private

Maximum number of nonlinear solver iterations at one solution.

If zero, this is the default of 4.

Definition at line 154 of file IdasIntegrator.h.

◆ m_maxNonlinConvFails

int m_maxNonlinConvFails = -1
private

Maximum number of nonlinear convergence failures.

Definition at line 157 of file IdasIntegrator.h.

◆ m_setSuppressAlg

bool m_setSuppressAlg = false
private

If true, the algebraic variables don't contribute to error tolerances.

Definition at line 160 of file IdasIntegrator.h.

◆ m_init_step

double m_init_step = 1e-14
private

Initial IDA step size.

Definition at line 163 of file IdasIntegrator.h.


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