Cantera  2.0
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
BEulerInt Class Reference

Wrapper class for 'beuler' integrator We derive the class from the class Integrator. More...

#include <BEulerInt.h>

Inheritance diagram for BEulerInt:
[legend]
Collaboration diagram for BEulerInt:
[legend]

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
 
ResidJacEvalm_func
 
double * m_rowScales
 
double * m_colScales
 
GeneralMatrixtdjac_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
 

Detailed Description

Wrapper class for 'beuler' integrator We derive the class from the class Integrator.

Definition at line 61 of file BEulerInt.h.

Constructor & Destructor Documentation

BEulerInt ( )

The default constructor doesn't take an argument.

Definition at line 48 of file BEulerInt.cpp.

~BEulerInt ( )
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.

Member Function Documentation

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

Set or reset the number of equations.

Set error tolerances.

Parameters
reltolscalar relative tolerance
nNumber of equations
abstolarray 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.

void setTolerances ( double  reltol,
double  abstol 
)
virtual

Set error tolerances.

Parameters
reltolscalar relative tolerance
abstolscalar absolute tolerance

Reimplemented from Integrator.

Definition at line 143 of file BEulerInt.cpp.

References BEulerInt::m_abstols, BEulerInt::m_itol, and BEulerInt::m_reltol.

void setProblemType ( int  probtype)
virtual

Set the problem type.

Parameters
probtypeType of the problem

Reimplemented from Integrator.

Definition at line 150 of file BEulerInt.cpp.

References BEulerInt::m_jacFormMethod.

double step ( double  tout)
virtual
virtual double& solution ( size_t  k)
inlinevirtual

The current value of the solution of equation k.

Reimplemented from Integrator.

Definition at line 78 of file BEulerInt.h.

double* solution ( )
inlinevirtual

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

Reimplemented from Integrator.

Definition at line 81 of file BEulerInt.h.

int nEquations ( ) const
inlinevirtual

The number of equations.

Reimplemented from Integrator.

Definition at line 84 of file BEulerInt.h.

References BEulerInt::m_neq.

int nEvals ( ) const
virtual

The number of function evaluations.

Reimplemented from Integrator.

Definition at line 349 of file BEulerInt.cpp.

References BEulerInt::m_nfe.

void setIterator ( IterType  t)
virtual

Set the linear iterator.

Reimplemented from Integrator.

Definition at line 211 of file BEulerInt.cpp.

References BEulerInt::m_iter.

double soln_error_norm ( const double * const  delta_y,
bool  printLargest = false 
)
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.

void internalMalloc ( )
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.

void calc_y_pred ( int  order)
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().

void calc_ydot ( int  order,
double *  y_curr,
double *  ydot_curr 
)
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().

double time_error_norm ( )
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().

double time_step_control ( int  m_order,
double  time_error_factor 
)
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().

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 
)
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

Parameters
y_comm[]Contains the input solution. On output y_comm[] contains the converged solution
ydot_commContains the input derivative solution. On output y_comm[] contains the converged derivative solution
CJInverse of the time step
time_currCurrent value of the time
jacJacobian
num_newt_itsnumber of newton iterations
num_linear_solvesnumber of linear solves
num_backtracksnumber of backtracs
loglevelLog 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().

void doNewtonSolve ( double  time_curr,
double *  y_curr,
double *  ydot_curr,
double *  delta_y,
GeneralMatrix jac,
int  loglevel 
)
protected
double boundStep ( const double *const  y,
const double *const  step0,
int  loglevel 
)
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

Parameters
yCurrent value of the solution
step0Current raw step change in y[]
loglevelLog level. This routine produces output if loglevel is greater than one
Returns
Returns the damping coefficient

Definition at line 1910 of file BEulerInt.cpp.

References BEulerInt::m_ewt, BEulerInt::m_neq, ckr::max(), and ckr::min().

virtual void setSensitivityTolerances ( doublereal  reltol,
doublereal  abstol 
)
inlinevirtualinherited

Set the sensitivity error tolerances.

Parameters
reltolscalar relative tolerance
abstolscalar absolute tolerance

Definition at line 94 of file Integrator.h.

virtual void initialize ( doublereal  t0,
FuncEval func 
)
inlinevirtualinherited

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 in CVodeInt.

Definition at line 111 of file Integrator.h.

Referenced by ImplicitSurfChem::initialize(), and ImplicitSurfChem::integrate().

