Cantera
2.2.1
|
Class that calculates the solution to a nonlinear system. More...
#include <NonlinearSolver.h>
Public Member Functions | |
NonlinearSolver (ResidJacEval *func) | |
Default constructor. More... | |
NonlinearSolver (const NonlinearSolver &right) | |
Copy Constructor. More... | |
virtual | ~NonlinearSolver () |
Destructor. More... | |
NonlinearSolver & | operator= (const NonlinearSolver &right) |
Assignment operator. More... | |
void | createSolnWeights (const doublereal *const y) |
Create solution weights for convergence criteria. More... | |
doublereal | solnErrorNorm (const doublereal *const delta_y, const char *title=0, int printLargest=0, const doublereal dampFactor=1.0) const |
L2 norm of the delta of the solution vector. More... | |
doublereal | residErrorNorm (const doublereal *const resid, const char *title=0, const int printLargest=0, const doublereal *const y=0) const |
L2 norm of the residual of the equation system. More... | |
int | doResidualCalc (const doublereal time_curr, const int typeCalc, const doublereal *const y_curr, const doublereal *const ydot_curr, const ResidEval_Type_Enum evalType=Base_ResidEval) const |
Compute the current residual. More... | |
int | doNewtonSolve (const doublereal time_curr, const doublereal *const y_curr, const doublereal *const ydot_curr, doublereal *const delta_y, GeneralMatrix &jac) |
Compute the undamped Newton step. More... | |
int | doAffineNewtonSolve (const doublereal *const y_curr, const doublereal *const ydot_curr, doublereal *const delta_y, GeneralMatrix &jac) |
Compute the Newton step, either by direct Newton's or by solving a close problem that is represented by a Hessian. More... | |
doublereal | trustRegionLength () const |
Calculate the length of the current trust region in terms of the solution error norm. More... | |
void | setDefaultDeltaBoundsMagnitudes () |
Set default deulta bounds amounts. More... | |
void | adjustUpStepMinimums () |
Adjust the step minimums. More... | |
void | setDeltaBoundsMagnitudes (const doublereal *const deltaBoundsMagnitudes) |
Set the delta Bounds magnitudes by hand. More... | |
doublereal | boundStep (const doublereal *const y, const doublereal *const step0) |
Bound the step. More... | |
void | setBoundsConstraints (const doublereal *const y_low_bounds, const doublereal *const y_high_bounds) |
Set bounds constraints for all variables in the problem. More... | |
std::vector< doublereal > & | lowBoundsConstraintVector () |
Return an editable vector of the low bounds constraints. More... | |
std::vector< doublereal > & | highBoundsConstraintVector () |
Return an editable vector of the high bounds constraints. More... | |
void | calc_ydot (const int order, const doublereal *const y_curr, doublereal *const ydot_curr) const |
Internal function to calculate the time derivative of the solution at the new step. More... | |
int | beuler_jac (GeneralMatrix &J, doublereal *const f, doublereal time_curr, doublereal CJ, doublereal *const y, doublereal *const ydot, int num_newt_its) |
Function called to evaluate the Jacobian matrix and the current residual vector at the current time step. More... | |
doublereal | filterNewStep (const doublereal timeCurrent, const doublereal *const ybase, doublereal *const step0) |
Apply a filtering process to the step. More... | |
doublereal | filterNewSolution (const doublereal timeCurrent, doublereal *const y_current, doublereal *const ydot_current) |
Apply a filter to the solution. More... | |
doublereal | deltaBoundStep (const doublereal *const y, const doublereal *const step0) |
Return the factor by which the undamped Newton step 'step0' must be multiplied in order to keep the update within the bounds of an accurate Jacobian. More... | |
int | dampStep (const doublereal time_curr, const doublereal *const y_n_curr, const doublereal *const ydot_n_curr, doublereal *const step_1, doublereal *const y_n_1, doublereal *const ydot_n_1, doublereal *step_2, doublereal &stepNorm_2, GeneralMatrix &jac, bool writetitle, int &num_backtracks) |
Find a damping coefficient through a look-ahead mechanism. More... | |
int | solve_nonlinear_problem (int SolnType, doublereal *const y_comm, doublereal *const ydot_comm, doublereal CJ, doublereal time_curr, GeneralMatrix &jac, int &num_newt_its, int &num_linear_solves, int &num_backtracks, int loglevelInput) |
Find the solution to F(X) = 0 by damped Newton iteration. More... | |
virtual void | setPreviousTimeStep (const std::vector< doublereal > &y_nm1, const std::vector< doublereal > &ydot_nm1) |
Set the values for the previous time step. More... | |
void | setColumnScaling (bool useColScaling, const double *const scaleFactors=0) |
Set the column scaling that are used for the inversion of the matrix. More... | |
void | setRowScaling (bool useRowScaling) |
Set the rowscaling that are used for the inversion of the matrix. More... | |
void | scaleMatrix (GeneralMatrix &jac, doublereal *const y_comm, doublereal *const ydot_comm, doublereal time_curr, int num_newt_its) |
Scale the matrix. More... | |
void | print_solnDelta_norm_contrib (const doublereal *const step_1, const char *const stepNorm_1, const doublereal *const step_2, const char *const stepNorm_2, const char *const title, const doublereal *const y_n_curr, const doublereal *const y_n_1, doublereal damp, size_t num_entries) |
Print solution norm contribution. More... | |
void | computeResidWts () |
Compute the Residual Weights. More... | |
void | getResidWts (doublereal *const residWts) const |
Return the residual weights. More... | |
int | convergenceCheck (int dampCode, doublereal s1) |
Check to see if the nonlinear problem has converged. More... | |
void | setAtol (const doublereal *const atol) |
Set the absolute tolerances for the solution variables. More... | |
void | setRtol (const doublereal rtol) |
Set the relative tolerances for the solution variables. More... | |
void | setResidualTols (double residRtol, double *residAtol, int residNormHandling=2) |
Set the relative and absolute tolerances for the Residual norm comparisons, if used. More... | |
void | setMaxNewtIts (const int maxNewtIts) |
Set the value of the maximum # of Newton iterations. More... | |
void | calcSolnToResNormVector () |
Calculate the scaling factor for translating residual norms into solution norms. More... | |
doublereal | doCauchyPointSolve (GeneralMatrix &jac) |
Calculate the steepest descent direction and the Cauchy Point where the quadratic formulation of the nonlinear problem expects a minimum along the descent direction. More... | |
void | descentComparison (doublereal time_curr, doublereal *ydot0, doublereal *ydot1, int &numTrials) |
This is a utility routine that can be used to print out the rates of the initial residual decline. More... | |
void | setupDoubleDogleg () |
Setup the parameters for the double dog leg. More... | |
int | lambdaToLeg (const doublereal lambda, doublereal &alpha) const |
Change the global lambda coordinate into the (leg,alpha) coordinate for the double dogleg. More... | |
int | calcTrustIntersection (doublereal trustVal, doublereal &lambda, doublereal &alpha) const |
Given a trust distance, this routine calculates the intersection of the this distance with the double dogleg curve. More... | |
void | initializeTrustRegion () |
Initialize the size of the trust vector. More... | |
void | setTrustRegionInitializationMethod (int method, doublereal factor) |
Set Trust region initialization strategy. More... | |
int | dampDogLeg (const doublereal time_curr, const doublereal *y_n_curr, const doublereal *ydot_n_curr, std::vector< doublereal > &step_1, doublereal *const y_n_1, doublereal *const ydot_n_1, doublereal &stepNorm_1, doublereal &stepNorm_2, GeneralMatrix &jac, int &num_backtracks) |
Damp using the dog leg approach. More... | |
int | decideStep (const doublereal time_curr, int leg, doublereal alpha, const doublereal *const y_n_curr, const doublereal *const ydot_n_curr, const std::vector< doublereal > &step_1, const doublereal *const y_n_1, const doublereal *const ydot_n_1, doublereal trustDeltaOld) |
Decide whether the current step is acceptable and adjust the trust region size. More... | |
doublereal | expectedResidLeg (int leg, doublereal alpha) const |
Calculated the expected residual along the double dogleg curve. More... | |
void | residualComparisonLeg (const doublereal time_curr, const doublereal *const ydot0, int &legBest, doublereal &alphaBest) const |
Here we print out the residual at various points along the double dogleg, comparing against the quadratic model in a table format. More... | |
void | setPrintLvl (int printLvl) |
Set the print level from the nonlinear solver. More... | |
void | setSolverScheme (int doDogLeg, int doAffineSolve) |
Parameter to turn on solution solver schemes. More... | |
Public Attributes | |
int | m_min_newt_its |
Minimum number of Newton iterations to use. More... | |
Static Public Attributes | |
static bool | s_TurnOffTiming |
Turn off printing of time. More... | |
static bool | s_print_NumJac |
Turn on or off printing of the Jacobian. More... | |
static bool | s_print_DogLeg |
Turn on extra printing of dogleg information. More... | |
static bool | s_doBothSolvesAndCompare |
Turn on solving both the Newton and Hessian systems and comparing the results. More... | |
static bool | s_alwaysAssumeNewtonGood |
This toggle turns off the use of the Hessian when it is warranted by the condition number. More... | |
Protected Member Functions | |
void | readjustTrustVector () |
Readjust the trust region vectors. More... | |
void | fillDogLegStep (int leg, doublereal alpha, std::vector< doublereal > &deltaX) const |
Fill a dogleg solution step vector. More... | |
doublereal | calcTrustDistance (std::vector< doublereal > const &deltaX) const |
Calculate the trust distance of a step in the solution variables. More... | |
Private Member Functions | |
void | calcColumnScales () |
Set the column scaling vector at the current time. More... | |
Private Attributes | |
ResidJacEval * | m_func |
Pointer to the residual and Jacobian evaluator for the function. More... | |
int | solnType_ |
Solution type. More... | |
size_t | neq_ |
Local copy of the number of equations. More... | |
std::vector< doublereal > | m_ewt |
Soln error weights. More... | |
int | m_manualDeltaStepSet |
Boolean indicating whether a manual delta bounds has been input. More... | |
std::vector< doublereal > | m_deltaStepMinimum |
Soln Delta bounds magnitudes. More... | |
std::vector< doublereal > | m_deltaStepMaximum |
Value of the delta step magnitudes. More... | |
std::vector< doublereal > | m_y_n_curr |
Vector containing the current solution vector within the nonlinear solver. More... | |
std::vector< doublereal > | m_ydot_n_curr |
Vector containing the time derivative of the current solution vector within the nonlinear solver (where applicable) More... | |
std::vector< doublereal > | m_y_nm1 |
Vector containing the solution at the previous time step. More... | |
std::vector< doublereal > | m_ydot_nm1 |
Vector containing the solution derivative at the previous time step. More... | |
std::vector< doublereal > | m_y_n_trial |
Vector containing the solution at the new point that is to be considered. More... | |
std::vector< doublereal > | m_ydot_trial |
Value of the solution time derivative at the new point that is to be considered. More... | |
std::vector< doublereal > | m_step_1 |
Value of the step to be taken in the solution. More... | |
std::vector< doublereal > | m_colScales |
Vector of column scaling factors. More... | |
std::vector< doublereal > | m_rowScales |
Weights for normalizing the values of the residuals. More... | |
std::vector< doublereal > | m_rowWtScales |
Weights for normalizing the values of the residuals. More... | |
std::vector< doublereal > | m_resid |
Value of the residual for the nonlinear problem. More... | |
std::vector< doublereal > | m_wksp |
Workspace of length neq_. More... | |
std::vector< doublereal > | m_wksp_2 |
Workspace of length neq_. More... | |
std::vector< doublereal > | m_residWts |
Vector of residual weights. More... | |
doublereal | m_normResid_0 |
Norm of the residual at the start of each nonlinear iteration. More... | |
doublereal | m_normResid_Bound |
Norm of the residual after it has been bounded. More... | |
doublereal | m_normResid_1 |
Norm of the residual at the end of the first leg of the current iteration. More... | |
doublereal | m_normResid_full |
Norm of the residual at the end of the first leg of the current iteration. More... | |
doublereal | m_normDeltaSoln_Newton |
Norm of the solution update created by the iteration in its raw, undamped form, using the solution norm. More... | |
doublereal | m_normDeltaSoln_CP |
Norm of the distance to the Cauchy point using the solution norm. More... | |
doublereal | m_normResidTrial |
Norm of the residual for a trial calculation which may or may not be used. More... | |
doublereal | m_normResidPoints [15] |
Vector of the norm. More... | |
bool | m_resid_scaled |
Boolean indicating whether we should scale the residual. More... | |
std::vector< doublereal > | m_y_high_bounds |
Bounds vector for each species. More... | |
std::vector< doublereal > | m_y_low_bounds |
Lower bounds vector for each species. More... | |
doublereal | m_dampBound |
Damping factor imposed by hard bounds and by delta bounds. More... | |
doublereal | m_dampRes |
Additional damping factor due to bounds on the residual and solution norms. More... | |
doublereal | delta_t_n |
Delta t for the current step. More... | |
int | m_nfe |
Counter for the total number of function evaluations. More... | |
int | m_colScaling |
The type of column scaling used in the matrix inversion of the problem. More... | |
int | m_rowScaling |
int indicating whether row scaling is turned on (1) or not (0) More... | |
int | m_numTotalLinearSolves |
Total number of linear solves taken by the solver object. More... | |
int | m_numLocalLinearSolves |
Number of local linear solves done during the current iteration. More... | |
int | m_numTotalNewtIts |
Total number of Newton iterations. More... | |
int | maxNewtIts_ |
Maximum number of Newton iterations. More... | |
int | m_jacFormMethod |
Jacobian formation method. More... | |
int | m_nJacEval |
Number of Jacobian evaluations. More... | |
doublereal | time_n |
Current system time. More... | |
int | m_matrixConditioning |
Boolean indicating matrix conditioning. More... | |
int | m_order |
Order of the time step method = 1. More... | |
doublereal | rtol_ |
value of the relative tolerance to use in solving the equation set More... | |
doublereal | atolBase_ |
Base value of the absolute tolerance. More... | |
std::vector< doublereal > | atolk_ |
absolute tolerance in the solution unknown More... | |
std::vector< doublereal > | userResidAtol_ |
absolute tolerance in the unscaled solution unknowns More... | |
doublereal | userResidRtol_ |
absolute tolerance in the unscaled solution unknowns More... | |
int | checkUserResidualTols_ |
Check the residual tolerances explicitly against user input. More... | |
int | m_print_flag |
Determines the level of printing for each time step. More... | |
doublereal | m_ScaleSolnNormToResNorm |
Scale factor for turning residual norms into solution norms. More... | |
Cantera::GeneralMatrix * | jacCopyPtr_ |
Copy of the Jacobian that doesn't get overwritten when the inverse is determined. More... | |
Cantera::GeneralMatrix * | HessianPtr_ |
Hessian. More... | |
std::vector< doublereal > | deltaX_CP_ |
Steepest descent direction. This is also the distance to the Cauchy Point. More... | |
std::vector< doublereal > | deltaX_Newton_ |
Newton Step - This is the Newton step determined from the straight Jacobian. More... | |
doublereal | residNorm2Cauchy_ |
Expected value of the residual norm at the Cauchy point if the quadratic model were valid. More... | |
int | dogLegID_ |
Current leg. More... | |
doublereal | dogLegAlpha_ |
Current Alpha param along the leg. More... | |
doublereal | RJd_norm_ |
Residual dot Jd norm. More... | |
doublereal | lambdaStar_ |
Value of lambdaStar_ which is used to calculate the Cauchy point. More... | |
std::vector< doublereal > | Jd_ |
Jacobian times the steepest descent direction in the normalized coordinates. More... | |
std::vector< doublereal > | deltaX_trust_ |
Vector of trust region values. More... | |
doublereal | norm_deltaX_trust_ |
Current norm of the vector deltaX_trust_ in terms of the solution norm. More... | |
doublereal | trustDelta_ |
Current value of trust radius. More... | |
int | trustRegionInitializationMethod_ |
Method for handling the trust region initialization. More... | |
doublereal | trustRegionInitializationFactor_ |
Factor used to set the initial trust region. More... | |
doublereal | Nuu_ |
Relative distance down the Newton step that the second dogleg starts. More... | |
doublereal | dist_R0_ |
Distance of the zeroeth leg of the dogleg in terms of the solution norm. More... | |
doublereal | dist_R1_ |
Distance of the first leg of the dogleg in terms of the solution norm. More... | |
doublereal | dist_R2_ |
Distance of the second leg of the dogleg in terms of the solution norm. More... | |
doublereal | dist_Total_ |
Distance of the sum of all legs of the doglegs in terms of the solution norm. More... | |
doublereal | JdJd_norm_ |
Dot product of the Jd_ variable defined above with itself. More... | |
doublereal | normTrust_Newton_ |
Norm of the Newton Step wrt trust region. More... | |
doublereal | normTrust_CP_ |
Norm of the Cauchy Step direction wrt trust region. More... | |
int | doDogLeg_ |
General toggle for turning on dog leg damping. More... | |
int | doAffineSolve_ |
General toggle for turning on Affine solve with Hessian. More... | |
doublereal | m_conditionNumber |
Condition number of the matrix. More... | |
doublereal | CurrentTrustFactor_ |
Factor indicating how much trust region has been changed this iteration - output variable. More... | |
doublereal | NextTrustFactor_ |
Factor indicating how much trust region has been changed next iteration - output variable. More... | |
bool | ResidWtsReevaluated_ |
Boolean indicating that the residual weights have been reevaluated this iteration - output variable. More... | |
doublereal | ResidDecreaseSDExp_ |
Expected DResid_dS for the steepest descent path - output variable. More... | |
doublereal | ResidDecreaseSD_ |
Actual DResid_dS for the steepest descent path - output variable. More... | |
doublereal | ResidDecreaseNewtExp_ |
Expected DResid_dS for the Newton path - output variable. More... | |
doublereal | ResidDecreaseNewt_ |
Actual DResid_dS for the Newton path - output variable. More... | |
Class that calculates the solution to a nonlinear system.
This is a small nonlinear solver that can solve highly nonlinear problems that must use a dense matrix to relax the system.
Newton's method is used.
Damping is used extensively when relaxing the system.
The basic idea is that we predict a direction that is parameterized by an overall coordinate value, beta, from zero to one, This may or may not be the same as the value, damp, depending upon whether the direction is straight.
TIME STEP TYPE
The code solves a nonlinear problem. Frequently the nonlinear problem is created from time-dependent residual. Whenever you change the solution vector, you are also changing the derivative of the solution vector. Therefore, the code has the option of altering ydot, a vector of time derivatives of the solution in tandem with the solution vector and then feeding a residual and Jacobian routine with the time derivatives as well as the solution. The code has support for a backwards euler method and a second order Adams-Bashforth or Trapezoidal Rule.
In order to use these methods, the solver must be initialized with delta_t and m_y_nm1[i] to specify the conditions at the previous time step. For second order methods, the time derivative at t_nm1 must also be supplied, m_ydot_nm1[i]. Then the solution type NSOLN_TYPE_TIME_DEPENDENT may be used to solve the problem.
For steady state problem whose residual doesn't have a solution time derivative in it, you should use the NSOLN_TYPE_STEADY_STATE problem type.
We have a NSOLN_TYPE_PSEUDO_TIME_DEPENDENT defined. However, this is not implemented yet. This would be a pseudo time dependent calculation, where an optional time derivative could be added in order to help equilibrate a nonlinear steady state system. The time transient is not important in and of itself. Many physical systems have a time dependence to them that provides a natural way to relax the nonlinear system.
MATRIX SCALING
Definition at line 126 of file NonlinearSolver.h.
NonlinearSolver | ( | ResidJacEval * | func | ) |
Default constructor.
func | Residual and Jacobian evaluator function object |
Definition at line 58 of file NonlinearSolver.cpp.
References NonlinearSolver::atolBase_, NonlinearSolver::atolk_, NonlinearSolver::deltaX_CP_, NonlinearSolver::deltaX_Newton_, NonlinearSolver::deltaX_trust_, NonlinearSolver::Jd_, NonlinearSolver::m_colScales, NonlinearSolver::m_deltaStepMaximum, NonlinearSolver::m_deltaStepMinimum, NonlinearSolver::m_ewt, NonlinearSolver::m_func, NonlinearSolver::m_resid, NonlinearSolver::m_residWts, NonlinearSolver::m_rowScales, NonlinearSolver::m_rowWtScales, NonlinearSolver::m_step_1, NonlinearSolver::m_wksp, NonlinearSolver::m_wksp_2, NonlinearSolver::m_y_high_bounds, NonlinearSolver::m_y_low_bounds, NonlinearSolver::m_y_n_curr, NonlinearSolver::m_y_n_trial, NonlinearSolver::m_y_nm1, NonlinearSolver::m_ydot_n_curr, NonlinearSolver::m_ydot_nm1, NonlinearSolver::m_ydot_trial, NonlinearSolver::neq_, ResidJacEval::nEquations(), NonlinearSolver::rtol_, and Cantera::warn_deprecated().
NonlinearSolver | ( | const NonlinearSolver & | right | ) |
Copy Constructor.
Definition at line 159 of file NonlinearSolver.cpp.
References NonlinearSolver::operator=().
|
virtual |
Destructor.
Definition at line 224 of file NonlinearSolver.cpp.
References NonlinearSolver::HessianPtr_, and NonlinearSolver::jacCopyPtr_.
NonlinearSolver & operator= | ( | const NonlinearSolver & | right | ) |
Assignment operator.
Definition at line 230 of file NonlinearSolver.cpp.
References NonlinearSolver::atolBase_, NonlinearSolver::atolk_, NonlinearSolver::checkUserResidualTols_, NonlinearSolver::CurrentTrustFactor_, NonlinearSolver::delta_t_n, NonlinearSolver::deltaX_CP_, NonlinearSolver::deltaX_Newton_, NonlinearSolver::deltaX_trust_, NonlinearSolver::dist_R0_, NonlinearSolver::dist_R1_, NonlinearSolver::dist_R2_, NonlinearSolver::dist_Total_, NonlinearSolver::doAffineSolve_, NonlinearSolver::doDogLeg_, NonlinearSolver::dogLegAlpha_, NonlinearSolver::dogLegID_, ResidJacEval::duplMyselfAsResidJacEval(), NonlinearSolver::HessianPtr_, NonlinearSolver::jacCopyPtr_, NonlinearSolver::Jd_, NonlinearSolver::JdJd_norm_, NonlinearSolver::lambdaStar_, NonlinearSolver::m_colScales, NonlinearSolver::m_colScaling, NonlinearSolver::m_dampBound, NonlinearSolver::m_dampRes, NonlinearSolver::m_deltaStepMinimum, NonlinearSolver::m_ewt, NonlinearSolver::m_func, NonlinearSolver::m_jacFormMethod, NonlinearSolver::m_manualDeltaStepSet, NonlinearSolver::m_matrixConditioning, NonlinearSolver::m_min_newt_its, NonlinearSolver::m_nfe, NonlinearSolver::m_nJacEval, NonlinearSolver::m_normDeltaSoln_CP, NonlinearSolver::m_normDeltaSoln_Newton, NonlinearSolver::m_normResid_0, NonlinearSolver::m_normResid_1, NonlinearSolver::m_normResid_Bound, NonlinearSolver::m_normResidTrial, NonlinearSolver::m_numTotalLinearSolves, NonlinearSolver::m_numTotalNewtIts, NonlinearSolver::m_order, NonlinearSolver::m_print_flag, NonlinearSolver::m_resid, NonlinearSolver::m_resid_scaled, NonlinearSolver::m_residWts, NonlinearSolver::m_rowScales, NonlinearSolver::m_rowScaling, NonlinearSolver::m_rowWtScales, NonlinearSolver::m_ScaleSolnNormToResNorm, NonlinearSolver::m_step_1, NonlinearSolver::m_wksp, NonlinearSolver::m_wksp_2, NonlinearSolver::m_y_high_bounds, NonlinearSolver::m_y_low_bounds, NonlinearSolver::m_y_n_curr, NonlinearSolver::m_y_n_trial, NonlinearSolver::m_y_nm1, NonlinearSolver::m_ydot_n_curr, NonlinearSolver::m_ydot_nm1, NonlinearSolver::m_ydot_trial, NonlinearSolver::maxNewtIts_, NonlinearSolver::neq_, NonlinearSolver::NextTrustFactor_, NonlinearSolver::norm_deltaX_trust_, NonlinearSolver::normTrust_CP_, NonlinearSolver::normTrust_Newton_, NonlinearSolver::Nuu_, NonlinearSolver::ResidDecreaseNewt_, NonlinearSolver::ResidDecreaseNewtExp_, NonlinearSolver::ResidDecreaseSD_, NonlinearSolver::ResidDecreaseSDExp_, NonlinearSolver::residNorm2Cauchy_, NonlinearSolver::ResidWtsReevaluated_, NonlinearSolver::RJd_norm_, NonlinearSolver::rtol_, NonlinearSolver::solnType_, NonlinearSolver::time_n, NonlinearSolver::trustDelta_, NonlinearSolver::trustRegionInitializationFactor_, NonlinearSolver::trustRegionInitializationMethod_, NonlinearSolver::userResidAtol_, and NonlinearSolver::userResidRtol_.
Referenced by NonlinearSolver::NonlinearSolver().
void createSolnWeights | ( | const doublereal *const | y | ) |
Create solution weights for convergence criteria.
We create soln weights from the following formula
wt[i] = rtol * abs(y[i]) + atol[i]
The program always assumes that atol is specific to the solution component
y | vector of the current solution values |
Definition at line 331 of file NonlinearSolver.cpp.
References NonlinearSolver::atolk_, NonlinearSolver::m_ewt, NonlinearSolver::neq_, and NonlinearSolver::rtol_.
Referenced by NonlinearSolver::solve_nonlinear_problem().
doublereal solnErrorNorm | ( | const doublereal *const | delta_y, |
const char * | title = 0 , |
||
int | printLargest = 0 , |
||
const doublereal | dampFactor = 1.0 |
||
) | const |
L2 norm of the delta of the solution vector.
calculate the norm of the solution vector. This will involve the column scaling of the matrix
The third argument has a default of false. However, if true, then a table of the largest values is printed out to standard output.
delta_y | Vector to take the norm of |
title | Optional title to be printed out |
printLargest | int indicating how many specific lines should be printed out |
dampFactor | Current value of the damping factor. Defaults to 1. only used for printout out a table. |
Definition at line 366 of file NonlinearSolver.cpp.
References Cantera::error(), Cantera::int2str(), NonlinearSolver::m_ewt, NonlinearSolver::m_print_flag, NonlinearSolver::m_y_n_curr, NonlinearSolver::neq_, and Cantera::npos.
Referenced by NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::descentComparison(), NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::readjustTrustVector(), NonlinearSolver::residualComparisonLeg(), NonlinearSolver::setupDoubleDogleg(), NonlinearSolver::solve_nonlinear_problem(), and NonlinearSolver::trustRegionLength().
doublereal residErrorNorm | ( | const doublereal *const | resid, |
const char * | title = 0 , |
||
const int | printLargest = 0 , |
||
const doublereal *const | y = 0 |
||
) | const |
L2 norm of the residual of the equation system.
Calculate the norm of the residual vector. This may involve using the row sum scaling from the matrix problem.
The third argument has a default of false. However, if true, then a table of the largest values is printed out to standard output.
resid | Vector of the residuals |
title | Optional title to be printed out |
printLargest | Number of specific entries to be printed |
y | Current value of y - only used for printouts |
Definition at line 443 of file NonlinearSolver.cpp.
References Cantera::checkFinite(), Cantera::error(), Cantera::int2str(), NonlinearSolver::m_print_flag, NonlinearSolver::m_residWts, NonlinearSolver::neq_, and Cantera::npos.
Referenced by NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::descentComparison(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::solve_nonlinear_problem().
int doResidualCalc | ( | const doublereal | time_curr, |
const int | typeCalc, | ||
const doublereal *const | y_curr, | ||
const doublereal *const | ydot_curr, | ||
const ResidEval_Type_Enum | evalType = Base_ResidEval |
||
) | const |
Compute the current residual.
The current value of the residual is stored in the internal work array m_resid, which is defined as mutable
time_curr | Value of the time |
typeCalc | Type of the calculation |
y_curr | Current value of the solution vector |
ydot_curr | Current value of the time derivative of the solution vector |
evalType | Base evaluation type. Defaults to Base_ResidEval |
Definition at line 564 of file NonlinearSolver.cpp.
References DATA_PTR, NonlinearSolver::delta_t_n, ResidJacEval::evalResidNJ(), NonlinearSolver::m_func, NonlinearSolver::m_nfe, NonlinearSolver::m_resid, and NonlinearSolver::m_resid_scaled.
Referenced by NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::descentComparison(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::solve_nonlinear_problem().
int doNewtonSolve | ( | const doublereal | time_curr, |
const doublereal *const | y_curr, | ||
const doublereal *const | ydot_curr, | ||
doublereal *const | delta_y, | ||
GeneralMatrix & | jac | ||
) |
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. The Jacobian is not recomputed.
A factored Jacobian is reused, if available. If a factored Jacobian is not available, then the Jacobian is factored. Before factoring, the Jacobian is row and column-scaled. Column scaling is not recomputed. The row scales are recomputed here, after column scaling has been implemented.
time_curr | Current value of the time |
y_curr | Current value of the solution |
ydot_curr | Current value of the solution derivative. |
delta_y | return value of the raw change in y |
jac | Jacobian |
Definition at line 776 of file NonlinearSolver.cpp.
References DATA_PTR, NonlinearSolver::m_colScales, NonlinearSolver::m_colScaling, NonlinearSolver::m_numLocalLinearSolves, NonlinearSolver::m_numTotalLinearSolves, NonlinearSolver::m_resid, NonlinearSolver::m_resid_scaled, NonlinearSolver::m_rowScales, NonlinearSolver::m_rowScaling, NonlinearSolver::neq_, and GeneralMatrix::solve().
Referenced by NonlinearSolver::dampStep(), and NonlinearSolver::solve_nonlinear_problem().
int doAffineNewtonSolve | ( | const doublereal *const | y_curr, |
const doublereal *const | ydot_curr, | ||
doublereal *const | delta_y, | ||
GeneralMatrix & | jac | ||
) |
Compute the Newton step, either by direct Newton's or by solving a close problem that is represented by a Hessian.
This is algorith A.6.5.1 in Dennis / Schnabel
Compute the QR decomposition
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. The Jacobian is not recomputed.
A factored Jacobian is reused, if available. If a factored Jacobian is not available, then the Jacobian is factored. Before factoring, the Jacobian is row and column-scaled. Column scaling is not recomputed. The row scales are recomputed here, after column scaling has been implemented.
y_curr | Current value of the solution |
ydot_curr | Current value of the solution derivative. |
delta_y | return value of the raw change in y |
jac | Jacobian |
internal m_resid Stored residual is used as input
Definition at line 859 of file NonlinearSolver.cpp.
References GeneralMatrix::begin(), DATA_PTR, NonlinearSolver::doDogLeg_, GeneralMatrix::duplMyselfAsGeneralMatrix(), GeneralMatrix::factor(), GeneralMatrix::factorAlgorithm(), GeneralMatrix::factored(), GeneralMatrix::factorQR(), NonlinearSolver::HessianPtr_, Cantera::int2str(), NonlinearSolver::jacCopyPtr_, NonlinearSolver::m_colScales, NonlinearSolver::m_colScaling, NonlinearSolver::m_conditionNumber, NonlinearSolver::m_numLocalLinearSolves, NonlinearSolver::m_numTotalLinearSolves, NonlinearSolver::m_print_flag, NonlinearSolver::m_resid, NonlinearSolver::m_resid_scaled, NonlinearSolver::m_rowScales, NonlinearSolver::m_rowScaling, GeneralMatrix::matrixType_, NonlinearSolver::neq_, GeneralMatrix::oneNorm(), GeneralMatrix::rcond(), GeneralMatrix::rcondQR(), NonlinearSolver::s_alwaysAssumeNewtonGood, NonlinearSolver::s_doBothSolvesAndCompare, NonlinearSolver::solnErrorNorm(), GeneralMatrix::solve(), GeneralMatrix::useFactorAlgorithm(), and GeneralMatrix::zero().
Referenced by NonlinearSolver::solve_nonlinear_problem().
doublereal trustRegionLength | ( | ) | const |
Calculate the length of the current trust region in terms of the solution error norm.
We carry out a norm of deltaX_trust_ first. Then, we multiply that value by trustDelta_
Definition at line 1723 of file NonlinearSolver.cpp.
References DATA_PTR, NonlinearSolver::deltaX_trust_, NonlinearSolver::norm_deltaX_trust_, NonlinearSolver::solnErrorNorm(), and NonlinearSolver::trustDelta_.
Referenced by NonlinearSolver::dampDogLeg(), and NonlinearSolver::decideStep().
void setDefaultDeltaBoundsMagnitudes | ( | ) |
Set default deulta bounds amounts.
Delta bounds are set to 0.01 for all unknowns arbitrarily and capriciously Then, for each call to the nonlinear solver Then, they are increased to 1000 x atol then, they are increased to 0.1 fab(y[i])
Definition at line 1729 of file NonlinearSolver.cpp.
References NonlinearSolver::atolk_, NonlinearSolver::m_deltaStepMinimum, NonlinearSolver::m_y_n_curr, and NonlinearSolver::neq_.
Referenced by NonlinearSolver::solve_nonlinear_problem().
void adjustUpStepMinimums | ( | ) |
Adjust the step minimums.
Definition at line 1737 of file NonlinearSolver.cpp.
References NonlinearSolver::deltaX_trust_, NonlinearSolver::m_deltaStepMinimum, NonlinearSolver::neq_, and NonlinearSolver::trustDelta_.
Referenced by NonlinearSolver::decideStep().
void setDeltaBoundsMagnitudes | ( | const doublereal *const | deltaBoundsMagnitudes | ) |
Set the delta Bounds magnitudes by hand.
deltaBoundsMagnitudes | set the deltaBoundsMagnitude vector |
Definition at line 1748 of file NonlinearSolver.cpp.
References NonlinearSolver::m_deltaStepMinimum, NonlinearSolver::m_manualDeltaStepSet, and NonlinearSolver::neq_.
|
protected |
Readjust the trust region vectors.
The trust region is made up of the trust region vector calculation and the trustDelta_ value We periodically recalculate the trustVector_ values so that they renormalize to the correct length. We change the trustDelta_ values regularly
The trust region calculate is based on
|| delta_x dot 1/trustDeltaX_ || <= trustDelta_
Definition at line 1854 of file NonlinearSolver.cpp.
References DATA_PTR, NonlinearSolver::deltaX_trust_, NonlinearSolver::doDogLeg_, NonlinearSolver::m_deltaStepMaximum, NonlinearSolver::m_deltaStepMinimum, NonlinearSolver::m_ewt, NonlinearSolver::m_print_flag, NonlinearSolver::m_y_n_curr, NonlinearSolver::neq_, NonlinearSolver::norm_deltaX_trust_, NonlinearSolver::solnErrorNorm(), and NonlinearSolver::trustDelta_.
Referenced by NonlinearSolver::initializeTrustRegion(), and NonlinearSolver::solve_nonlinear_problem().
|
protected |
Fill a dogleg solution step vector.
Previously, we have filled up deltaX_Newton_[], deltaX_CP_[], and Nuu_, so that this routine is straightforward.
leg | Leg of the dog leg you are on (0, 1, or 2) |
alpha | Relative length along the dog length that you are on. |
deltaX | Vector to be filled up |
Definition at line 1957 of file NonlinearSolver.cpp.
References NonlinearSolver::deltaX_CP_, NonlinearSolver::deltaX_Newton_, NonlinearSolver::neq_, and NonlinearSolver::Nuu_.
Referenced by NonlinearSolver::dampDogLeg().
|
protected |
Calculate the trust distance of a step in the solution variables.
The trust distance is defined as the length of the step according to the norm wrt to the trust region. We calculate the trust distance by the following method.
trustDist = || delta_x dot 1/trustDeltaX_ ||
deltaX | Current value of deltaX |
Definition at line 1974 of file NonlinearSolver.cpp.
References NonlinearSolver::deltaX_trust_, NonlinearSolver::neq_, and NonlinearSolver::trustDelta_.
Referenced by NonlinearSolver::initializeTrustRegion(), NonlinearSolver::setupDoubleDogleg(), and NonlinearSolver::solve_nonlinear_problem().
doublereal boundStep | ( | const doublereal *const | y, |
const doublereal *const | step0 | ||
) |
Bound the step.
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 solution value of the old step |
step0 | Proposed step change in the solution |
Definition at line 2024 of file NonlinearSolver.cpp.
References NonlinearSolver::deltaBoundStep(), Cantera::int2str(), NonlinearSolver::m_print_flag, NonlinearSolver::m_y_high_bounds, NonlinearSolver::m_y_low_bounds, NonlinearSolver::neq_, and Cantera::npos.
Referenced by NonlinearSolver::dampDogLeg(), and NonlinearSolver::dampStep().
void setBoundsConstraints | ( | const doublereal *const | y_low_bounds, |
const doublereal *const | y_high_bounds | ||
) |
Set bounds constraints for all variables in the problem.
y_low_bounds | Vector of lower bounds |
y_high_bounds | Vector of high bounds |
Definition at line 341 of file NonlinearSolver.cpp.
References NonlinearSolver::m_y_high_bounds, NonlinearSolver::m_y_low_bounds, and NonlinearSolver::neq_.
std::vector< doublereal > & lowBoundsConstraintVector | ( | ) |
Return an editable vector of the low bounds constraints.
Definition at line 356 of file NonlinearSolver.cpp.
References NonlinearSolver::m_y_low_bounds.
std::vector< doublereal > & highBoundsConstraintVector | ( | ) |
Return an editable vector of the high bounds constraints.
Definition at line 361 of file NonlinearSolver.cpp.
References NonlinearSolver::m_y_high_bounds.
void calc_ydot | ( | const int | order, |
const doublereal *const | y_curr, | ||
doublereal *const | ydot_curr | ||
) | const |
Internal function to calculate the time derivative of the solution at the new step.
Previously, the user must have supplied information about the previous time step for this routine to work as intended.
order | of the BDF method |
y_curr | current value of the solution |
ydot_curr | Calculated value of the solution derivative that is consistent with y_curr |
Definition at line 3471 of file NonlinearSolver.cpp.
References NonlinearSolver::delta_t_n, NonlinearSolver::m_y_nm1, NonlinearSolver::m_ydot_nm1, and NonlinearSolver::neq_.
Referenced by NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::solve_nonlinear_problem().
int beuler_jac | ( | GeneralMatrix & | J, |
doublereal *const | f, | ||
doublereal | time_curr, | ||
doublereal | CJ, | ||
doublereal *const | y, | ||
doublereal *const | ydot, | ||
int | num_newt_its | ||
) |
Function called to evaluate the Jacobian matrix and the current residual vector 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. |
time_curr | Current time |
CJ | inverse of the value of deltaT |
y | value of the solution vector |
ydot | value of the time derivative of the solution vector |
num_newt_its | Number of Newton iterations |
Definition at line 3199 of file NonlinearSolver.cpp.
References ResidJacEval::calcDeltaSolnVariables(), GeneralMatrix::checkColumns(), Cantera::checkFinite(), GeneralMatrix::checkRows(), GeneralMatrix::clearFactorFlag(), DATA_PTR, NonlinearSolver::delta_t_n, ResidJacEval::evalJacobian(), ResidJacEval::evalResidNJ(), Cantera::int2str(), Cantera::JacBase_ResidEval, NonlinearSolver::jacCopyPtr_, Cantera::JacDelta_ResidEval, NonlinearSolver::m_ewt, NonlinearSolver::m_func, NonlinearSolver::m_jacFormMethod, NonlinearSolver::m_nfe, NonlinearSolver::m_nJacEval, NonlinearSolver::m_print_flag, NonlinearSolver::m_wksp, GeneralMatrix::matrixType_, NonlinearSolver::neq_, GeneralMatrix::nRowsAndStruct(), NSOLN_JAC_ANAL, NSOLN_TYPE_STEADY_STATE, GeneralMatrix::ptrColumn(), NonlinearSolver::s_print_NumJac, NonlinearSolver::solnType_, and Cantera::subtractRD().
Referenced by NonlinearSolver::solve_nonlinear_problem().
doublereal filterNewStep | ( | const doublereal | timeCurrent, |
const doublereal *const | ybase, | ||
doublereal *const | step0 | ||
) |
Apply a filtering process to the step.
timeCurrent | Current value of the time |
ybase | current value of the solution |
step0 | Proposed step change in the solution |
Definition at line 3499 of file NonlinearSolver.cpp.
References ResidJacEval::filterNewStep(), and NonlinearSolver::m_func.
Referenced by NonlinearSolver::solve_nonlinear_problem().
doublereal filterNewSolution | ( | const doublereal | timeCurrent, |
doublereal *const | y_current, | ||
doublereal *const | ydot_current | ||
) |
Apply a filter to the solution.
timeCurrent | Current value of the time |
y_current | current value of the solution |
ydot_current | Current value of the solution derivative. |
Definition at line 3505 of file NonlinearSolver.cpp.
References ResidJacEval::filterSolnPrediction(), and NonlinearSolver::m_func.
Referenced by NonlinearSolver::solve_nonlinear_problem().
double deltaBoundStep | ( | const doublereal *const | y, |
const doublereal *const | step0 | ||
) |
Return the factor by which the undamped Newton step 'step0' must be multiplied in order to keep the update within the bounds of an accurate Jacobian.
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 1.5 Maximum decrease in variable in any one Newton iteration: factor of 2
y | Initial value of the solution vector |
step0 | initial proposed step size |
Definition at line 1757 of file NonlinearSolver.cpp.
References Cantera::int2str(), NonlinearSolver::m_deltaStepMinimum, NonlinearSolver::m_print_flag, and NonlinearSolver::neq_.
Referenced by NonlinearSolver::boundStep().
int dampStep | ( | const doublereal | time_curr, |
const doublereal *const | y_n_curr, | ||
const doublereal *const | ydot_n_curr, | ||
doublereal *const | step_1, | ||
doublereal *const | y_n_1, | ||
doublereal *const | ydot_n_1, | ||
doublereal * | step_2, | ||
doublereal & | stepNorm_2, | ||
GeneralMatrix & | jac, | ||
bool | writetitle, | ||
int & | num_backtracks | ||
) |
Find a damping coefficient through a look-ahead mechanism.
On entry, step_1 must contain an undamped Newton step for the solution y_n_curr. This method attempts to find a damping coefficient such that all components stay in bounds, and the next undamped step would have a norm smaller than that of step_1. If successful, the new solution after taking the damped step is returned in y_n_1, and the undamped step at y_n_1 is returned in step_2.
time_curr | Current physical time |
y_n_curr | Base value of the solution before any steps are taken |
ydot_n_curr | Base value of the time derivative of the solution |
step_1 | Initial step suggested. |
y_n_1 | Value of y1, the suggested solution after damping |
ydot_n_1 | Value of the time derivative of the solution at y_n_1 |
step_2 | Value of the step change from y_n_1 to y_n_2 |
stepNorm_2 | norm of the step change in going from y_n_1 to y_n_2 |
jac | Jacobian |
writetitle | Write a title line |
num_backtracks | Number of backtracks taken |
Definition at line 2074 of file NonlinearSolver.cpp.
References Cantera::Base_LaggedSolutionComponents, NonlinearSolver::boundStep(), NonlinearSolver::calc_ydot(), Cantera::DampFactor, DATA_PTR, NonlinearSolver::doNewtonSolve(), NonlinearSolver::doResidualCalc(), NonlinearSolver::m_dampBound, NonlinearSolver::m_dampRes, NonlinearSolver::m_normResid_0, NonlinearSolver::m_normResid_1, NonlinearSolver::m_normResid_Bound, NonlinearSolver::m_normResidTrial, NonlinearSolver::m_order, NonlinearSolver::m_print_flag, NonlinearSolver::m_resid, NonlinearSolver::m_wksp, Cantera::NDAMP, NonlinearSolver::neq_, NSOLN_RETN_FAIL_DAMPSTEP, NSOLN_TYPE_STEADY_STATE, NonlinearSolver::print_solnDelta_norm_contrib(), NonlinearSolver::residErrorNorm(), NonlinearSolver::solnErrorNorm(), and NonlinearSolver::solnType_.
Referenced by NonlinearSolver::solve_nonlinear_problem().
int solve_nonlinear_problem | ( | int | SolnType, |
doublereal *const | y_comm, | ||
doublereal *const | ydot_comm, | ||
doublereal | CJ, | ||
doublereal | time_curr, | ||
GeneralMatrix & | jac, | ||
int & | num_newt_its, | ||
int & | num_linear_solves, | ||
int & | num_backtracks, | ||
int | loglevelInput | ||
) |
Find the solution to F(X) = 0 by damped Newton iteration.
On entry, x0 contains an initial estimate of the solution. On successful return, x1 contains the converged solution.
SolnType = TRANSIENT -> we will assume we are relaxing a transient equation system for now. Will make it more general later, if an application comes up.
SolnType | Solution type |
y_comm | Initial value of the solution. On return this is the converged value of the solution |
ydot_comm | Initial value of the solution derivative. On return this is the converged value of the solution derivative. |
CJ | Inverse of the value of deltaT |
time_curr | Current value of the time |
jac | Matrix that will be used to store the Jacobian |
num_newt_its | Number of Newton iterations taken |
num_linear_solves | Number of linear solves taken |
num_backtracks | Number of backtracking steps taken |
loglevelInput | Input log level determines the amount of printing. |
Definition at line 2591 of file NonlinearSolver.cpp.
References Cantera::Base_ShowSolution, NonlinearSolver::beuler_jac(), NonlinearSolver::calc_ydot(), NonlinearSolver::calcColumnScales(), NonlinearSolver::calcTrustDistance(), NonlinearSolver::convergenceCheck(), NonlinearSolver::createSolnWeights(), NonlinearSolver::CurrentTrustFactor_, NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), DATA_PTR, NonlinearSolver::deltaX_Newton_, NonlinearSolver::deltaX_trust_, NonlinearSolver::descentComparison(), NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doAffineSolve_, NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::doDogLeg_, NonlinearSolver::dogLegAlpha_, NonlinearSolver::dogLegID_, NonlinearSolver::doNewtonSolve(), NonlinearSolver::doResidualCalc(), GeneralMatrix::duplMyselfAsGeneralMatrix(), NonlinearSolver::filterNewSolution(), NonlinearSolver::filterNewStep(), NonlinearSolver::initializeTrustRegion(), NonlinearSolver::jacCopyPtr_, NonlinearSolver::m_conditionNumber, NonlinearSolver::m_dampBound, NonlinearSolver::m_dampRes, NonlinearSolver::m_manualDeltaStepSet, NonlinearSolver::m_min_newt_its, NonlinearSolver::m_normDeltaSoln_CP, NonlinearSolver::m_normDeltaSoln_Newton, NonlinearSolver::m_normResid_0, NonlinearSolver::m_normResid_1, NonlinearSolver::m_normResid_Bound, NonlinearSolver::m_normResid_full, NonlinearSolver::m_normResidTrial, NonlinearSolver::m_numLocalLinearSolves, NonlinearSolver::m_numTotalLinearSolves, NonlinearSolver::m_numTotalNewtIts, NonlinearSolver::m_order, NonlinearSolver::m_print_flag, NonlinearSolver::m_resid, NonlinearSolver::m_step_1, NonlinearSolver::m_wksp_2, NonlinearSolver::m_y_n_curr, NonlinearSolver::m_y_n_trial, NonlinearSolver::m_ydot_n_curr, NonlinearSolver::m_ydot_trial, NonlinearSolver::maxNewtIts_, NonlinearSolver::neq_, NonlinearSolver::NextTrustFactor_, NonlinearSolver::norm_deltaX_trust_, NSOLN_RETN_CONTINUE, NSOLN_RETN_JACOBIANFORMATIONERROR, NSOLN_RETN_MATRIXINVERSIONERROR, NSOLN_RETN_MAXIMUMITERATIONSEXCEEDED, NSOLN_RETN_RESIDUALFORMATIONERROR, NSOLN_RETN_SUCCESS, NSOLN_TYPE_STEADY_STATE, NonlinearSolver::readjustTrustVector(), NonlinearSolver::ResidDecreaseNewt_, NonlinearSolver::ResidDecreaseNewtExp_, NonlinearSolver::ResidDecreaseSD_, NonlinearSolver::ResidDecreaseSDExp_, NonlinearSolver::residErrorNorm(), NonlinearSolver::residualComparisonLeg(), NonlinearSolver::ResidWtsReevaluated_, NonlinearSolver::s_print_DogLeg, NonlinearSolver::s_TurnOffTiming, NonlinearSolver::scaleMatrix(), clockWC::secondsWC(), NonlinearSolver::setDefaultDeltaBoundsMagnitudes(), NonlinearSolver::setupDoubleDogleg(), NonlinearSolver::solnErrorNorm(), NonlinearSolver::solnType_, NonlinearSolver::time_n, NonlinearSolver::trustDelta_, NonlinearSolver::trustRegionInitializationFactor_, and NonlinearSolver::trustRegionInitializationMethod_.
|
virtual |
Set the values for the previous time step.
We set the values for the previous time step here. These are used in the nonlinear solve because they affect the calculation of ydot.
y_nm1 | Value of the solution vector at the previous time step |
ydot_nm1 | Value of the solution vector derivative at the previous time step |
Definition at line 3105 of file NonlinearSolver.cpp.
References NonlinearSolver::m_y_nm1, and NonlinearSolver::m_ydot_nm1.
|
private |
Set the column scaling vector at the current time.
Definition at line 548 of file NonlinearSolver.cpp.
References ResidJacEval::calcSolnScales(), DATA_PTR, NonlinearSolver::m_colScales, NonlinearSolver::m_colScaling, NonlinearSolver::m_ewt, NonlinearSolver::m_func, NonlinearSolver::m_y_n_curr, NonlinearSolver::m_y_nm1, NonlinearSolver::neq_, and NonlinearSolver::time_n.
Referenced by NonlinearSolver::solve_nonlinear_problem().
void setColumnScaling | ( | bool | useColScaling, |
const double *const | scaleFactors = 0 |
||
) |
Set the column scaling that are used for the inversion of the matrix.
There are three ways to do this.
The first method is to set the bool useColScaling to true, leaving the scaling factors unset. Then, the column scales will be set to the solution error weighting factors. This has the effect of ensuring that all delta variables will have the same order of magnitude at convergence end.
The second way is the explicitly set the column factors in the second argument of this function call.
The final way to input the scales is to override the ResidJacEval member function call,
calcSolnScales(double time_n, const double *m_y_n_curr, const double *m_y_nm1, double *m_colScales)
Overriding this function call will trump all other ways to specify the column scaling factors.
useColScaling | Turn this on if you want to use column scaling in the calculations |
scaleFactors | A vector of doubles that specifies the column factors. |
Definition at line 524 of file NonlinearSolver.cpp.
References NonlinearSolver::m_colScales, NonlinearSolver::m_colScaling, and NonlinearSolver::neq_.
void setRowScaling | ( | bool | useRowScaling | ) |
Set the rowscaling that are used for the inversion of the matrix.
Row scaling is set here. Right now the row scaling is set internally in the code.
useRowScaling | Turn row scaling on or off. |
Definition at line 543 of file NonlinearSolver.cpp.
References NonlinearSolver::m_rowScaling.
void scaleMatrix | ( | GeneralMatrix & | jac, |
doublereal *const | y_comm, | ||
doublereal *const | ydot_comm, | ||
doublereal | time_curr, | ||
int | num_newt_its | ||
) |
Scale the matrix.
jac | Jacobian |
y_comm | Current value of the solution vector |
ydot_comm | Current value of the time derivative of the solution vector |
time_curr | current value of the time |
num_newt_its | Current value of the number of newt its |
Definition at line 573 of file NonlinearSolver.cpp.
References GeneralMatrix::begin(), Cantera::checkFinite(), NonlinearSolver::computeResidWts(), GeneralMatrix::factored(), NonlinearSolver::m_colScales, NonlinearSolver::m_colScaling, NonlinearSolver::m_ewt, NonlinearSolver::m_rowScales, NonlinearSolver::m_rowScaling, NonlinearSolver::m_rowWtScales, GeneralMatrix::matrixType_, NonlinearSolver::neq_, GeneralMatrix::nRowsAndStruct(), and GeneralMatrix::ptrColumn().
Referenced by NonlinearSolver::solve_nonlinear_problem().
void print_solnDelta_norm_contrib | ( | const doublereal *const | step_1, |
const char *const | stepNorm_1, | ||
const doublereal *const | step_2, | ||
const char *const | stepNorm_2, | ||
const char *const | title, | ||
const doublereal *const | y_n_curr, | ||
const doublereal *const | y_n_1, | ||
doublereal | damp, | ||
size_t | num_entries | ||
) |
Print solution norm contribution.
Prints out the most important entries to the update to the solution vector for the current step
step_1 | Raw update vector for the current nonlinear step |
stepNorm_1 | Norm of the vector step_1 |
step_2 | Raw update vector for the next solution value based on the old matrix |
stepNorm_2 | Norm of the vector step_2 |
title | title of the printout |
y_n_curr | Old value of the solution |
y_n_1 | New value of the solution after damping corrections |
damp | Value of the damping factor |
num_entries | Number of entries to print out |
Definition at line 3112 of file NonlinearSolver.cpp.
References Cantera::error(), Cantera::int2str(), NonlinearSolver::m_ewt, NonlinearSolver::neq_, and Cantera::npos.
Referenced by NonlinearSolver::dampStep().
void computeResidWts | ( | ) |
Compute the Residual Weights.
The residual weights are defined here to be equal to the inverse of the row scaling factors used to row scale the matrix, after column scaling is used. They are multiplied by 10-3 because the column weights are also multiplied by that same quantity.
The basic idea is that a change in the solution vector on the order of the convergence tolerance multiplied by [RJC] which is of order one after row scaling should give you the relative weight of the row. Values of the residual for that row can then be normalized by the value of this weight. When the tolerance in delta x is achieved, the tolerance in the residual is also achieved.
Definition at line 3511 of file NonlinearSolver.cpp.
References NonlinearSolver::atolBase_, Cantera::checkFinite(), NonlinearSolver::checkUserResidualTols_, NonlinearSolver::m_residWts, NonlinearSolver::m_rowWtScales, NonlinearSolver::m_ScaleSolnNormToResNorm, NonlinearSolver::neq_, NonlinearSolver::ResidWtsReevaluated_, NonlinearSolver::rtol_, NonlinearSolver::userResidAtol_, and NonlinearSolver::userResidRtol_.
Referenced by NonlinearSolver::calcSolnToResNormVector(), and NonlinearSolver::scaleMatrix().
void getResidWts | ( | doublereal *const | residWts | ) | const |
Return the residual weights.
residWts | Vector of length neq_ |
Definition at line 3543 of file NonlinearSolver.cpp.
References NonlinearSolver::m_residWts, and NonlinearSolver::neq_.
int convergenceCheck | ( | int | dampCode, |
doublereal | s1 | ||
) |
Check to see if the nonlinear problem has converged.
dampCode | Code from the damping routine |
s1 | Value of the norm of the step change |
Definition at line 3550 of file NonlinearSolver.cpp.
References NonlinearSolver::m_dampBound, NonlinearSolver::m_dampRes, NonlinearSolver::m_normDeltaSoln_Newton, and NonlinearSolver::m_normResid_full.
Referenced by NonlinearSolver::solve_nonlinear_problem().
void setAtol | ( | const doublereal *const | atol | ) |
Set the absolute tolerances for the solution variables.
Set the absolute tolerances used in the calculation
atol | Vector of length neq_ that contains the tolerances to be used for the solution variables |
Definition at line 3601 of file NonlinearSolver.cpp.
References NonlinearSolver::atolk_, and NonlinearSolver::neq_.
void setRtol | ( | const doublereal | rtol | ) |
Set the relative tolerances for the solution variables.
Set the relative tolerances used in the calculation
rtol | single double |
Definition at line 3612 of file NonlinearSolver.cpp.
References NonlinearSolver::rtol_.
void setResidualTols | ( | double | residRtol, |
double * | residAtol, | ||
int | residNormHandling = 2 |
||
) |
Set the relative and absolute tolerances for the Residual norm comparisons, if used.
Residual norms are used to calculate convergence within the nonlinear solver, since these are the norms that are associated with convergence proofs, especially for ill-conditioned systems. Usually the residual weights for each row are calculated by the program such that they correlate with the convergence requirements on the solution variables input by the user using the routines setAtol() and setRtol(). The residual weights are essentially calculated from the value
residWeightNorm[i] = m_ScaleSolnNormToResNorm * sum_j ( fabs(A_i,j) ewt(j))
The factor, m_ScaleSolnNormToResNorm, is computed periodically to ensure that the solution norms and the residual norms are converging at the same time and thus accounts for some-illconditioning issues but not all.
With this routine the user can override or add to the residual weighting norm evaluation by specifying their own vector of residual absolute and relative tolerances.
The user specified tolerance for the residual is given by the following quantity
residWeightNorm[i] = residAtol[i] + residRtol * m_rowWtScales[i] / neq
residRtol | scalar residual relative tolerance |
residAtol | vector of residual absolute tolerances |
residNormHandling | Parameter that sets the default handling of the residual norms 0 The residual weighting vector is calculated to make sure that the solution norms are roughly 1 when the residual norm is roughly 1. This is the default if this routine is not called. 1 Use the user residual norm specified by the parameters in this routine 2 Use the minimum value of the residual weights calculated by method 1 and 2. This is the default if this routine is called and this parameter isn't specified. |
Definition at line 3621 of file NonlinearSolver.cpp.
References NonlinearSolver::checkUserResidualTols_, NonlinearSolver::neq_, NonlinearSolver::userResidAtol_, and NonlinearSolver::userResidRtol_.
void setMaxNewtIts | ( | const int | maxNewtIts | ) |
Set the value of the maximum # of Newton iterations.
maxNewtIts | Maximum number of Newton iterations |
void calcSolnToResNormVector | ( | ) |
Calculate the scaling factor for translating residual norms into solution norms.
This routine calls computeResidWts() a couple of times in the calculation of m_ScaleSolnNormToResNorm. A more sophisticated routine may do more with signs to get a better value. Perhaps, a series of calculations with different signs attached may be in order. Then, m_ScaleSolnNormToResNorm would be calculated as the minimum of a series of calculations.
Definition at line 719 of file NonlinearSolver.cpp.
References NonlinearSolver::atolBase_, NonlinearSolver::checkUserResidualTols_, NonlinearSolver::computeResidWts(), Cantera::error(), GeneralMatrix::factored(), NonlinearSolver::jacCopyPtr_, NonlinearSolver::m_ewt, NonlinearSolver::m_residWts, NonlinearSolver::m_rowWtScales, NonlinearSolver::m_ScaleSolnNormToResNorm, NonlinearSolver::m_wksp, NonlinearSolver::neq_, NonlinearSolver::userResidAtol_, and NonlinearSolver::userResidRtol_.
doublereal doCauchyPointSolve | ( | GeneralMatrix & | jac | ) |
Calculate the steepest descent direction and the Cauchy Point where the quadratic formulation of the nonlinear problem expects a minimum along the descent direction.
jac | Jacobian matrix: must be unfactored. |
Definition at line 1204 of file NonlinearSolver.cpp.
References Cantera::checkFinite(), DATA_PTR, NonlinearSolver::deltaX_CP_, NonlinearSolver::doDogLeg_, NonlinearSolver::Jd_, NonlinearSolver::JdJd_norm_, NonlinearSolver::lambdaStar_, NonlinearSolver::m_colScales, NonlinearSolver::m_colScaling, NonlinearSolver::m_ewt, NonlinearSolver::m_normResid_0, NonlinearSolver::m_print_flag, NonlinearSolver::m_resid, NonlinearSolver::m_residWts, NonlinearSolver::m_rowScales, NonlinearSolver::m_rowScaling, NonlinearSolver::neq_, NonlinearSolver::residNorm2Cauchy_, NonlinearSolver::RJd_norm_, NonlinearSolver::s_print_DogLeg, and NonlinearSolver::solnErrorNorm().
Referenced by NonlinearSolver::solve_nonlinear_problem().
void descentComparison | ( | doublereal | time_curr, |
doublereal * | ydot0, | ||
doublereal * | ydot1, | ||
int & | numTrials | ||
) |
This is a utility routine that can be used to print out the rates of the initial residual decline.
The residual**2 decline for various directions is printed out. The rate of decline of the square of the residuals multiplied by the number of equations along each direction is printed out This quantity can be directly related to the theory, and may be calculated from derivatives at the original point.
( (r)**2 * neq - (r0)**2 * neq ) / distance
What's printed out:
The theoretical linearized residual decline The actual residual decline in the steepest descent direction determined by numerical differencing The actual residual decline in the Newton direction determined by numerical differencing
This routine doesn't need to be called for the solution of the nonlinear problem.
time_curr | Current time |
ydot0 | INPUT Current value of the derivative of the solution vector |
ydot1 | INPUT Time derivatives of solution at the conditions which are evaluated for success |
numTrials | OUTPUT Counter for the number of residual evaluations |
Definition at line 1331 of file NonlinearSolver.cpp.
References Cantera::Base_LaggedSolutionComponents, DATA_PTR, NonlinearSolver::deltaX_CP_, NonlinearSolver::deltaX_Newton_, NonlinearSolver::doDogLeg_, NonlinearSolver::doResidualCalc(), NonlinearSolver::lambdaStar_, NonlinearSolver::m_normResid_0, NonlinearSolver::m_print_flag, NonlinearSolver::m_resid, NonlinearSolver::m_wksp, NonlinearSolver::m_y_n_curr, NonlinearSolver::neq_, NSOLN_TYPE_STEADY_STATE, NonlinearSolver::ResidDecreaseNewt_, NonlinearSolver::ResidDecreaseNewtExp_, NonlinearSolver::ResidDecreaseSD_, NonlinearSolver::ResidDecreaseSDExp_, NonlinearSolver::residErrorNorm(), NonlinearSolver::RJd_norm_, NonlinearSolver::s_print_DogLeg, NonlinearSolver::solnErrorNorm(), and NonlinearSolver::solnType_.
Referenced by NonlinearSolver::solve_nonlinear_problem().
void setupDoubleDogleg | ( | ) |
Setup the parameters for the double dog leg.
The calls to the doCauchySolve() and doNewtonSolve() routines are done at the main level. This routine comes after those calls. We calculate the point Nuu_ here, the distances of the dog-legs, and the norms of the CP and Newton points in terms of the trust vectors.
Definition at line 1454 of file NonlinearSolver.cpp.
References NonlinearSolver::calcTrustDistance(), DATA_PTR, NonlinearSolver::deltaX_CP_, NonlinearSolver::deltaX_Newton_, NonlinearSolver::dist_R0_, NonlinearSolver::dist_R1_, NonlinearSolver::dist_R2_, NonlinearSolver::dist_Total_, NonlinearSolver::expectedResidLeg(), NonlinearSolver::m_normDeltaSoln_CP, NonlinearSolver::m_normDeltaSoln_Newton, NonlinearSolver::m_normResid_0, NonlinearSolver::m_wksp, NonlinearSolver::neq_, NonlinearSolver::normTrust_CP_, NonlinearSolver::normTrust_Newton_, NonlinearSolver::Nuu_, and NonlinearSolver::solnErrorNorm().
Referenced by NonlinearSolver::solve_nonlinear_problem().
int lambdaToLeg | ( | const doublereal | lambda, |
doublereal & | alpha | ||
) | const |
Change the global lambda coordinate into the (leg,alpha) coordinate for the double dogleg.
lambda | Global value of the distance along the double dogleg |
alpha | relative value along the particular leg |
Definition at line 1500 of file NonlinearSolver.cpp.
References NonlinearSolver::dist_R0_, NonlinearSolver::dist_R1_, NonlinearSolver::dist_R2_, and NonlinearSolver::dist_Total_.
int calcTrustIntersection | ( | doublereal | trustVal, |
doublereal & | lambda, | ||
doublereal & | alpha | ||
) | const |
Given a trust distance, this routine calculates the intersection of the this distance with the double dogleg curve.
trustVal | (INPUT) Value of the trust distance |
lambda | (OUTPUT) Returns the internal coordinate of the double dogleg |
alpha | (OUTPUT) Returns the relative distance along the appropriate leg |
Definition at line 1985 of file NonlinearSolver.cpp.
References NonlinearSolver::deltaX_CP_, NonlinearSolver::deltaX_Newton_, NonlinearSolver::deltaX_trust_, NonlinearSolver::dist_R0_, NonlinearSolver::dist_R1_, NonlinearSolver::dist_R2_, NonlinearSolver::dist_Total_, NonlinearSolver::neq_, NonlinearSolver::normTrust_CP_, NonlinearSolver::normTrust_Newton_, and NonlinearSolver::Nuu_.
Referenced by NonlinearSolver::dampDogLeg().
void initializeTrustRegion | ( | ) |
Initialize the size of the trust vector.
The algorithm we use is to set it equal to the length of the Distance to the Cauchy point.
Definition at line 1914 of file NonlinearSolver.cpp.
References NonlinearSolver::calcTrustDistance(), NonlinearSolver::deltaX_CP_, NonlinearSolver::deltaX_Newton_, NonlinearSolver::deltaX_trust_, NonlinearSolver::doDogLeg_, NonlinearSolver::m_ewt, NonlinearSolver::m_normDeltaSoln_CP, NonlinearSolver::m_normDeltaSoln_Newton, NonlinearSolver::m_print_flag, NonlinearSolver::neq_, NonlinearSolver::readjustTrustVector(), NonlinearSolver::trustDelta_, NonlinearSolver::trustRegionInitializationFactor_, and NonlinearSolver::trustRegionInitializationMethod_.
Referenced by NonlinearSolver::solve_nonlinear_problem().
void setTrustRegionInitializationMethod | ( | int | method, |
doublereal | factor | ||
) |
Set Trust region initialization strategy.
The default is use method 2 with a factor of 1. Then, on subsequent invocations of solve_nonlinear_problem() the strategy flips to method 0.
method | Method to set the strategy 0 No strategy - Use the previous strategy 1 Factor of the solution error weights 2 Factor of the first Cauchy Point distance 3 Factor of the first Newton step distance |
factor | Factor to use in combination with the method |
int dampDogLeg | ( | const doublereal | time_curr, |
const doublereal * | y_n_curr, | ||
const doublereal * | ydot_n_curr, | ||
std::vector< doublereal > & | step_1, | ||
doublereal *const | y_n_1, | ||
doublereal *const | ydot_n_1, | ||
doublereal & | stepNorm_1, | ||
doublereal & | stepNorm_2, | ||
GeneralMatrix & | jac, | ||
int & | num_backtracks | ||
) |
Damp using the dog leg approach.
time_curr | INPUT Current value of the time |
y_n_curr | INPUT Current value of the solution vector |
ydot_n_curr | INPUT Current value of the derivative of the solution vector |
step_1 | INPUT First trial step for the first iteration |
y_n_1 | INPUT First trial value of the solution vector |
ydot_n_1 | INPUT First trial value of the derivative of the solution vector |
stepNorm_1 | OUTPUT Norm of the vector step_1 |
stepNorm_2 | OUTPUT Estimated norm of the vector step_2 |
jac | INPUT Jacobian |
num_backtracks | OUTPUT number of backtracks taken in the current damping step |
Definition at line 2282 of file NonlinearSolver.cpp.
References NonlinearSolver::boundStep(), NonlinearSolver::calc_ydot(), NonlinearSolver::calcTrustIntersection(), DATA_PTR, NonlinearSolver::decideStep(), NonlinearSolver::dogLegAlpha_, NonlinearSolver::dogLegID_, NonlinearSolver::fillDogLegStep(), NonlinearSolver::m_dampBound, NonlinearSolver::m_dampRes, NonlinearSolver::m_normResid_0, NonlinearSolver::m_normResid_1, NonlinearSolver::m_normResid_Bound, NonlinearSolver::m_normResidTrial, NonlinearSolver::m_order, NonlinearSolver::m_print_flag, NonlinearSolver::m_wksp, Cantera::NDAMP, NonlinearSolver::neq_, NonlinearSolver::normTrust_Newton_, NSOLN_RETN_FAIL_DAMPSTEP, NSOLN_TYPE_STEADY_STATE, NonlinearSolver::solnErrorNorm(), NonlinearSolver::solnType_, NonlinearSolver::trustDelta_, and NonlinearSolver::trustRegionLength().
Referenced by NonlinearSolver::solve_nonlinear_problem().
int decideStep | ( | const doublereal | time_curr, |
int | leg, | ||
doublereal | alpha, | ||
const doublereal *const | y_n_curr, | ||
const doublereal *const | ydot_n_curr, | ||
const std::vector< doublereal > & | step_1, | ||
const doublereal *const | y_n_1, | ||
const doublereal *const | ydot_n_1, | ||
doublereal | trustDeltaOld | ||
) |
Decide whether the current step is acceptable and adjust the trust region size.
This is an extension of algorithm 6.4.5 of Dennis and Schnabel.
Here we decide whether to accept the current step At the end of the calculation a new estimate of the trust region is calculated
time_curr | INPUT Current value of the time |
leg | INPUT Leg of the dogleg that we are on |
alpha | INPUT Distance down that leg that we are on |
y_n_curr | INPUT Current value of the solution vector |
ydot_n_curr | INPUT Current value of the derivative of the solution vector |
step_1 | INPUT Trial step |
y_n_1 | OUTPUT Solution values at the conditions which are evaluated for success |
ydot_n_1 | OUTPUT Time derivatives of solution at the conditions which are evaluated for success |
trustDeltaOld | INPUT Value of the trust length at the old conditions |
Definition at line 2431 of file NonlinearSolver.cpp.
References NonlinearSolver::adjustUpStepMinimums(), Cantera::Base_LaggedSolutionComponents, NonlinearSolver::CurrentTrustFactor_, DATA_PTR, NonlinearSolver::deltaX_CP_, NonlinearSolver::dogLegAlpha_, NonlinearSolver::dogLegID_, NonlinearSolver::doResidualCalc(), NonlinearSolver::expectedResidLeg(), NonlinearSolver::lambdaStar_, NonlinearSolver::m_dampBound, NonlinearSolver::m_normDeltaSoln_Newton, NonlinearSolver::m_normResid_0, NonlinearSolver::m_normResid_1, NonlinearSolver::m_normResidTrial, NonlinearSolver::m_print_flag, NonlinearSolver::m_resid, NonlinearSolver::neq_, NonlinearSolver::NextTrustFactor_, NSOLN_TYPE_STEADY_STATE, NonlinearSolver::residErrorNorm(), NonlinearSolver::RJd_norm_, NonlinearSolver::rtol_, NonlinearSolver::solnErrorNorm(), NonlinearSolver::solnType_, NonlinearSolver::trustDelta_, and NonlinearSolver::trustRegionLength().
Referenced by NonlinearSolver::dampDogLeg().
doublereal expectedResidLeg | ( | int | leg, |
doublereal | alpha | ||
) | const |
Calculated the expected residual along the double dogleg curve.
leg | 0, 1, or 2 representing the curves of the dogleg |
alpha | Relative distance along the particular curve. |
Definition at line 1514 of file NonlinearSolver.cpp.
References NonlinearSolver::lambdaStar_, NonlinearSolver::m_normResid_0, NonlinearSolver::neq_, NonlinearSolver::Nuu_, and NonlinearSolver::RJd_norm_.
Referenced by NonlinearSolver::decideStep(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::setupDoubleDogleg().
void residualComparisonLeg | ( | const doublereal | time_curr, |
const doublereal *const | ydot0, | ||
int & | legBest, | ||
doublereal & | alphaBest | ||
) | const |
Here we print out the residual at various points along the double dogleg, comparing against the quadratic model in a table format.
time_curr | INPUT current time |
ydot0 | INPUT Current value of the derivative of the solution vector for non-time dependent determinations |
legBest | OUTPUT leg of the dogleg that gives the lowest residual |
alphaBest | OUTPUT distance along dogleg for best result. |
Definition at line 1567 of file NonlinearSolver.cpp.
References Cantera::Base_LaggedSolutionComponents, NonlinearSolver::calc_ydot(), DATA_PTR, NonlinearSolver::deltaX_CP_, NonlinearSolver::deltaX_Newton_, NonlinearSolver::doDogLeg_, NonlinearSolver::doResidualCalc(), NonlinearSolver::expectedResidLeg(), NonlinearSolver::m_order, NonlinearSolver::m_print_flag, NonlinearSolver::m_resid, NonlinearSolver::m_wksp, NonlinearSolver::m_wksp_2, NonlinearSolver::m_y_n_curr, NonlinearSolver::neq_, NSOLN_TYPE_STEADY_STATE, NonlinearSolver::Nuu_, NonlinearSolver::residErrorNorm(), NonlinearSolver::s_print_DogLeg, NonlinearSolver::solnErrorNorm(), and NonlinearSolver::solnType_.
Referenced by NonlinearSolver::solve_nonlinear_problem().
void setPrintLvl | ( | int | printLvl | ) |
Set the print level from the nonlinear solver.
0 -> absolutely nothing is printed for a single time step. 1 -> One line summary per solve_nonlinear call 2 -> short description, points of interest: Table of nonlinear solve - one line per iteration 3 -> Table is included -> More printing per nonlinear iteration (default) that occurs during the table 4 -> Summaries of the nonlinear solve iteration as they are occurring -> table no longer printed 5 -> Algorithm information on the nonlinear iterates are printed out 6 -> Additional info on the nonlinear iterates are printed out 7 -> Additional info on the linear solve is printed out. 8 -> Info on a per iterate of the linear solve is printed out.
printLvl | integer value |
Definition at line 3642 of file NonlinearSolver.cpp.
References NonlinearSolver::m_print_flag.
void setSolverScheme | ( | int | doDogLeg, |
int | doAffineSolve | ||
) |
Parameter to turn on solution solver schemes.
doDogLeg | Parameter to turn on the double dog leg scheme Default is to always use a damping scheme in the Newton Direction. When this is nonzero, a model trust region approach is used using a double dog leg with the steepest descent direction used for small step sizes. |
doAffineSolve | Parameter to turn on or off the solution of the system using a Hessian if the matrix has a bad condition number. |
Definition at line 350 of file NonlinearSolver.cpp.
References NonlinearSolver::doAffineSolve_, and NonlinearSolver::doDogLeg_.
|
private |
Pointer to the residual and Jacobian evaluator for the function.
See ResidJacEval.h for an evaluator.
Definition at line 899 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac(), NonlinearSolver::calcColumnScales(), NonlinearSolver::doResidualCalc(), NonlinearSolver::filterNewSolution(), NonlinearSolver::filterNewStep(), NonlinearSolver::NonlinearSolver(), and NonlinearSolver::operator=().
|
private |
Solution type.
Definition at line 902 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac(), NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::descentComparison(), NonlinearSolver::operator=(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Local copy of the number of equations.
Definition at line 905 of file NonlinearSolver.h.
Referenced by NonlinearSolver::adjustUpStepMinimums(), NonlinearSolver::beuler_jac(), NonlinearSolver::boundStep(), NonlinearSolver::calc_ydot(), NonlinearSolver::calcColumnScales(), NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::calcTrustDistance(), NonlinearSolver::calcTrustIntersection(), NonlinearSolver::computeResidWts(), NonlinearSolver::createSolnWeights(), NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::deltaBoundStep(), NonlinearSolver::descentComparison(), NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::doNewtonSolve(), NonlinearSolver::expectedResidLeg(), NonlinearSolver::fillDogLegStep(), NonlinearSolver::getResidWts(), NonlinearSolver::initializeTrustRegion(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::print_solnDelta_norm_contrib(), NonlinearSolver::readjustTrustVector(), NonlinearSolver::residErrorNorm(), NonlinearSolver::residualComparisonLeg(), NonlinearSolver::scaleMatrix(), NonlinearSolver::setAtol(), NonlinearSolver::setBoundsConstraints(), NonlinearSolver::setColumnScaling(), NonlinearSolver::setDefaultDeltaBoundsMagnitudes(), NonlinearSolver::setDeltaBoundsMagnitudes(), NonlinearSolver::setResidualTols(), NonlinearSolver::setupDoubleDogleg(), NonlinearSolver::solnErrorNorm(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Soln error weights.
Definition at line 908 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac(), NonlinearSolver::calcColumnScales(), NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::createSolnWeights(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::initializeTrustRegion(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::print_solnDelta_norm_contrib(), NonlinearSolver::readjustTrustVector(), NonlinearSolver::scaleMatrix(), and NonlinearSolver::solnErrorNorm().
|
private |
Boolean indicating whether a manual delta bounds has been input.
Definition at line 911 of file NonlinearSolver.h.
Referenced by NonlinearSolver::operator=(), NonlinearSolver::setDeltaBoundsMagnitudes(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Soln Delta bounds magnitudes.
Definition at line 914 of file NonlinearSolver.h.
Referenced by NonlinearSolver::adjustUpStepMinimums(), NonlinearSolver::deltaBoundStep(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::readjustTrustVector(), NonlinearSolver::setDefaultDeltaBoundsMagnitudes(), and NonlinearSolver::setDeltaBoundsMagnitudes().
|
private |
Value of the delta step magnitudes.
Definition at line 917 of file NonlinearSolver.h.
Referenced by NonlinearSolver::NonlinearSolver(), and NonlinearSolver::readjustTrustVector().
|
private |
Vector containing the current solution vector within the nonlinear solver.
Definition at line 920 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcColumnScales(), NonlinearSolver::descentComparison(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::readjustTrustVector(), NonlinearSolver::residualComparisonLeg(), NonlinearSolver::setDefaultDeltaBoundsMagnitudes(), NonlinearSolver::solnErrorNorm(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Vector containing the time derivative of the current solution vector within the nonlinear solver (where applicable)
Definition at line 924 of file NonlinearSolver.h.
Referenced by NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Vector containing the solution at the previous time step.
Definition at line 927 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calc_ydot(), NonlinearSolver::calcColumnScales(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::setPreviousTimeStep().
|
private |
Vector containing the solution derivative at the previous time step.
Definition at line 930 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calc_ydot(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::setPreviousTimeStep().
|
private |
Vector containing the solution at the new point that is to be considered.
Definition at line 933 of file NonlinearSolver.h.
Referenced by NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Value of the solution time derivative at the new point that is to be considered.
Definition at line 936 of file NonlinearSolver.h.
Referenced by NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Value of the step to be taken in the solution.
Definition at line 939 of file NonlinearSolver.h.
Referenced by NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Vector of column scaling factors.
Definition at line 942 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcColumnScales(), NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::doNewtonSolve(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::scaleMatrix(), and NonlinearSolver::setColumnScaling().
|
private |
Weights for normalizing the values of the residuals.
These are computed if row scaling, m_rowScaling, is turned on. They are calculated currently as the sum of the absolute values of the rows of the Jacobian.
Definition at line 949 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::doNewtonSolve(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::scaleMatrix().
|
private |
Weights for normalizing the values of the residuals.
They are calculated as the sum of the absolute values of the Jacobian multiplied by the solution weight function. This is carried out in scaleMatrix().
Definition at line 957 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::computeResidWts(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::scaleMatrix().
|
mutableprivate |
Value of the residual for the nonlinear problem.
Definition at line 960 of file NonlinearSolver.h.
Referenced by NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::descentComparison(), NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::doNewtonSolve(), NonlinearSolver::doResidualCalc(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::solve_nonlinear_problem().
|
mutableprivate |
Workspace of length neq_.
Definition at line 963 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac(), NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::descentComparison(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::setupDoubleDogleg().
|
mutableprivate |
Workspace of length neq_.
Definition at line 966 of file NonlinearSolver.h.
Referenced by NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Vector of residual weights.
These are used to establish useful and informative weighted norms of the residual vector.
Definition at line 975 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::computeResidWts(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::getResidWts(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::residErrorNorm().
|
private |
Norm of the residual at the start of each nonlinear iteration.
Definition at line 978 of file NonlinearSolver.h.
Referenced by NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::descentComparison(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::expectedResidLeg(), NonlinearSolver::operator=(), NonlinearSolver::setupDoubleDogleg(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Norm of the residual after it has been bounded.
Definition at line 981 of file NonlinearSolver.h.
Referenced by NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Norm of the residual at the end of the first leg of the current iteration.
Definition at line 984 of file NonlinearSolver.h.
Referenced by NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Norm of the residual at the end of the first leg of the current iteration.
Definition at line 987 of file NonlinearSolver.h.
Referenced by NonlinearSolver::convergenceCheck(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Norm of the solution update created by the iteration in its raw, undamped form, using the solution norm.
Definition at line 990 of file NonlinearSolver.h.
Referenced by NonlinearSolver::convergenceCheck(), NonlinearSolver::decideStep(), NonlinearSolver::initializeTrustRegion(), NonlinearSolver::operator=(), NonlinearSolver::setupDoubleDogleg(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Norm of the distance to the Cauchy point using the solution norm.
Definition at line 993 of file NonlinearSolver.h.
Referenced by NonlinearSolver::initializeTrustRegion(), NonlinearSolver::operator=(), NonlinearSolver::setupDoubleDogleg(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Norm of the residual for a trial calculation which may or may not be used.
Definition at line 996 of file NonlinearSolver.h.
Referenced by NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Vector of the norm.
Definition at line 999 of file NonlinearSolver.h.
|
mutableprivate |
Boolean indicating whether we should scale the residual.
Definition at line 1002 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doNewtonSolve(), NonlinearSolver::doResidualCalc(), and NonlinearSolver::operator=().
|
private |
Bounds vector for each species.
Definition at line 1009 of file NonlinearSolver.h.
Referenced by NonlinearSolver::boundStep(), NonlinearSolver::highBoundsConstraintVector(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::setBoundsConstraints().
|
private |
Lower bounds vector for each species.
Definition at line 1012 of file NonlinearSolver.h.
Referenced by NonlinearSolver::boundStep(), NonlinearSolver::lowBoundsConstraintVector(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::setBoundsConstraints().
|
private |
Damping factor imposed by hard bounds and by delta bounds.
Definition at line 1015 of file NonlinearSolver.h.
Referenced by NonlinearSolver::convergenceCheck(), NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Additional damping factor due to bounds on the residual and solution norms.
Definition at line 1018 of file NonlinearSolver.h.
Referenced by NonlinearSolver::convergenceCheck(), NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Delta t for the current step.
Definition at line 1021 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac(), NonlinearSolver::calc_ydot(), NonlinearSolver::doResidualCalc(), and NonlinearSolver::operator=().
|
mutableprivate |
Counter for the total number of function evaluations.
Definition at line 1024 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac(), NonlinearSolver::doResidualCalc(), and NonlinearSolver::operator=().
|
private |
The type of column scaling used in the matrix inversion of the problem.
If 1 then colScaling = m_ewt[] If 2 then colScaling = user set if 0 then colScaling = 1.0
Definition at line 1036 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcColumnScales(), NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::doNewtonSolve(), NonlinearSolver::operator=(), NonlinearSolver::scaleMatrix(), and NonlinearSolver::setColumnScaling().
|
private |
int indicating whether row scaling is turned on (1) or not (0)
Definition at line 1039 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::doNewtonSolve(), NonlinearSolver::operator=(), NonlinearSolver::scaleMatrix(), and NonlinearSolver::setRowScaling().
|
private |
Total number of linear solves taken by the solver object.
Definition at line 1042 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doNewtonSolve(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Number of local linear solves done during the current iteration.
Definition at line 1045 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doNewtonSolve(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Total number of Newton iterations.
Definition at line 1048 of file NonlinearSolver.h.
Referenced by NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
int m_min_newt_its |
Minimum number of Newton iterations to use.
Definition at line 1052 of file NonlinearSolver.h.
Referenced by NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Maximum number of Newton iterations.
Definition at line 1056 of file NonlinearSolver.h.
Referenced by NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Jacobian formation method.
1 = numerical (default) 2 = analytical
Definition at line 1063 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac(), and NonlinearSolver::operator=().
|
private |
Number of Jacobian evaluations.
Definition at line 1066 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac(), and NonlinearSolver::operator=().
|
private |
Current system time.
Note, we assume even for steady state problems that the residual is a function of a system time.
Definition at line 1073 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcColumnScales(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Boolean indicating matrix conditioning.
Definition at line 1076 of file NonlinearSolver.h.
Referenced by NonlinearSolver::operator=().
|
private |
Order of the time step method = 1.
Definition at line 1079 of file NonlinearSolver.h.
Referenced by NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::operator=(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
value of the relative tolerance to use in solving the equation set
Definition at line 1082 of file NonlinearSolver.h.
Referenced by NonlinearSolver::computeResidWts(), NonlinearSolver::createSolnWeights(), NonlinearSolver::decideStep(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), and NonlinearSolver::setRtol().
|
private |
Base value of the absolute tolerance.
Definition at line 1085 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::computeResidWts(), NonlinearSolver::NonlinearSolver(), and NonlinearSolver::operator=().
|
private |
absolute tolerance in the solution unknown
This is used to evaluating the weighting factor
Definition at line 1091 of file NonlinearSolver.h.
Referenced by NonlinearSolver::createSolnWeights(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::setAtol(), and NonlinearSolver::setDefaultDeltaBoundsMagnitudes().
|
private |
absolute tolerance in the unscaled solution unknowns
Definition at line 1094 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::computeResidWts(), NonlinearSolver::operator=(), and NonlinearSolver::setResidualTols().
|
private |
absolute tolerance in the unscaled solution unknowns
Definition at line 1097 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::computeResidWts(), NonlinearSolver::operator=(), and NonlinearSolver::setResidualTols().
|
private |
Check the residual tolerances explicitly against user input.
0 Don't calculate residual weights from residual tolerance inputs 1 Calculate residual weights from residual tolerance inputs only 2 Calculate residual weights from a minimum of the solution error weights process and the direct residual tolerance inputs
Definition at line 1105 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::computeResidWts(), NonlinearSolver::operator=(), and NonlinearSolver::setResidualTols().
|
private |
Determines the level of printing for each time step.
0 -> absolutely nothing is printed for a single time step. 1 -> One line summary per solve_nonlinear call 2 -> short description, points of interest: Table of nonlinear solve - one line per iteration 3 -> Table is included -> More printing per nonlinear iteration (default) that occurs during the table 4 -> Summaries of the nonlinear solve iteration as they are occurring -> table no longer printed Base_ShowSolution Residual called for residual printing at the end of convergence. 5 -> Algorithm information on the nonlinear iterates are printed out 6 -> Additional info on the nonlinear iterates are printed out Base_ShowSolution Residual called for residual printing at the end of each step. 7 -> Additional info on the linear solve is printed out. 8 -> Info on a per iterate of the linear solve is printed out.
Definition at line 1121 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac(), NonlinearSolver::boundStep(), NonlinearSolver::dampDogLeg(), NonlinearSolver::dampStep(), NonlinearSolver::decideStep(), NonlinearSolver::deltaBoundStep(), NonlinearSolver::descentComparison(), NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::initializeTrustRegion(), NonlinearSolver::operator=(), NonlinearSolver::readjustTrustVector(), NonlinearSolver::residErrorNorm(), NonlinearSolver::residualComparisonLeg(), NonlinearSolver::setPrintLvl(), NonlinearSolver::solnErrorNorm(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Scale factor for turning residual norms into solution norms.
Definition at line 1124 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::computeResidWts(), and NonlinearSolver::operator=().
|
private |
Copy of the Jacobian that doesn't get overwritten when the inverse is determined.
The Jacobian stored here is the raw matrix, before any row or column scaling is carried out
Definition at line 1130 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac(), NonlinearSolver::calcSolnToResNormVector(), NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::operator=(), NonlinearSolver::solve_nonlinear_problem(), and NonlinearSolver::~NonlinearSolver().
|
private |
Hessian.
Definition at line 1133 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::operator=(), and NonlinearSolver::~NonlinearSolver().
|
private |
Steepest descent direction. This is also the distance to the Cauchy Point.
Definition at line 1140 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcTrustIntersection(), NonlinearSolver::decideStep(), NonlinearSolver::descentComparison(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::fillDogLegStep(), NonlinearSolver::initializeTrustRegion(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::setupDoubleDogleg().
|
private |
Newton Step - This is the Newton step determined from the straight Jacobian.
Definition at line 1146 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcTrustIntersection(), NonlinearSolver::descentComparison(), NonlinearSolver::fillDogLegStep(), NonlinearSolver::initializeTrustRegion(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::residualComparisonLeg(), NonlinearSolver::setupDoubleDogleg(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Expected value of the residual norm at the Cauchy point if the quadratic model were valid.
Definition at line 1150 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doCauchyPointSolve(), and NonlinearSolver::operator=().
|
private |
Current leg.
Definition at line 1153 of file NonlinearSolver.h.
Referenced by NonlinearSolver::dampDogLeg(), NonlinearSolver::decideStep(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Current Alpha param along the leg.
Definition at line 1156 of file NonlinearSolver.h.
Referenced by NonlinearSolver::dampDogLeg(), NonlinearSolver::decideStep(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Residual dot Jd norm.
This is equal to R_hat dot J_hat d_y_descent
Definition at line 1162 of file NonlinearSolver.h.
Referenced by NonlinearSolver::decideStep(), NonlinearSolver::descentComparison(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::expectedResidLeg(), and NonlinearSolver::operator=().
|
private |
Value of lambdaStar_ which is used to calculate the Cauchy point.
Definition at line 1165 of file NonlinearSolver.h.
Referenced by NonlinearSolver::decideStep(), NonlinearSolver::descentComparison(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::expectedResidLeg(), and NonlinearSolver::operator=().
|
private |
Jacobian times the steepest descent direction in the normalized coordinates.
This is equal to [ Jhat d^y_{descent} ] in the notes, Eqn. 18.
Definition at line 1171 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::NonlinearSolver(), and NonlinearSolver::operator=().
|
private |
Vector of trust region values.
Definition at line 1174 of file NonlinearSolver.h.
Referenced by NonlinearSolver::adjustUpStepMinimums(), NonlinearSolver::calcTrustDistance(), NonlinearSolver::calcTrustIntersection(), NonlinearSolver::initializeTrustRegion(), NonlinearSolver::NonlinearSolver(), NonlinearSolver::operator=(), NonlinearSolver::readjustTrustVector(), NonlinearSolver::solve_nonlinear_problem(), and NonlinearSolver::trustRegionLength().
|
mutableprivate |
Current norm of the vector deltaX_trust_ in terms of the solution norm.
Definition at line 1177 of file NonlinearSolver.h.
Referenced by NonlinearSolver::operator=(), NonlinearSolver::readjustTrustVector(), NonlinearSolver::solve_nonlinear_problem(), and NonlinearSolver::trustRegionLength().
|
private |
Current value of trust radius.
This is used with deltaX_trust_ to calculate the max step size.
Definition at line 1181 of file NonlinearSolver.h.
Referenced by NonlinearSolver::adjustUpStepMinimums(), NonlinearSolver::calcTrustDistance(), NonlinearSolver::dampDogLeg(), NonlinearSolver::decideStep(), NonlinearSolver::initializeTrustRegion(), NonlinearSolver::operator=(), NonlinearSolver::readjustTrustVector(), NonlinearSolver::solve_nonlinear_problem(), and NonlinearSolver::trustRegionLength().
|
private |
Method for handling the trust region initialization.
Then, on subsequent invocations of solve_nonlinear_problem() the strategy flips to method 0.
method Method to set the strategy 0 No strategy - Use the previous strategy 1 Factor of the solution error weights 2 Factor of the first Cauchy Point distance 3 Factor of the first Newton step distance
Definition at line 1193 of file NonlinearSolver.h.
Referenced by NonlinearSolver::initializeTrustRegion(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Factor used to set the initial trust region.
Definition at line 1196 of file NonlinearSolver.h.
Referenced by NonlinearSolver::initializeTrustRegion(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Relative distance down the Newton step that the second dogleg starts.
Definition at line 1199 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcTrustIntersection(), NonlinearSolver::expectedResidLeg(), NonlinearSolver::fillDogLegStep(), NonlinearSolver::operator=(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::setupDoubleDogleg().
|
private |
Distance of the zeroeth leg of the dogleg in terms of the solution norm.
Definition at line 1202 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcTrustIntersection(), NonlinearSolver::lambdaToLeg(), NonlinearSolver::operator=(), and NonlinearSolver::setupDoubleDogleg().
|
private |
Distance of the first leg of the dogleg in terms of the solution norm.
Definition at line 1205 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcTrustIntersection(), NonlinearSolver::lambdaToLeg(), NonlinearSolver::operator=(), and NonlinearSolver::setupDoubleDogleg().
|
private |
Distance of the second leg of the dogleg in terms of the solution norm.
Definition at line 1208 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcTrustIntersection(), NonlinearSolver::lambdaToLeg(), NonlinearSolver::operator=(), and NonlinearSolver::setupDoubleDogleg().
|
private |
Distance of the sum of all legs of the doglegs in terms of the solution norm.
Definition at line 1211 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcTrustIntersection(), NonlinearSolver::lambdaToLeg(), NonlinearSolver::operator=(), and NonlinearSolver::setupDoubleDogleg().
|
private |
Dot product of the Jd_ variable defined above with itself.
Definition at line 1214 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doCauchyPointSolve(), and NonlinearSolver::operator=().
|
private |
Norm of the Newton Step wrt trust region.
Definition at line 1217 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcTrustIntersection(), NonlinearSolver::dampDogLeg(), NonlinearSolver::operator=(), and NonlinearSolver::setupDoubleDogleg().
|
private |
Norm of the Cauchy Step direction wrt trust region.
Definition at line 1220 of file NonlinearSolver.h.
Referenced by NonlinearSolver::calcTrustIntersection(), NonlinearSolver::operator=(), and NonlinearSolver::setupDoubleDogleg().
|
private |
General toggle for turning on dog leg damping.
Definition at line 1223 of file NonlinearSolver.h.
Referenced by NonlinearSolver::descentComparison(), NonlinearSolver::doAffineNewtonSolve(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::initializeTrustRegion(), NonlinearSolver::operator=(), NonlinearSolver::readjustTrustVector(), NonlinearSolver::residualComparisonLeg(), NonlinearSolver::setSolverScheme(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
General toggle for turning on Affine solve with Hessian.
Definition at line 1226 of file NonlinearSolver.h.
Referenced by NonlinearSolver::operator=(), NonlinearSolver::setSolverScheme(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Condition number of the matrix.
Definition at line 1229 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doAffineNewtonSolve(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Factor indicating how much trust region has been changed this iteration - output variable.
Definition at line 1232 of file NonlinearSolver.h.
Referenced by NonlinearSolver::decideStep(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Factor indicating how much trust region has been changed next iteration - output variable.
Definition at line 1235 of file NonlinearSolver.h.
Referenced by NonlinearSolver::decideStep(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Boolean indicating that the residual weights have been reevaluated this iteration - output variable.
Definition at line 1238 of file NonlinearSolver.h.
Referenced by NonlinearSolver::computeResidWts(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Expected DResid_dS for the steepest descent path - output variable.
Definition at line 1241 of file NonlinearSolver.h.
Referenced by NonlinearSolver::descentComparison(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Actual DResid_dS for the steepest descent path - output variable.
Definition at line 1244 of file NonlinearSolver.h.
Referenced by NonlinearSolver::descentComparison(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Expected DResid_dS for the Newton path - output variable.
Definition at line 1247 of file NonlinearSolver.h.
Referenced by NonlinearSolver::descentComparison(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
private |
Actual DResid_dS for the Newton path - output variable.
Definition at line 1250 of file NonlinearSolver.h.
Referenced by NonlinearSolver::descentComparison(), NonlinearSolver::operator=(), and NonlinearSolver::solve_nonlinear_problem().
|
static |
Turn off printing of time.
Necessary to do for test suites
Definition at line 1261 of file NonlinearSolver.h.
Referenced by NonlinearSolver::solve_nonlinear_problem().
|
static |
Turn on or off printing of the Jacobian.
Definition at line 1264 of file NonlinearSolver.h.
Referenced by NonlinearSolver::beuler_jac().
|
static |
Turn on extra printing of dogleg information.
Definition at line 1267 of file NonlinearSolver.h.
Referenced by NonlinearSolver::descentComparison(), NonlinearSolver::doCauchyPointSolve(), NonlinearSolver::residualComparisonLeg(), and NonlinearSolver::solve_nonlinear_problem().
|
static |
Turn on solving both the Newton and Hessian systems and comparing the results.
This is off by default
Definition at line 1273 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doAffineNewtonSolve().
|
static |
This toggle turns off the use of the Hessian when it is warranted by the condition number.
This is a debugging option.
Definition at line 1279 of file NonlinearSolver.h.
Referenced by NonlinearSolver::doAffineNewtonSolve().