Cantera
2.0
|
Wrapper class for 'beuler' integrator We derive the class from the class Integrator. More...
#include <BEulerInt.h>
Public Member Functions | |
BEulerInt () | |
The default constructor doesn't take an argument. | |
virtual | ~BEulerInt () |
Destructor. | |
virtual void | setTolerances (double reltol, size_t n, double *abstol) |
Set or reset the number of equations. | |
virtual void | setTolerances (double reltol, double abstol) |
Set error tolerances. | |
virtual void | setProblemType (int probtype) |
Set the problem type. | |
virtual void | initializeRJE (double t0, ResidJacEval &func) |
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. | |
virtual void | setSolnWeights () |
virtual double & | solution (size_t k) |
The current value of the solution of equation k. | |
double * | solution () |
The current value of the solution of the system of equations. | |
int | nEquations () const |
The number of equations. | |
virtual int | nEvals () const |
The number of function evaluations. | |
virtual void | setMethodBEMT (BEulerMethodType t) |
virtual void | setIterator (IterType t) |
Set the linear iterator. | |
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) |
void | setNonLinOptions (int min_newt_its=0, bool matrixConditioning=false, bool colScaling=false, bool rowScaling=true) |
virtual void | setPrintFlag (int print_flag) |
virtual void | setColumnScales () |
virtual double | soln_error_norm (const double *const, bool printLargest=false) |
calculate the solution error norm | |
virtual void | setInitialTimeStep (double delta_t) |
void | beuler_jac (GeneralMatrix &, double *const, double, double, double *const, double *const, int) |
virtual void | setSensitivityTolerances (doublereal reltol, doublereal abstol) |
Set the sensitivity error tolerances. | |
virtual void | initialize (doublereal t0, FuncEval &func) |
Initialize the integrator for a new problem. | |
virtual void | reinitialize (doublereal t0, FuncEval &func) |
virtual void | integrate (doublereal tout) |
Integrate the system of equations. | |
virtual void | setMaxOrder (int n) |
Set the maximum integration order that will be used. | |
virtual void | setMethod (MethodType t) |
Set the solution method. | |
virtual void | setMaxStepSize (double hmax) |
Set the maximum step size. | |
virtual void | setMinStepSize (double hmin) |
Set the minimum step size. | |
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. | |
void | calc_y_pred (int) |
Internal function to calculate the predicted solution at a time step. | |
void | calc_ydot (int, double *, double *) |
Internal function to calculate the time derivative at the new step. | |
double | time_error_norm () |
Internal function to calculate the time step truncation error for a predictor corrector time step. | |
double | time_step_control (int m_order, double time_error_factor) |
Internal function to calculate the time step for the next step based on the time-truncation error on the current time step. | |
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. | |
void | doNewtonSolve (double, double *, double *, double *, GeneralMatrix &, int) |
Compute the undamped Newton step. | |
double | boundStep (const double *const y, const double *const step0, int loglevel) |
Bound the Newton step while relaxing the solution. | |
int | dampStep (double, const double *, const double *, const double *, double *, double *, double *, double &, GeneralMatrix &, int &, bool, int &) |
void | computeResidWts (GeneralMatrix &jac) |
double | filterNewStep (double, double *, double *) |
double | getPrintTime (double time_current) |
Protected Attributes | |
IterType | m_iter |
IterType is used to specify how the nonlinear equations are to be relaxed at each time step. | |
BEulerMethodType | m_method |
MethodType is used to specify how the time step is to be chosen. | |
int | m_jacFormMethod |
m_jacFormMethod determines how a matrix is formed. | |
bool | m_rowScaling |
m_rowScaling is a boolean. | |
bool | m_colScaling |
m_colScaling is a boolean. | |
bool | m_matrixConditioning |
m_matrixConditioning is a boolean. | |
int | m_itol |
If m_itol =1 then each component has an individual value of atol. | |
double | m_reltol |
Relative time truncation error tolerances. | |
double | m_abstols |
Absolute time truncation error tolerances, when uniform for all variables. | |
double * | m_abstol |
Vector of absolute time truncation error tolerance when not uniform for all variables. | |
double * | m_ewt |
Error Weights. | |
double | m_hmax |
Maximum step size. | |
int | m_maxord |
Maximum integration order. | |
int | m_order |
Current integration order. | |
int | m_time_step_num |
Time step number. | |
int | m_time_step_attempts |
int | m_max_time_step_attempts |
Max time steps allowed. | |
int | m_numInitialConstantDeltaTSteps |
Number of initial time steps to take where the time truncation error tolerances are not checked. | |
int | m_failure_counter |
Failure Counter -> keeps track of the number of consequetive failures. | |
int | m_min_newt_its |
Minimum Number of Newton Iterations per nonlinear step default = 0. | |
int | m_printSolnStepInterval |
Step Interval at which to print out the solution default = 1; If set to zero, there is no printout. | |
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. | |
int | m_printSolnFirstSteps |
Number of initial steps that the solution is printed out. | |
bool | m_dumpJacobians |
Dump Jacobians to disk default false. | |
int | m_neq |
Number of equations in the ode integrator. | |
double * | m_y_n |
double * | m_y_nm1 |
double * | m_y_pred_n |
double * | m_ydot_n |
double * | m_ydot_nm1 |
double | m_t0 |
Initial time at the start of the integration. | |
double | m_time_final |
Final time. | |
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. | |
double * | m_resid |
double * | m_residWts |
double * | m_wksp |
ResidJacEval * | m_func |
double * | m_rowScales |
double * | m_colScales |
GeneralMatrix * | tdjac_ptr |
Pointer to the jacobian representing the time dependent problem. | |
int | m_print_flag |
Determines the level of printing for each time step. | |
int | m_nfe |
Number of function evaluations. | |
int | m_nJacEval |
Number of Jacobian Evaluations and factorization steps (they are the same) | |
int | m_numTotalNewtIts |
Number of total newton iterations. | |
int | m_numTotalLinearSolves |
Total number of linear iterations. | |
int | m_numTotalConvFails |
Total number of convergence failures. | |
int | m_numTotalTruncFails |
Total Number of time truncation error failures. | |
int | num_failures |
Wrapper class for 'beuler' integrator We derive the class from the class Integrator.
Definition at line 61 of file BEulerInt.h.
BEulerInt | ( | ) |
The default constructor doesn't take an argument.
Definition at line 48 of file BEulerInt.cpp.
|
virtual |
Destructor.
Definition at line 110 of file BEulerInt.cpp.
References BEulerInt::m_abstol, BEulerInt::m_ewt, mdp::mdp_safe_free(), and BEulerInt::tdjac_ptr.
|
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 127 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 143 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 150 of file BEulerInt.cpp.
References BEulerInt::m_jacFormMethod.
|
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 1231 of file BEulerInt.cpp.
References BEulerInt::calc_y_pred(), BEulerInt::calc_ydot(), BEulerInt::delta_t_max, 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, ckr::max(), mdp::mdp_copy_dbl_1(), ckr::min(), BEulerInt::solve_nonlinear_problem(), BEulerInt::tdjac_ptr, BEulerInt::time_error_norm(), and BEulerInt::time_step_control().
|
inlinevirtual |
The current value of the solution of equation k.
Reimplemented from Integrator.
Definition at line 78 of file BEulerInt.h.
|
inlinevirtual |
The current value of the solution of the system of equations.
Reimplemented from Integrator.
Definition at line 81 of file BEulerInt.h.
|
inlinevirtual |
The number of equations.
Reimplemented from Integrator.
Definition at line 84 of file BEulerInt.h.
References BEulerInt::m_neq.
|
virtual |
The number of function evaluations.
Reimplemented from Integrator.
Definition at line 349 of file BEulerInt.cpp.
References BEulerInt::m_nfe.
|
virtual |
Set the linear iterator.
Reimplemented from Integrator.
Definition at line 211 of file BEulerInt.cpp.
References BEulerInt::m_iter.
|
virtual |
calculate the solution error norm
L2 Norm of a delta in the solution.
The second argument has a default of false. However, if true, then a table of the largest values is printed out to standard output.
Definition at line 1648 of file BEulerInt.cpp.
References Cantera::error(), BEulerInt::m_ewt, and BEulerInt::m_neq.
|
protected |
Internal routine that sets up the fixed length storage based on the size of the problem to solve.
Definition at line 361 of file BEulerInt.cpp.
References BEulerInt::m_colScaling, BEulerInt::m_ewt, BEulerInt::m_neq, BEulerInt::m_rowScaling, and BEulerInt::tdjac_ptr.
|
protected |
Internal function to calculate the predicted solution at a time step.
Definition at line 817 of file BEulerInt.cpp.
References ResidJacEval::filterSolnPrediction(), and BEulerInt::m_neq.
Referenced by BEulerInt::step().
|
protected |
Internal function to calculate the time derivative at the new step.
Definition at line 880 of file BEulerInt.cpp.
References BEulerInt::m_neq.
Referenced by BEulerInt::solve_nonlinear_problem(), and BEulerInt::step().
|
protected |
Internal function to calculate the time step truncation error for a predictor corrector time step.
Definition at line 925 of file BEulerInt.cpp.
References Cantera::error(), BEulerInt::m_ewt, BEulerInt::m_neq, and BEulerInt::m_print_flag.
Referenced by BEulerInt::step().
|
protected |
Internal function to calculate the time step for the next step based on the time-truncation error on the current time step.
Definition at line 1021 of file BEulerInt.cpp.
References BEulerInt::m_print_flag, ckr::max(), and ckr::min().
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 2165 of file BEulerInt.cpp.
References BEulerInt::calc_ydot(), BEulerInt::doNewtonSolve(), BEulerInt::m_min_newt_its, BEulerInt::m_numTotalLinearSolves, BEulerInt::m_numTotalNewtIts, BEulerInt::m_order, and mdp::mdp_copy_dbl_1().
Referenced by BEulerInt::step().
|
protected |
Compute the undamped Newton step.
The residual function is evaluated at x, but the Jacobian is not recomputed.
Definition at line 1711 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(), and GeneralMatrix::solve().
Referenced by 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 1910 of file BEulerInt.cpp.
References BEulerInt::m_ewt, BEulerInt::m_neq, ckr::max(), and ckr::min().
|
inlinevirtualinherited |
Set the sensitivity error tolerances.
reltol | scalar relative tolerance |
abstol | scalar absolute tolerance |
Definition at line 94 of file Integrator.h.
|
inlinevirtualinherited |
Initialize the integrator for a new problem.
Call after all options have been set.
t0 | initial time |
func | RHS evaluator object for system of equations. |
Reimplemented in CVodeInt.
Definition at line 111 of file Integrator.h.
Referenced by ImplicitSurfChem::initialize(), and ImplicitSurfChem::integrate().
|
inlinevirtualinherited |
Integrate the system of equations.
tout | Integrate to this time. Note that this is the absolute time value, not a time interval. |
Reimplemented in CVodeInt.
Definition at line 124 of file Integrator.h.
Referenced by ImplicitSurfChem::integrate(), and ImplicitSurfChem::integrate0().
|
inlinevirtualinherited |
Set the maximum integration order that will be used.
Reimplemented in CVodeInt.
Definition at line 163 of file Integrator.h.
|
inlinevirtualinherited |
Set the solution method.
Reimplemented in CVodeInt.
Definition at line 168 of file Integrator.h.
Referenced by ImplicitSurfChem::ImplicitSurfChem().
|
inlinevirtualinherited |
Set the maximum step size.
Reimplemented in CVodeInt.
Definition at line 178 of file Integrator.h.
Referenced by ImplicitSurfChem::integrate().
|
inlinevirtualinherited |
|
protected |
IterType is used to specify how the nonlinear equations are to be relaxed at each time step.
Definition at line 254 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 261 of file BEulerInt.h.
Referenced by BEulerInt::step().
|
protected |
m_jacFormMethod determines how a matrix is formed.
Definition at line 265 of file BEulerInt.h.
Referenced by 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 271 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::internalMalloc().
|
protected |
m_colScaling is a boolean.
If true, then column scaling is performed on each solution of the linear system.
Definition at line 276 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::internalMalloc().
|
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 283 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve().
|
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 288 of file BEulerInt.h.
Referenced by BEulerInt::setTolerances().
|
protected |
Relative time truncation error tolerances.
Definition at line 292 of file BEulerInt.h.
Referenced by BEulerInt::setTolerances().
|
protected |
Absolute time truncation error tolerances, when uniform for all variables.
Definition at line 297 of file BEulerInt.h.
Referenced by BEulerInt::setTolerances().
|
protected |
Vector of absolute time truncation error tolerance when not uniform for all variables.
Definition at line 302 of file BEulerInt.h.
Referenced by BEulerInt::setTolerances(), and BEulerInt::~BEulerInt().
|
protected |
Error Weights.
This is a surprisingly important quantity.
Definition at line 306 of file BEulerInt.h.
Referenced by BEulerInt::boundStep(), BEulerInt::internalMalloc(), BEulerInt::soln_error_norm(), BEulerInt::time_error_norm(), and BEulerInt::~BEulerInt().
|
protected |
Maximum step size.
Definition at line 309 of file BEulerInt.h.
|
protected |
Maximum integration order.
Definition at line 313 of file BEulerInt.h.
|
protected |
Current integration order.
Definition at line 317 of file BEulerInt.h.
Referenced by BEulerInt::solve_nonlinear_problem(), and BEulerInt::step().
|
protected |
Time step number.
Definition at line 321 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::step().
|
protected |
Max time steps allowed.
Definition at line 326 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 332 of file BEulerInt.h.
Referenced by BEulerInt::step().
|
protected |
Failure Counter -> keeps track of the number of consequetive failures.
Definition at line 337 of file BEulerInt.h.
Referenced by BEulerInt::step().
|
protected |
Minimum Number of Newton Iterations per nonlinear step default = 0.
Definition at line 342 of file BEulerInt.h.
Referenced by 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 351 of file BEulerInt.h.
|
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 358 of file BEulerInt.h.
|
protected |
Number of initial steps that the solution is printed out.
default = 0
Definition at line 365 of file BEulerInt.h.
|
protected |
Dump Jacobians to disk default false.
Definition at line 371 of file BEulerInt.h.
|
protected |
Number of equations in the ode integrator.
Definition at line 379 of file BEulerInt.h.
Referenced by BEulerInt::boundStep(), BEulerInt::calc_y_pred(), BEulerInt::calc_ydot(), BEulerInt::doNewtonSolve(), BEulerInt::internalMalloc(), BEulerInt::nEquations(), BEulerInt::setTolerances(), BEulerInt::soln_error_norm(), and BEulerInt::time_error_norm().
|
protected |
Initial time at the start of the integration.
Definition at line 391 of file BEulerInt.h.
|
protected |
Final time.
Definition at line 395 of file BEulerInt.h.
|
protected |
Maximum permissible time step.
Definition at line 409 of file BEulerInt.h.
Referenced by BEulerInt::step().
|
protected |
Pointer to the jacobian representing the time dependent problem.
Definition at line 423 of file BEulerInt.h.
Referenced by BEulerInt::internalMalloc(), BEulerInt::step(), and BEulerInt::~BEulerInt().
|
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 433 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 440 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::nEvals().
|
protected |
Number of Jacobian Evaluations and factorization steps (they are the same)
Definition at line 445 of file BEulerInt.h.
|
protected |
Number of total newton iterations.
Definition at line 449 of file BEulerInt.h.
Referenced by BEulerInt::solve_nonlinear_problem().
|
protected |
Total number of linear iterations.
Definition at line 453 of file BEulerInt.h.
Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::solve_nonlinear_problem().
|
protected |
Total number of convergence failures.
Definition at line 457 of file BEulerInt.h.
Referenced by BEulerInt::step().
|
protected |
Total Number of time truncation error failures.
Definition at line 461 of file BEulerInt.h.
Referenced by BEulerInt::step().