Cantera
2.2.1
|
#include <BEulerInt.h>
Public Member Functions | |
BEulerInt () | |
virtual void | setTolerances (double reltol, size_t n, double *abstol) |
Set or reset the number of equations. More... | |
virtual void | setTolerances (double reltol, double abstol) |
Set error tolerances. More... | |
virtual void | setProblemType (int probtype) |
Set the problem type. More... | |
virtual void | initializeRJE (double t0, ResidJacEval &func) |
Find the initial conditions for y and ydot. More... | |
virtual void | reinitializeRJE (double t0, ResidJacEval &func) |
virtual double | integrateRJE (double tout, double tinit=0.0) |
virtual doublereal | step (double tout) |
Integrate the system of equations. More... | |
virtual void | setSolnWeights () |
Set the solution weights. More... | |
virtual double & | solution (size_t k) |
The current value of the solution of equation k. More... | |
double * | solution () |
The current value of the solution of the system of equations. More... | |
int | nEquations () const |
The number of equations. More... | |
virtual int | nEvals () const |
Return the total number of function evaluations. More... | |
virtual void | setMethodBEMT (BEulerMethodType t) |
virtual void | setIterator (IterType t) |
Set the linear iterator. More... | |
virtual void | setMaxStep (double hmax) |
virtual void | setMaxNumTimeSteps (int) |
virtual void | setNumInitialConstantDeltaTSteps (int) |
void | print_solnDelta_norm_contrib (const double *const soln0, const char *const s0, const double *const soln1, const char *const s1, const char *const title, const double *const y0, const double *const y1, double damp, int num_entries) |
virtual void | setPrintSolnOptions (int printSolnStepInterval, int printSolnNumberToTout, int printSolnFirstSteps=0, bool dumpJacobians=false) |
This routine controls when the solution is printed. More... | |
void | setNonLinOptions (int min_newt_its=0, bool matrixConditioning=false, bool colScaling=false, bool rowScaling=true) |
Set the options for the nonlinear method. More... | |
virtual void | setPrintFlag (int print_flag) |
virtual void | setColumnScales () |
Set the column scaling vector at the current time. More... | |
virtual double | soln_error_norm (const double *const, bool printLargest=false) |
Calculate the solution error norm. More... | |
virtual void | setInitialTimeStep (double delta_t) |
void | beuler_jac (GeneralMatrix &J, double *const f, double, double, double *const, double *const, int) |
Public Member Functions inherited from Integrator | |
Integrator () | |
Default Constructor. More... | |
virtual | ~Integrator () |
Destructor. More... | |
virtual void | setSensitivityTolerances (doublereal reltol, doublereal abstol) |
Set the sensitivity error tolerances. More... | |
virtual void | initialize (doublereal t0, FuncEval &func) |
Initialize the integrator for a new problem. More... | |
virtual void | reinitialize (doublereal t0, FuncEval &func) |
virtual void | integrate (doublereal tout) |
Integrate the system of equations. More... | |
virtual void | setMaxOrder (int n) |
Set the maximum integration order that will be used. More... | |
virtual void | setMethod (MethodType t) |
Set the solution method. More... | |
virtual void | setMaxStepSize (double hmax) |
Set the maximum step size. More... | |
virtual void | setMinStepSize (double hmin) |
Set the minimum step size. More... | |
virtual void | setMaxErrTestFails (int n) |
Set the maximum permissible number of error test failures. More... | |
virtual void | setMaxSteps (int nmax) |
virtual void | setBandwidth (int N_Upper, int N_Lower) |
virtual int | nSensParams () |
virtual double | sensitivity (size_t k, size_t p) |
Protected Member Functions | |
void | internalMalloc () |
Internal routine that sets up the fixed length storage based on the size of the problem to solve. More... | |
void | calc_y_pred (int) |
void | calc_ydot (int, double *, double *) |
double | time_error_norm () |
double | time_step_control (int m_order, double time_error_factor) |
int | solve_nonlinear_problem (double *const y_comm, double *const ydot_comm, double CJ, double time_curr, GeneralMatrix &jac, int &num_newt_its, int &num_linear_solves, int &num_backtracks, int loglevel) |
Solve a nonlinear system. More... | |
void | doNewtonSolve (double, double *, double *, double *, GeneralMatrix &, int) |
Compute the undamped Newton step. More... | |
double | boundStep (const double *const y, const double *const step0, int loglevel) |
Bound the Newton step while relaxing the solution. More... | |
int | dampStep (double, const double *, const double *, const double *, double *, double *, double *, double &, GeneralMatrix &, int &, bool, int &) |
void | computeResidWts (GeneralMatrix &jac) |
Compute Residual Weights. More... | |
double | filterNewStep (double, double *, double *) |
Filter a new step. More... | |
double | getPrintTime (double time_current) |
Get the next time to print out. More... | |
Protected Attributes | |
IterType | m_iter |
IterType is used to specify how the nonlinear equations are to be relaxed at each time step. More... | |
BEulerMethodType | m_method |
MethodType is used to specify how the time step is to be chosen. More... | |
int | m_jacFormMethod |
m_jacFormMethod determines how a matrix is formed. More... | |
bool | m_rowScaling |
m_rowScaling is a boolean. More... | |
bool | m_colScaling |
m_colScaling is a boolean. More... | |
bool | m_matrixConditioning |
m_matrixConditioning is a boolean. More... | |
int | m_itol |
If m_itol =1 then each component has an individual value of atol. More... | |
double | m_reltol |
Relative time truncation error tolerances. More... | |
double | m_abstols |
Absolute time truncation error tolerances, when uniform for all variables. More... | |
vector_fp | m_abstol |
Vector of absolute time truncation error tolerance when not uniform for all variables. More... | |
vector_fp | m_ewt |
Error Weights. This is a surprisingly important quantity. More... | |
double | m_hmax |
Maximum step size. More... | |
int | m_maxord |
Maximum integration order. More... | |
int | m_order |
Current integration order. More... | |
int | m_time_step_num |
Time step number. More... | |
int | m_time_step_attempts |
int | m_max_time_step_attempts |
Max time steps allowed. More... | |
int | m_numInitialConstantDeltaTSteps |
Number of initial time steps to take where the time truncation error tolerances are not checked. More... | |
int | m_failure_counter |
Failure Counter -> keeps track of the number of consecutive failures. More... | |
int | m_min_newt_its |
Minimum Number of Newton Iterations per nonlinear step. default = 0. More... | |
int | m_printSolnStepInterval |
Step Interval at which to print out the solution default = 1; If set to zero, there is no printout. More... | |
int | m_printSolnNumberToTout |
Number of evenly spaced printouts of the solution If zero, there is no printout from this option default 1 If set to zero there is no printout. More... | |
int | m_printSolnFirstSteps |
Number of initial steps that the solution is printed out. default = 0. More... | |
bool | m_dumpJacobians |
Dump Jacobians to disk. default false. More... | |
int | m_neq |
Number of equations in the ode integrator. More... | |
vector_fp | m_y_n |
vector_fp | m_y_nm1 |
vector_fp | m_y_pred_n |
vector_fp | m_ydot_n |
vector_fp | m_ydot_nm1 |
double | m_t0 |
Initial time at the start of the integration. More... | |
double | m_time_final |
Final time. More... | |
double | time_n |
double | time_nm1 |
double | time_nm2 |
double | delta_t_n |
double | delta_t_nm1 |
double | delta_t_nm2 |
double | delta_t_np1 |
double | delta_t_max |
Maximum permissible time step. More... | |
vector_fp | m_resid |
vector_fp | m_residWts |
vector_fp | m_wksp |
ResidJacEval * | m_func |
vector_fp | m_rowScales |
vector_fp | m_colScales |
GeneralMatrix * | tdjac_ptr |
Pointer to the Jacobian representing the time dependent problem. More... | |
int | m_print_flag |
Determines the level of printing for each time step. More... | |
int | m_nfe |
Number of function evaluations. More... | |
int | m_nJacEval |
Number of Jacobian Evaluations and factorization steps (they are the same) More... | |
int | m_numTotalNewtIts |
Number of total Newton iterations. More... | |
int | m_numTotalLinearSolves |
Total number of linear iterations. More... | |
int | m_numTotalConvFails |
Total number of convergence failures. More... | |
int | m_numTotalTruncFails |
Total Number of time truncation error failures. More... | |
int | num_failures |
Wrapper class for 'beuler' integrator We derive the class from the class Integrator
Definition at line 56 of file BEulerInt.h.
BEulerInt | ( | ) |
Constructor. Default settings: dense Jacobian, no user-supplied Jacobian function, Newton iteration.
Definition at line 26 of file BEulerInt.cpp.
References Cantera::warn_deprecated().
|
virtual |
Set or reset the number of equations.
Set error tolerances.
reltol | scalar relative tolerance |
n | Number of equations |
abstol | array of N absolute tolerance values |
Reimplemented from Integrator.
Definition at line 78 of file BEulerInt.cpp.
References BEulerInt::m_abstol, BEulerInt::m_itol, BEulerInt::m_neq, and BEulerInt::m_reltol.
|
virtual |
Set error tolerances.
reltol | scalar relative tolerance |
abstol | scalar absolute tolerance |
Reimplemented from Integrator.
Definition at line 92 of file BEulerInt.cpp.
References BEulerInt::m_abstols, BEulerInt::m_itol, and BEulerInt::m_reltol.
|
virtual |
Set the problem type.
probtype | Type of the problem |
Reimplemented from Integrator.
Definition at line 99 of file BEulerInt.cpp.
References BEulerInt::m_jacFormMethod.
|
virtual |
Find the initial conditions for y and ydot.
Definition at line 165 of file BEulerInt.cpp.
References ResidJacEval::getInitialConditions(), BEulerInt::internalMalloc(), BEulerInt::m_neq, BEulerInt::m_t0, and ResidJacEval::nEquations().
|
virtual |
Integrate the system of equations.
tout | integrate to this time. Note that this is the absolute time value, not a time interval. |
Reimplemented from Integrator.
Definition at line 888 of file BEulerInt.cpp.
References BEulerInt::calc_y_pred(), BEulerInt::calc_ydot(), BEulerInt::delta_t_max, BEulerInt::filterNewStep(), BEulerInt::m_failure_counter, BEulerInt::m_method, BEulerInt::m_numInitialConstantDeltaTSteps, BEulerInt::m_numTotalConvFails, BEulerInt::m_numTotalTruncFails, BEulerInt::m_order, BEulerInt::m_print_flag, BEulerInt::m_time_step_num, BEulerInt::setSolnWeights(), BEulerInt::solve_nonlinear_problem(), BEulerInt::tdjac_ptr, BEulerInt::time_error_norm(), and BEulerInt::time_step_control().
|
virtual |
Set the solution weights.
This is a very important routine as it affects quite a few operations involving convergence.
Definition at line 250 of file BEulerInt.cpp.
References BEulerInt::m_abstol, BEulerInt::m_abstols, BEulerInt::m_ewt, BEulerInt::m_itol, BEulerInt::m_neq, and BEulerInt::m_reltol.
Referenced by BEulerInt::step().
|
inlinevirtual |
The current value of the solution of equation k.
Reimplemented from Integrator.
Definition at line 84 of file BEulerInt.h.
|
inlinevirtual |
The current value of the solution of the system of equations.
Reimplemented from Integrator.
Definition at line 87 of file BEulerInt.h.
|
inlinevirtual |
The number of equations.
Reimplemented from Integrator.
Definition at line 90 of file BEulerInt.h.
References BEulerInt::m_neq.
|
virtual |
Return the total number of function evaluations.
Reimplemented from Integrator.
Definition at line 225 of file BEulerInt.cpp.
References BEulerInt::m_nfe.
|
virtual |
Set the linear iterator.
Reimplemented from Integrator.
Definition at line 135 of file BEulerInt.cpp.
References BEulerInt::m_iter.
|
virtual |
This routine controls when the solution is printed.
printSolnStepInterval | If greater than 0, then the soln is printed every printSolnStepInterval steps. |
printSolnNumberToTout | The solution is printed at regular invervals a total of "printSolnNumberToTout" times. |
printSolnFirstSteps | The solution is printed out the first "printSolnFirstSteps" steps. After these steps the other parameters determine the printing. default = 0 |
dumpJacobians | Dump Jacobians to disk. |
Definition at line 124 of file BEulerInt.cpp.
References BEulerInt::m_dumpJacobians, BEulerInt::m_printSolnFirstSteps, BEulerInt::m_printSolnNumberToTout, and BEulerInt::m_printSolnStepInterval.
void setNonLinOptions | ( | int | min_newt_its = 0 , |
bool | matrixConditioning = false , |
||
bool | colScaling = false , |
||
bool | rowScaling = true |
||
) |
Set the options for the nonlinear method.
Defaults are set in the .h file. These are the defaults: min_newt_its = 0 matrixConditioning = false colScaling = false rowScaling = true
Definition at line 140 of file BEulerInt.cpp.
References BEulerInt::m_colScaling, BEulerInt::m_matrixConditioning, BEulerInt::m_min_newt_its, BEulerInt::m_neq, and BEulerInt::m_rowScaling.
|
virtual |
Set the column scaling vector at the current time.
Definition at line 266 of file BEulerInt.cpp.
References ResidJacEval::calcSolnScales().
Referenced by BEulerInt::doNewtonSolve().
|
virtual |
Calculate the solution error norm.
if printLargest is true, then a table of the largest values is printed to standard output.
Definition at line 1295 of file BEulerInt.cpp.
References Cantera::error(), BEulerInt::m_ewt, and BEulerInt::m_neq.
Referenced by BEulerInt::dampStep().
void beuler_jac | ( | GeneralMatrix & | J, |
double *const | f, | ||
double | time_curr, | ||
double | CJ, | ||
double * const | y, | ||
double * const | ydot, | ||
int | num_newt_its | ||
) |
Function called by to evaluate the Jacobian matrix and the current residual at the current time step.
J | = Jacobian matrix to be filled in |
f | = Right hand side. This routine returns the current value of the RHS (output), so that it does not have to be computed again. |
Definition at line 483 of file BEulerInt.cpp.
References ResidJacEval::calcDeltaSolnVariables(), GeneralMatrix::clearFactorFlag(), ResidJacEval::evalJacobian(), ResidJacEval::evalResidNJ(), Cantera::JacBase_ResidEval, Cantera::JacDelta_ResidEval, BEulerInt::m_ewt, BEulerInt::m_jacFormMethod, BEulerInt::m_neq, BEulerInt::m_nfe, BEulerInt::m_nJacEval, GeneralMatrix::ptrColumn(), and Cantera::subtractRD().
Referenced by BEulerInt::solve_nonlinear_problem().
|
protected |
Internal routine that sets up the fixed length storage based on the size of the problem to solve.
Definition at line 230 of file BEulerInt.cpp.
References BEulerInt::m_colScaling, BEulerInt::m_ewt, BEulerInt::m_neq, BEulerInt::m_rowScaling, and BEulerInt::tdjac_ptr.
Referenced by BEulerInt::initializeRJE().
|
protected |
Function to calculate the predicted solution vector, m_y_pred_n for the (n+1)th time step. This routine can be used by a first order - forward Euler / backward Euler predictor / corrector method or for a second order Adams-Bashforth / Trapezoidal Rule predictor / corrector method. See Nachos documentation Sand86-1816 and Gresho, Lee, Sani LLNL report UCRL - 83282 for more information.
on input:
N - number of unknowns order - indicates order of method = 1 -> first order forward Euler/backward Euler predictor/corrector = 2 -> second order Adams-Bashforth/Trapezoidal Rule predictor/corrector
delta_t_n - magnitude of time step at time n (i.e., = t_n+1 - t_n) delta_t_nm1 - magnitude of time step at time n - 1 (i.e., = t_n - t_n-1) y_n[] - solution vector at time n y_dot_n[] - acceleration vector from the predictor at time n y_dot_nm1[] - acceleration vector from the predictor at time n - 1
on output:
m_y_pred_n[] - predicted solution vector at time n + 1
Definition at line 590 of file BEulerInt.cpp.
References ResidJacEval::filterSolnPrediction(), and BEulerInt::m_neq.
Referenced by BEulerInt::step().
|
protected |
Function to calculate the acceleration vector ydot for the first or second order predictor/corrector time integrator. This routine can be called by a first order - forward Euler / backward Euler predictor / corrector or for a second order Adams - Bashforth / Trapezoidal Rule predictor / corrector. See Nachos documentation Sand86-1816 and Gresho, Lee, Sani LLNL report UCRL - 83282 for more information.
on input:
N - number of local unknowns on the processor This is equal to internal plus border unknowns. order - indicates order of method = 1 -> first order forward Euler/backward Euler predictor/corrector = 2 -> second order Adams-Bashforth/Trapezoidal Rule predictor/corrector
delta_t_n - Magnitude of the current time step at time n (i.e., = t_n - t_n-1) y_curr[] - Current Solution vector at time n y_nm1[] - Solution vector at time n-1 ydot_nm1[] - Acceleration vector at time n-1
on output:
ydot_curr[] - Current acceleration vector at time n
Note we use the current attribute to denote the possibility that y_curr[] may not be equal to m_y_n[] during the nonlinear solve because we may be using a look-ahead scheme.
Definition at line 618 of file BEulerInt.cpp.
References BEulerInt::m_neq.
Referenced by BEulerInt::dampStep(), BEulerInt::solve_nonlinear_problem(), and BEulerInt::step().
|
protected |
Calculates the time step truncation error estimate from a very simple formula based on Gresho et al. This routine can be called for a first order - forward Euler/backward Euler predictor/ corrector and for a second order Adams- Bashforth/Trapezoidal Rule predictor/corrector. See Nachos documentation Sand86-1816 and Gresho, Lee, LLNL report UCRL - 83282 for more information.
on input:
abs_error - Generic absolute error tolerance rel_error - Generic relative error tolerance x_coor[] - Solution vector from the implicit corrector x_pred_n[] - Solution vector from the explicit predictor
on output:
delta_t_n - Magnitude of next time step at time t_n+1 delta_t_nm1 - Magnitude of previous time step at time t_n
Definition at line 639 of file BEulerInt.cpp.
References Cantera::error(), BEulerInt::m_ewt, BEulerInt::m_neq, and BEulerInt::m_print_flag.
Referenced by BEulerInt::step().
|
protected |
Time step control function for the selection of the time step size based on a desired accuracy of time integration and on an estimate of the relative error of the time integration process. This routine can be called for a first order - forward Euler/backward Euler predictor/ corrector and for a second order Adams- Bashforth/Trapezoidal Rule predictor/corrector. See Nachos documentation Sand86-1816 and Gresho, Lee, Sani LLNL report UCRL - 83282 for more information.
on input:
order - indicates order of method = 1 -> first order forward Euler/backward Euler predictor/corrector = 2 -> second order forward Adams-Bashforth/Trapezoidal rule predictor/corrector
delta_t_n - Magnitude of time step at time t_n delta_t_nm1 - Magnitude of time step at time t_n-1 rel_error - Generic relative error tolerance time_error_factor - Estimated value of the time step truncation error factor. This value is a ratio of the computed error norms. The premultiplying constants and the power are not yet applied to normalize the predictor/corrector ratio. (see output value)
on output:
return - delta_t for the next time step If delta_t is negative, then the current time step is rejected because the time-step truncation error is too large. The return value will contain the negative of the recommended next time step.
time_error_factor - This output value is normalized so that values greater than one indicate the current time integration error is greater than the user specified magnitude.
Definition at line 693 of file BEulerInt.cpp.
References BEulerInt::m_print_flag.
Referenced by BEulerInt::step().
|
protected |
Solve a nonlinear system.
Find the solution to F(X, xprime) = 0 by damped Newton iteration. On entry, y_comm[] contains an initial estimate of the solution and ydot_comm[] contains an estimate of the derivative. On successful return, y_comm[] contains the converged solution and ydot_comm[] contains the derivative
y_comm[] | Contains the input solution. On output y_comm[] contains the converged solution |
ydot_comm | Contains the input derivative solution. On output y_comm[] contains the converged derivative solution |
CJ | Inverse of the time step |
time_curr | Current value of the time |
jac | Jacobian |
num_newt_its | number of Newton iterations |
num_linear_solves | number of linear solves |
num_backtracks | number of backtracs |
loglevel | Log level |
Definition at line 1731 of file BEulerInt.cpp.
References BEulerInt::beuler_jac(), BEulerInt::calc_ydot(), BEulerInt::dampStep(), BEulerInt::doNewtonSolve(), BEulerInt::filterNewStep(), BEulerInt::m_min_newt_its, BEulerInt::m_neq, BEulerInt::m_numTotalLinearSolves, BEulerInt::m_numTotalNewtIts, and BEulerInt::m_order.
Referenced by BEulerInt::step().
|
protected |
Compute the undamped Newton step.
The residual function is evaluated at the current time, t_n, at the current values of the solution vector, m_y_n, and the solution time derivative, m_ydot_n, but the Jacobian is not recomputed.
Definition at line 1349 of file BEulerInt.cpp.
References Cantera::Base_ResidEval, GeneralMatrix::begin(), ResidJacEval::evalResidNJ(), GeneralMatrix::factored(), BEulerInt::m_colScaling, BEulerInt::m_matrixConditioning, BEulerInt::m_neq, BEulerInt::m_nfe, BEulerInt::m_numTotalLinearSolves, BEulerInt::m_rowScaling, BEulerInt::m_time_step_num, ResidJacEval::matrixConditioning(), ResidJacEval::nEquations(), BEulerInt::setColumnScales(), and GeneralMatrix::solve().
Referenced by BEulerInt::dampStep(), and BEulerInt::solve_nonlinear_problem().
|
protected |
Bound the Newton step while relaxing the solution.
Return the factor by which the undamped Newton step 'step0' must be multiplied in order to keep all solution components in all domains between their specified lower and upper bounds. Other bounds may be applied here as well.
Currently the bounds are hard coded into this routine:
Minimum value for all variables: - 0.01 * m_ewt[i] Maximum value = none.
Thus, this means that all solution components are expected to be numerical greater than zero in the limit of time step truncation errors going to zero.
Delta bounds: The idea behind these is that the Jacobian couldn't possibly be representative if the variable is changed by a lot. (true for nonlinear systems, false for linear systems) Maximum increase in variable in any one Newton iteration: factor of 2 Maximum decrease in variable in any one Newton iteration: factor of 5
y | Current value of the solution |
step0 | Current raw step change in y[] |
loglevel | Log level. This routine produces output if loglevel is greater than one |
Definition at line 1515 of file BEulerInt.cpp.
References BEulerInt::m_ewt, and BEulerInt::m_neq.
Referenced by BEulerInt::dampStep().
|
protected |
On entry, step0 must contain an undamped Newton step for the solution x0. This method attempts to find a damping coefficient such that the next undamped step would have a norm smaller than that of step0. If successful, the new solution after taking the damped step is returned in y1, and the undamped step at y1 is returned in step1.
Definition at line 1582 of file BEulerInt.cpp.
References BEulerInt::boundStep(), BEulerInt::calc_ydot(), Cantera::checkFinite(), Cantera::DampFactor, BEulerInt::doNewtonSolve(), BEulerInt::m_neq, BEulerInt::m_order, Cantera::NDAMP, and BEulerInt::soln_error_norm().
Referenced by BEulerInt::solve_nonlinear_problem().
|
protected |
Compute Residual Weights.
Definition at line 271 of file BEulerInt.cpp.
References GeneralMatrix::begin(), BEulerInt::m_ewt, and BEulerInt::m_neq.
|
protected |
Filter a new step.
Definition at line 295 of file BEulerInt.cpp.
Referenced by BEulerInt::solve_nonlinear_problem(), and BEulerInt::step().
|
protected |
Get the next time to print out.
Definition at line 210 of file BEulerInt.cpp.
References BEulerInt::m_printSolnNumberToTout, BEulerInt::m_t0, and BEulerInt::m_time_final.
|
protected |
IterType is used to specify how the nonlinear equations are to be relaxed at each time step.
Definition at line 394 of file BEulerInt.h.
Referenced by BEulerInt::setIterator().
|
protected |
MethodType is used to specify how the time step is to be chosen.
Currently, there are two choices, one is a fixed step method while the other is based on a predictor-corrector algorithm and a time-step truncation error tolerance.
Definition at line 401 of file BEulerInt.h.
Referenced by BEulerInt::step().
|
protected |
m_jacFormMethod determines how a matrix is formed.
Definition at line 405 of file BEulerInt.h.
Referenced by BEulerInt::beuler_jac(), and BEulerInt::setProblemType().
|
protected |
m_rowScaling is a boolean.
If true then row sum scaling of the Jacobian matrix is carried out when solving the linear systems.
Definition at line 411 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve(), BEulerInt::internalMalloc(), and BEulerInt::setNonLinOptions().
|
protected |
m_colScaling is a boolean.
If true, then column scaling is performed on each solution of the linear system.
Definition at line 416 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve(), BEulerInt::internalMalloc(), and BEulerInt::setNonLinOptions().
|
protected |
m_matrixConditioning is a boolean.
If true, then the Jacobian and every RHS is multiplied by the inverse of a matrix that is suppose to reduce the condition number of the matrix. This is done before row scaling.
Definition at line 423 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::setNonLinOptions().
|
protected |
If m_itol =1 then each component has an individual value of atol.
If m_itol = 0, the all atols are equal.
Definition at line 428 of file BEulerInt.h.
Referenced by BEulerInt::setSolnWeights(), and BEulerInt::setTolerances().
|
protected |
Relative time truncation error tolerances.
Definition at line 431 of file BEulerInt.h.
Referenced by BEulerInt::setSolnWeights(), and BEulerInt::setTolerances().
|
protected |
Absolute time truncation error tolerances, when uniform for all variables.
Definition at line 437 of file BEulerInt.h.
Referenced by BEulerInt::setSolnWeights(), and BEulerInt::setTolerances().
|
protected |
Vector of absolute time truncation error tolerance when not uniform for all variables.
Definition at line 442 of file BEulerInt.h.
Referenced by BEulerInt::setSolnWeights(), and BEulerInt::setTolerances().
|
protected |
Error Weights. This is a surprisingly important quantity.
Definition at line 445 of file BEulerInt.h.
Referenced by BEulerInt::beuler_jac(), BEulerInt::boundStep(), BEulerInt::computeResidWts(), BEulerInt::internalMalloc(), BEulerInt::setSolnWeights(), BEulerInt::soln_error_norm(), and BEulerInt::time_error_norm().
|
protected |
Maximum step size.
Definition at line 448 of file BEulerInt.h.
|
protected |
Maximum integration order.
Definition at line 451 of file BEulerInt.h.
|
protected |
Current integration order.
Definition at line 454 of file BEulerInt.h.
Referenced by BEulerInt::dampStep(), BEulerInt::solve_nonlinear_problem(), and BEulerInt::step().
|
protected |
Time step number.
Definition at line 457 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::step().
|
protected |
Max time steps allowed.
Definition at line 461 of file BEulerInt.h.
|
protected |
Number of initial time steps to take where the time truncation error tolerances are not checked.
Instead the delta T is uniform
Definition at line 466 of file BEulerInt.h.
Referenced by BEulerInt::step().
|
protected |
Failure Counter -> keeps track of the number of consecutive failures.
Definition at line 469 of file BEulerInt.h.
Referenced by BEulerInt::step().
|
protected |
Minimum Number of Newton Iterations per nonlinear step. default = 0.
Definition at line 472 of file BEulerInt.h.
Referenced by BEulerInt::setNonLinOptions(), and BEulerInt::solve_nonlinear_problem().
|
protected |
Step Interval at which to print out the solution default = 1; If set to zero, there is no printout.
Definition at line 481 of file BEulerInt.h.
Referenced by BEulerInt::setPrintSolnOptions().
|
protected |
Number of evenly spaced printouts of the solution If zero, there is no printout from this option default 1 If set to zero there is no printout.
Definition at line 488 of file BEulerInt.h.
Referenced by BEulerInt::getPrintTime(), and BEulerInt::setPrintSolnOptions().
|
protected |
Number of initial steps that the solution is printed out. default = 0.
Definition at line 491 of file BEulerInt.h.
Referenced by BEulerInt::setPrintSolnOptions().
|
protected |
Dump Jacobians to disk. default false.
Definition at line 494 of file BEulerInt.h.
Referenced by BEulerInt::setPrintSolnOptions().
|
protected |
Number of equations in the ode integrator.
Definition at line 501 of file BEulerInt.h.
Referenced by BEulerInt::beuler_jac(), BEulerInt::boundStep(), BEulerInt::calc_y_pred(), BEulerInt::calc_ydot(), BEulerInt::computeResidWts(), BEulerInt::dampStep(), BEulerInt::doNewtonSolve(), BEulerInt::initializeRJE(), BEulerInt::internalMalloc(), BEulerInt::nEquations(), BEulerInt::setNonLinOptions(), BEulerInt::setSolnWeights(), BEulerInt::setTolerances(), BEulerInt::soln_error_norm(), BEulerInt::solve_nonlinear_problem(), and BEulerInt::time_error_norm().
|
protected |
Initial time at the start of the integration.
Definition at line 512 of file BEulerInt.h.
Referenced by BEulerInt::getPrintTime(), and BEulerInt::initializeRJE().
|
protected |
|
protected |
Maximum permissible time step.
Definition at line 526 of file BEulerInt.h.
Referenced by BEulerInt::step().
|
protected |
Pointer to the Jacobian representing the time dependent problem.
Definition at line 536 of file BEulerInt.h.
Referenced by BEulerInt::internalMalloc(), and BEulerInt::step().
|
protected |
Determines the level of printing for each time step.
0 -> absolutely nothing is printed for a single time step. 1 -> One line summary per time step 2 -> short description, points of interest 3 -> Lots printed per time step (default)
Definition at line 547 of file BEulerInt.h.
Referenced by BEulerInt::step(), BEulerInt::time_error_norm(), and BEulerInt::time_step_control().
|
protected |
Number of function evaluations.
Definition at line 553 of file BEulerInt.h.
Referenced by BEulerInt::beuler_jac(), BEulerInt::doNewtonSolve(), and BEulerInt::nEvals().
|
protected |
Number of Jacobian Evaluations and factorization steps (they are the same)
Definition at line 559 of file BEulerInt.h.
Referenced by BEulerInt::beuler_jac().
|
protected |
Number of total Newton iterations.
Definition at line 562 of file BEulerInt.h.
Referenced by BEulerInt::solve_nonlinear_problem().
|
protected |
Total number of linear iterations.
Definition at line 565 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::solve_nonlinear_problem().
|
protected |
Total number of convergence failures.
Definition at line 568 of file BEulerInt.h.
Referenced by BEulerInt::step().
|
protected |
Total Number of time truncation error failures.
Definition at line 571 of file BEulerInt.h.
Referenced by BEulerInt::step().