virtual void integrate ( doublereal  tout)
inlinevirtualinherited

Integrate the system of equations.

Parameters
toutIntegrate 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().

virtual void setMaxOrder ( int  n)
inlinevirtualinherited

Set the maximum integration order that will be used.

Reimplemented in CVodeInt.

Definition at line 163 of file Integrator.h.

virtual void setMethod ( MethodType  t)
inlinevirtualinherited

Set the solution method.

Reimplemented in CVodeInt.

Definition at line 168 of file Integrator.h.

Referenced by ImplicitSurfChem::ImplicitSurfChem().

virtual void setMaxStepSize ( double  hmax)
inlinevirtualinherited

Set the maximum step size.

Reimplemented in CVodeInt.

Definition at line 178 of file Integrator.h.

Referenced by ImplicitSurfChem::integrate().

virtual void setMinStepSize ( double  hmin)
inlinevirtualinherited

Set the minimum step size.

Reimplemented in CVodeInt.

Definition at line 183 of file Integrator.h.

Member Data Documentation

IterType m_iter
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().

BEulerMethodType m_method
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().

int m_jacFormMethod
protected

m_jacFormMethod determines how a matrix is formed.

Definition at line 265 of file BEulerInt.h.

Referenced by BEulerInt::setProblemType().

bool m_rowScaling
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().

bool m_colScaling
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().

bool m_matrixConditioning
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().

int m_itol
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().

double m_reltol
protected

Relative time truncation error tolerances.

Definition at line 292 of file BEulerInt.h.

Referenced by BEulerInt::setTolerances().

double m_abstols
protected

Absolute time truncation error tolerances, when uniform for all variables.

Definition at line 297 of file BEulerInt.h.

Referenced by BEulerInt::setTolerances().

double* m_abstol
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().

double* m_ewt
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().

double m_hmax
protected

Maximum step size.

Definition at line 309 of file BEulerInt.h.

int m_maxord
protected

Maximum integration order.

Definition at line 313 of file BEulerInt.h.

int m_order
protected

Current integration order.

Definition at line 317 of file BEulerInt.h.

Referenced by BEulerInt::solve_nonlinear_problem(), and BEulerInt::step().

int m_time_step_num
protected

Time step number.

Definition at line 321 of file BEulerInt.h.

Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::step().

int m_max_time_step_attempts
protected

Max time steps allowed.

Definition at line 326 of file BEulerInt.h.

int m_numInitialConstantDeltaTSteps
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().

int m_failure_counter
protected

Failure Counter -> keeps track of the number of consequetive failures.

Definition at line 337 of file BEulerInt.h.

Referenced by BEulerInt::step().

int m_min_newt_its
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().

int m_printSolnStepInterval
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.

int m_printSolnNumberToTout
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.

int m_printSolnFirstSteps
protected

Number of initial steps that the solution is printed out.

default = 0

Definition at line 365 of file BEulerInt.h.

bool m_dumpJacobians
protected

Dump Jacobians to disk default false.

Definition at line 371 of file BEulerInt.h.

int m_neq
protected
double m_t0
protected

Initial time at the start of the integration.

Definition at line 391 of file BEulerInt.h.

double m_time_final
protected

Final time.

Definition at line 395 of file BEulerInt.h.

double delta_t_max
protected

Maximum permissible time step.

Definition at line 409 of file BEulerInt.h.

Referenced by BEulerInt::step().

GeneralMatrix* tdjac_ptr
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().

int m_print_flag
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().

int m_nfe
protected

Number of function evaluations.

Definition at line 440 of file BEulerInt.h.

Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::nEvals().

int m_nJacEval
protected

Number of Jacobian Evaluations and factorization steps (they are the same)

Definition at line 445 of file BEulerInt.h.

int m_numTotalNewtIts
protected

Number of total newton iterations.

Definition at line 449 of file BEulerInt.h.

Referenced by BEulerInt::solve_nonlinear_problem().

int m_numTotalLinearSolves
protected

Total number of linear iterations.

Definition at line 453 of file BEulerInt.h.

Referenced by BEulerInt::doNewtonSolve(), and BEulerInt::solve_nonlinear_problem().

int m_numTotalConvFails
protected

Total number of convergence failures.

Definition at line 457 of file BEulerInt.h.

Referenced by BEulerInt::step().

int m_numTotalTruncFails
protected

Total Number of time truncation error failures.

Definition at line 461 of file BEulerInt.h.

Referenced by BEulerInt::step().


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