10 #ifndef CT_RESIDJACEVAL_H 11 #define CT_RESIDJACEVAL_H 101 virtual int evalResidNJ(
const doublereal t,
const doublereal delta_t,
102 const doublereal*
const y,
103 const doublereal*
const ydot,
104 doublereal*
const resid,
107 const doublereal delta_x = 0.0);
109 virtual int eval(
const doublereal t,
const doublereal*
const y,
110 const doublereal*
const ydot,
111 doublereal*
const r);
113 virtual int getInitialConditions(
const doublereal t0, doublereal*
const y, doublereal*
const ydot);
125 virtual doublereal
filterNewStep(
const doublereal t,
const doublereal*
const ybase,
126 doublereal*
const step);
160 virtual int evalTimeTrackingEqns(
const doublereal t,
const doublereal delta_t,
const doublereal*
const y,
161 const doublereal*
const ydot);
177 const doublereal delta_t,
178 const doublereal*
const y,
179 const doublereal*
const ydot);
199 const doublereal*
const y,
200 const doublereal*
const ydot,
201 doublereal*
const delta_y,
202 const doublereal*
const solnWeights = 0);
214 virtual void calcSolnScales(
const doublereal t,
const doublereal*
const y,
215 const doublereal*
const y_old, doublereal*
const yScales);
230 virtual void user_out2(
const int ifunc,
const doublereal t,
231 const doublereal delta_t,
232 const doublereal*
const y,
233 const doublereal*
const ydot);
245 virtual void user_out(
const int ifunc,
const doublereal t,
247 const doublereal* ydot);
269 doublereal*
const rhs);
287 virtual int evalJacobian(
const doublereal t,
const doublereal delta_t, doublereal cj,
288 const doublereal*
const y,
const doublereal*
const ydot,
307 virtual int evalJacobianDP(
const doublereal t,
const doublereal delta_t, doublereal cj,
308 const doublereal*
const y,
309 const doublereal*
const ydot,
310 doublereal*
const* jacobianColPts,
311 doublereal*
const resid);
virtual int evalJacobianDP(const doublereal t, const doublereal delta_t, doublereal cj, const doublereal *const y, const doublereal *const ydot, doublereal *const *jacobianColPts, doublereal *const resid)
Calculate an analytical Jacobian and the residual at the current time and values. ...
virtual void user_out2(const int ifunc, const doublereal t, const doublereal delta_t, const doublereal *const y, const doublereal *const ydot)
This function may be used to create output at various points in the execution of an application...
Base residual calculation for the showSolution routine.
virtual int evalResidNJ(const doublereal t, const doublereal delta_t, const doublereal *const y, const doublereal *const ydot, doublereal *const resid, const ResidEval_Type_Enum evalType=Base_ResidEval, const int id_x=-1, const doublereal delta_x=0.0)
Evaluate the residual function.
virtual ResidJacEval * duplMyselfAsResidJacEval() const
Duplication routine for objects derived from residJacEval.
int neq_
Number of equations.
virtual int eval(const doublereal t, const doublereal *const y, const doublereal *const ydot, doublereal *const r)
Evaluate the residual function.
Base residual calculation for the Jacobian calculation.
virtual int getInitialConditions(const doublereal t0, doublereal *const y, doublereal *const ydot)
Fill in the initial conditions.
ResidJacEval(doublereal atol=1.0e-13)
Default constructor.
virtual void user_out(const int ifunc, const doublereal t, const doublereal *y, const doublereal *ydot)
This function may be used to create output at various points in the execution of an application...
ResidEval_Type_Enum
Differentiates the type of residual evaluations according to functionality.
Wrappers for the function evaluators for Nonlinear solvers and Time steppers.
virtual int nEquations() const
Return the number of equations in the equation system.
virtual int calcDeltaSolnVariables(const doublereal t, const doublereal *const y, const doublereal *const ydot, doublereal *const delta_y, const doublereal *const solnWeights=0)
Return a vector of delta y's for calculation of the numerical Jacobian.
virtual int matrixConditioning(doublereal *const matrix, const int nrows, doublereal *const rhs)
Multiply the matrix by another matrix that leads to better conditioning.
Virtual base class for DAE residual function evaluators.
virtual int evalJacobian(const doublereal t, const doublereal delta_t, doublereal cj, const doublereal *const y, const doublereal *const ydot, DenseMatrix &J, doublereal *const resid)
Calculate an analytical Jacobian and the residual at the current time and values. ...
doublereal m_atol
constant value of atol
virtual int evalTimeTrackingEqns(const doublereal t, const doublereal delta_t, const doublereal *const y, const doublereal *const ydot)
Evaluate the time tracking equations, if any.
virtual void calcSolnScales(const doublereal t, const doublereal *const y, const doublereal *const y_old, doublereal *const yScales)
Returns a vector of column scale factors that can be used to column scale Jacobians.
ResidJacEval & operator=(const ResidJacEval &right)
Headers for the DenseMatrix object, which deals with dense rectangular matrices and description of th...
Delta residual calculation for the Jacbobian calculation.
virtual doublereal filterSolnPrediction(const doublereal t, doublereal *const y)
Filter the solution predictions.
Namespace for the Cantera kernel.
virtual bool evalStoppingCritera(const doublereal t, const doublereal delta_t, const doublereal *const y, const doublereal *const ydot)
Evaluate any stopping criteria other than a final time limit.
void setAtol(doublereal atol)
Set a global value of the absolute tolerance.
virtual doublereal filterNewStep(const doublereal t, const doublereal *const ybase, doublereal *const step)
Filter the solution predictions.
Base residual calculation for the time-stepping function.
Base residual calculation containing any lagged components.
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...