Aristos  Development
 All Classes Functions Pages
Public Member Functions | List of all members
Aristos::SQPAlgo Class Reference

Implements a trust-region (TR) sequential quadratic programming (SQP) algorithm for equality-constrained optimization. More...

#include <Aristos_SQPAlgoDecl.hpp>

Public Member Functions

 SQPAlgo (DataPool &dat, Objective &obj, Constraints &constr, HessVec &hessvec, LagMult &lagmult, FeasStep &feasstep)
 
bool run (Vector &x, Vector &c, Vector &l, int &iter, int &iflag, Teuchos::ParameterList &parlist)
 Runs the trust-region SQP algorithm. More...
 
bool runMultiplierEst (const Vector &x, const Vector &g, Vector &l, double tol)
 Computes a Lagrange multiplier estimate. More...
 
bool runQuasiNormalStep (const Vector &x, const Vector &c, Vector &s, double delta, double tol)
 Computes the quasi-normal step. More...
 
bool runTangentialStep (const Vector &x, const Vector &g, const Vector &v, const Vector &l, Vector &s, double delta, int cgitmax, double cgtol, double tol, int &cgiter, int &iflag)
 Computes the tangential step. More...
 
bool runTangentialStepInx (const Vector &x, const Vector &g, const Vector &v, const Vector &l, Vector &s, double delta, int cgitmax, double cgtol, double fixedtol, bool istolfixed, bool fullortho, bool orthocheck, bool fcdcheck, int &cgiter, int &iflag)
 Computes the tangential step using the inexact CG algorithm with full orthogonalization and inexactness control. More...
 
bool runAcceptStep (Vector &x, const Vector &l, double f, const Vector &g, const Vector &c, const Vector &v, const Vector &w, double &delta, double &rho, double tol, int &iaccept)
 Checks whether the computed step is acceptable, and accordingly updates the trust-region radius and the merit function penalty parameter. More...
 
bool runAcceptStepInx (Vector &x, const Vector &l, double f, const Vector &g, const Vector &c, const Vector &v, Vector &w, double &delta, double &rho, bool addproj, double &projtol, double tol, int &iaccept)
 Checks whether the computed step is acceptable, and accordingly updates the trust-region radius and the merit function penalty parameter. Should be used whenever the SQP algorithm manages inexactness in linear system solves. More...
 
bool runDerivativeCheck (const Vector &x, const Vector &l, const Vector &dir)
 Derivative checks.. More...
 

Detailed Description

Implements a trust-region (TR) sequential quadratic programming (SQP) algorithm for equality-constrained optimization.

Aristos::SQPAlgo solves equality-constrained nonlinear optimization problems of the type:

\[ \min f(x) \]

\[ \mbox{subject to } c(x) = 0 \]

  where \form#6, \form#7,
  and \form#8 and \form#9 are Hilbert spaces. We utilize a trust-region
  SQP approach in which the quadratic subproblems are solved using a composite step strategy.
  The TR quadratic subproblems are of the form

\[ \min \frac{1}{2} \langle H_k s^x_k, s^x_k \rangle_{\mathcal{X}} + \langle \nabla_x \mathcal{L} (x_k, \lambda_k), s^x_k \rangle_{\mathcal{X}} \]

subject to

\[ c_x(x_k)s^x_k + c(x_k) = 0 \]

\[ \|s^x_k\|_{\mathcal{X}} \le \Delta_k \]

where $\mathcal{L} (x_k, \lambda_k)$ is the Lagrangian functional with multipliers $\lambda_k$, $H_k$ is a representation of the Hessian of the Lagrangian at $x_k$, and $\Delta_k$ is the trust-region radius at iteration $k$. The step $s_k$ consists of a quasi-normal step that moves the model toward feasibility and a tangential step that moves the model toward optimality while staying in the null space of the linearized constraints (at least in absence of inexact linear system solves). Based on this approach, the SQP algorithm consists of the following steps:

Member Function Documentation

bool Aristos::SQPAlgo::run ( Vector x,
Vector c,
Vector l,
int &  iter,
int &  iflag,
Teuchos::ParameterList &  parlist 
)

Runs the trust-region SQP algorithm.

Parameters
x[in,out] - Initial SQP iterate vector, final solution vector.
c[out] - The vector of constraint values.
l[out] - The vector of Lagrange multipliers.
iter[out] - Total number of SQP iterations.
iflag[out] - Return code: if 0, the algorithm has converged to desired tolerances within the maximum number of iterations; if 1, the maximum number of iterations has been exceeded; if 2, the trust-region radius has fallen below its allowed minimum.
parlist[in] - Parameter list.
Returns
true if the SQP run is successful.
Detailed Description:

The SQP algorithm consists of the following steps:

  • Lagrange multiplier estimation,
  • computation of the quasi-normal step,
  • computation of the tangential step,
  • decision on whether to accept the computed total step and modify the trust-region radius.

Please read the documentation provided for the corresponding member functions.

References Aristos::Vector::createVector(), Aristos::Vector::innerProd(), runAcceptStepInx(), runMultiplierEst(), runQuasiNormalStep(), and runTangentialStepInx().

bool Aristos::SQPAlgo::runAcceptStep ( Vector x,
const Vector l,
double  f,
const Vector g,
const Vector c,
const Vector v,
const Vector w,
double &  delta,
double &  rho,
double  tol,
int &  iaccept 
)

Checks whether the computed step is acceptable, and accordingly updates the trust-region radius and the merit function penalty parameter.

Parameters
x[in,out] - SQP iterate.
l[in] - Current Lagrange multiplier estimate.
f[in] - Current value of the objective.
g[in] - Current gradient of the objective vector.
c[in] - The vector of constraint values.
v[in] - Current quasi-normal step.
w[in] - Current tangential step.
s[out] - The computed tangential step vector.
delta[in,out] - Trust-region radius.
rho[in,out] - Penalty parameter.
tol[in] - Tolerance for inexact computations.
iaccept[out] - Integer return code: 0 - step rejected 1 - step accepted
Returns
true if runTangentialStep terminates successfully.
Detailed Description:

Global convergence in trust-region SQP methods is usually ensured through the use of a merit function, which helps us determine whether a step is acceptable and whether the trust-region radius needs to be modified. We use the augmented Lagrangian merit function

\begin{equation*} \phi(x,\lambda;\rho) = f(x) + \langle \lambda, c(x) \rangle_{\mathcal{Y}} + \rho \|c(x)\|^2_{\mathcal{Y}} = \mathcal{L}(x,\lambda) + \rho \|c(x)\|^2_{\mathcal{Y}}. \end{equation*}

Let $s_k^x=n_k+t_k$ be a trial step where $n_k$ is the computed quasi-normal step, and $t_k$ is the computed tangential step, and let $\lambda_{k+1} = \lambda_{k} + \Delta\lambda_k$ be an updated Lagrange multiplier. To measure the improvement in the merit function $\phi$, we compare the actual reduction and the predicted reduction in moving from the current iterate $x_k$ to the trial iterate $x_k+s_k^x$. The actual reduction is defined by

\begin{equation*} ared(s_k^x; \rho_k) = \phi(x_k,\lambda_k;\rho_k) - \phi(x_k+s_k,\lambda_{k+1};\rho_k), \end{equation*}

and the predicted reduction is given by

\begin{equation*} pred(s_k^x; \rho_k) = \phi(x_k,\lambda_k;\rho_k) - \widetilde{\phi}(x_k,\Delta\lambda_k;\rho_k), \end{equation*}

where

\begin{equation*} \widetilde{\phi}(x_k,\Delta\lambda_k;\rho_k) = \frac{1}{2} \langle H_k s_k^x, s_k^x \rangle_{\mathcal{X}} + \langle \nabla_x \mathcal{L}(x_k,\lambda_k), s_k^x \rangle_{\mathcal{X}} + \langle \Delta\lambda_k, c_x(x_k)s_k^x+c(x_k) \rangle_{\mathcal{Y}} + \rho_k \|c_x(x_k)s_k^x+c(x_k)\|^2_{\mathcal{Y}} . \end{equation*}

The following decision procedure is then used (we leave out the details of the penalty parameter update):

\begin{enumerate} \addtocounter{enumi}{-1} \item Assume $0 < \eta_1 < \eta_2 < 1$, $0 < \gamma_1 < \gamma_2 < 1$. \item Acceptance Test. \begin{enumerate} \item Update penalty parameter $\rho_k$ (details omitted). \item Compute actual reduction $ared(s_k^x; \rho_k)$ and predicted reduction $pred(s_k^x; \rho_k)$, and their ratio $\theta_k=\frac{ared(s_k^x; \rho_k)}{pred(s_k^x; \rho_k)}$. \item If $\theta_k \ge \eta_1$, set $x_{k+1} = x_k+s_k$, otherwise set $x_{k+1} = x_k$ and reset $\lambda_{k+1} = \lambda_k$. \end{enumerate} \item Trust-Region Radius Update. Set \[ \Delta_{k+1} \in \begin{cases} [ \Delta_k, \infty] & \mbox{ if } \; \theta_k \ge \eta_2, \\ [ \gamma_2 \Delta_k, \Delta_k ] & \mbox{ if } \; \theta_k \in [\eta_1, \eta_2), \\ [ \gamma_1 \Delta_k, \gamma_2 \Delta_k ] & \mbox{ if } \; \theta_k < \eta_1. \end{cases} \] \end{enumerate}

References Aristos::Vector::createVector(), Aristos::Vector::innerProd(), Aristos::Vector::linComb(), and runMultiplierEst().

bool Aristos::SQPAlgo::runAcceptStepInx ( Vector x,
const Vector l,
double  f,
const Vector g,
const Vector c,
const Vector v,
Vector w,
double &  delta,
double &  rho,
bool  addproj,
double &  projtol,
double  tol,
int &  iaccept 
)

Checks whether the computed step is acceptable, and accordingly updates the trust-region radius and the merit function penalty parameter. Should be used whenever the SQP algorithm manages inexactness in linear system solves.

Parameters
x[in,out] - SQP iterate.
l[in] - Current Lagrange multiplier estimate.
f[in] - Current value of the objective.
g[in] - Current gradient of the objective vector.
c[in] - The vector of constraint values.
v[in] - Current quasi-normal step.
w[in] - Current tangential step.
s[out] - The computed tangential step vector.
delta[in,out] - Trust-region radius.
rho[in,out] - Penalty parameter.
addproj[in]- if true, run an additional null-space projection
projtol[in,out]- Tolerance for inexact null-space projections.
tol[in] - Tolerance for inexact derivative, etc computations.
iaccept[out] - Integer return code: 0 - step rejected 1 - step accepted
Returns
true if runTangentialStepInx terminates successfully.
Detailed Description:

Global convergence in trust-region SQP methods is usually ensured through the use of a merit function, which helps us determine whether a step is acceptable and whether the trust-region radius needs to be modified. We use the augmented Lagrangian merit function

\begin{equation*} \phi(x,\lambda;\rho) = f(x) + \langle \lambda, c(x) \rangle_{\mathcal{Y}} + \rho \|c(x)\|^2_{\mathcal{Y}} = \mathcal{L}(x,\lambda) + \rho \|c(x)\|^2_{\mathcal{Y}}. \end{equation*}

Let $s_k^x=n_k+t_k$ be a trial step where $n_k$ is the computed quasi-normal step, and $t_k$ is the computed tangential step, and let $\lambda_{k+1} = \lambda_{k} + \Delta\lambda_k$ be an updated Lagrange multiplier. To measure the improvement in the merit function $\phi$, we compare the actual reduction and the predicted reduction in moving from the current iterate $x_k$ to the trial iterate $x_k+s_k^x$. The actual reduction is defined by

\begin{equation*} ared(s_k^x; \rho_k) = \phi(x_k,\lambda_k;\rho_k) - \phi(x_k+s_k,\lambda_{k+1};\rho_k), \end{equation*}

and the predicted reduction is given by

\begin{equation*} pred(s_k^x; \rho_k) = \phi(x_k,\lambda_k;\rho_k) - \widetilde{\phi}(x_k,\Delta\lambda_k;\rho_k), \end{equation*}

where

\begin{equation*} \widetilde{\phi}(x_k,\Delta\lambda_k;\rho_k) = \frac{1}{2} \langle H_k s_k^x, s_k^x \rangle_{\mathcal{X}} + \langle \nabla_x \mathcal{L}(x_k,\lambda_k), s_k^x \rangle_{\mathcal{X}} + \langle \Delta\lambda_k, c_x(x_k)s_k^x+c(x_k) \rangle_{\mathcal{Y}} + \rho_k \|c_x(x_k)s_k^x+c(x_k)\|^2_{\mathcal{Y}} . \end{equation*}

The following decision procedure is then used (we leave out the details of the penalty parameter update):

\begin{enumerate} \addtocounter{enumi}{-1} \item Assume $0 < \eta_1 < \eta_2 < 1$, $0 < \gamma_1 < \gamma_2 < 1$. \item Acceptance Test. \begin{enumerate} \item Update penalty parameter $\rho_k$ (details omitted). \item Compute actual reduction $ared(s_k^x; \rho_k)$ and predicted reduction $pred(s_k^x; \rho_k)$, and their ratio $\theta_k=\frac{ared(s_k^x; \rho_k)}{pred(s_k^x; \rho_k)}$. \item If $\theta_k \ge \eta_1$, set $x_{k+1} = x_k+s_k$, otherwise set $x_{k+1} = x_k$ and reset $\lambda_{k+1} = \lambda_k$. \end{enumerate} \item Trust-Region Radius Update. Set \[ \Delta_{k+1} \in \begin{cases} [ \Delta_k, \infty] & \mbox{ if } \; \theta_k \ge \eta_2, \\ [ \gamma_2 \Delta_k, \Delta_k ] & \mbox{ if } \; \theta_k \in [\eta_1, \eta_2), \\ [ \gamma_1 \Delta_k, \gamma_2 \Delta_k ] & \mbox{ if } \; \theta_k < \eta_1. \end{cases} \] \end{enumerate}

References Aristos::Vector::createVector(), Aristos::Vector::innerProd(), Aristos::Vector::linComb(), and runMultiplierEst().

Referenced by run().

bool Aristos::SQPAlgo::runDerivativeCheck ( const Vector x,
const Vector l,
const Vector dir 
)

Derivative checks..

Parameters
x[in] - Iterate vector.
l[in] - Lagrange multiplier vector.
dir[in] - Direction of differentiation.
Returns
true if runDerivativeCheck terminates successfully.
Detailed Description:

Checks validity of the computation of the following quantities:

  • objective gradient
  • constraint Jacobian
  • Hessian of the Lagrangian

It is best to run the checks once before the SQP run, with random input vectors.

References Aristos::Vector::createVector().

bool Aristos::SQPAlgo::runMultiplierEst ( const Vector x,
const Vector g,
Vector l,
double  tol 
)

Computes a Lagrange multiplier estimate.

Parameters
x[in] - Current SQP iterate.
g[in] - Gradient of the objective.
l[out] - The vector of Lagrange multipliers.
tol[in] - Tolerance for inexact computations.
Returns
true if runMultiplierEst is successful.
Detailed Description:

Calls a user-defined multiplier estimate function which is passed in through the data pool object. We might provide a default routine in the future.

A common practice in SQP methods is to approximate the Lagrange multipliers by a least-squares estimate based on the KKT stationarity conditions (eq:eqstationarity}). Thus at each iterate $x_k$, we solve a problem related to the least-squares problem

\begin{eqnarray*} \min_{\lambda} & & \| \nabla_x f(x_k) + c_x(x_k)^* \lambda \|_{\mathcal{X}}. \end{eqnarray*}

INEXACT VERSION:

It should be noted that there are many ways in which Lagrange multiplier estimation problems can be formulated, due to specific problem features and data representations provided by the user. The convergence theory for SQP methods poses rather loose requirements on the Lagrange multipliers. In fact, our SQP algorithm merely requires that the sequence of Lagrange multipliers be uniformly bounded.

This is guaranteed through our choice of the stopping criteria for the iterative linear system solver that computes the Lagrange multipliers.

Referenced by run(), runAcceptStep(), and runAcceptStepInx().

bool Aristos::SQPAlgo::runQuasiNormalStep ( const Vector x,
const Vector c,
Vector s,
double  delta,
double  tol 
)

Computes the quasi-normal step.

Parameters
x[in] - Current SQP iterate.
c[in] - The vector of constraint values.
s[out] - The computed quasi-normal step vector.
delta[in] - Scaled trust-region radius.
tol[in] - Tolerance for inexact computations.
Returns
true if runQuasiNormalStep is successful.
Detailed Description:

We ask that the quasi-normal step $n_k$ approximately solve a problem of the following type:

\begin{eqnarray*} \min & & \|c_x(x_k) n + c(x_k)\|^2_{\mathcal{Y}} \\ \mbox{s.t.} & & \|n\|_{\mathcal{X}} \le \zeta \Delta_k \end{eqnarray*}

The approximate solution is computed using the Powell dogleg method. The dogleg path is computed using the Cauchy point $n^{cp}_k = \alpha^{cp} c_x(x_k)^*c(x_k)$, where the step length $\alpha^{cp}$ is given by

\begin{equation*} \begin{cases} \frac{\| c_x(x_k)^*c(x_k) \|^2_{\mathcal{X}}}{\| c_x(x_k)c_x(x_k)^* c(x_k) \|^2_{\mathcal{Y}}} & \mbox{if } \frac{\| c_x(x_k)^*c(x_k) \|^3_{\mathcal{X}}}{\| c_x(x_k)c_x(x_k)^* c(x_k) \|^2_{\mathcal{Y}}} \le \zeta \Delta_k \\ \frac{\zeta \Delta_k}{\| c_x(x_k)^* c(x_k) \|_{\mathcal{X}}} & \mbox{otherwise}, \end{cases} \end{equation*}

and the feasibility step $n^N_k$ given, for example, as the minimum norm minimizer of the 2-norm of the linearized constraints (other choices are possible):

\begin{equation*} n^N_k = - c_x(x_k)^* \left( c_x(x_k)c(x_k)^* \right)^{-1} c(x_k). \end{equation*}

The intersection of the dogleg path (origin - Cauchy point - feasibility point) with the trust-region constraint gives the quasi-normal step.

INEXACT VERSION:

We require from the quasi-normal step $n_k$ to give as much decrease as $n^{cp}_k$, and

\begin{equation*} \|c(x_k)\|^2_{\mathcal{Y}} - \|c_x(x_k)n_k+c(x_k)\|^2_{\mathcal{Y}} \ge \sigma_1 \left( \|c(x_k)\|^2_{\mathcal{Y}} - \|c_x(x_k)n^{cp}_k+c(x_k)\|^2_{\mathcal{Y}} \right), \end{equation*}

for some $0 < \sigma_1 \le 1$ independent of $k$. Additionally, due to possible nonnormality of $n_k$ with respect to the tangent space of the constraints, the following condition must hold for $\kappa_1>0$ independent of $k$:

\begin{equation*} \|n_k\|_{\mathcal{X}} \le \kappa_1 \|c(x_k)\|_{\mathcal{Y}}. \end{equation*}

All of the above is guaranteed through our choice of the stopping criteria for an iterative linear system solver that computes the feasibility step $n^N_k$.

References Aristos::Vector::createVector(), Aristos::Vector::innerProd(), Aristos::Vector::linComb(), and Aristos::Vector::Set().

Referenced by run().

bool Aristos::SQPAlgo::runTangentialStep ( const Vector x,
const Vector g,
const Vector v,
const Vector l,
Vector s,
double  delta,
int  cgitmax,
double  cgtol,
double  tol,
int &  cgiter,
int &  iflag 
)

Computes the tangential step.

Parameters
x[in] - Current SQP iterate.
g[in] - Current gradient of the objective vector.
v[in] - Current quasi-normal step.
l[in] - Current Lagrange multiplier estimate.
s[out] - The computed tangential step vector.
delta[in] - Trust-region radius.
cgitmax[in] - Maximum number of conjugate gradient iterations (internal use only).
cgtol[in] - Stopping tolerance for the conjugate gradient method (internal use only).
tol[in] - Tolerance for inexact computations.
cgiter[out] - Required number of conjugate gradient iterations.
iflag[out] - Return flag: 0 - normal convergence 1 - negative curvature detected 2 - trust-region radius active 3 - maximum number of conjugate gradient iterations exceeded
Returns
true upon successful termination (not the same as convergence!) of runTangentialStep.
Detailed Description:

The standard requirements on the tangential step are that it lie in the tangent space of the constraints and that it move the quadratic model of the Lagrangian toward optimality. Thus we ask that the tangential step $t_k$ approximately solve a problem of the type

\begin{eqnarray*} \min & & \frac{1}{2} \langle H_k (t+n_k), t+n_k \rangle_{\mathcal{X}} + \langle \nabla_x \mathcal{L} (x_k, \lambda_k), t+n_k \rangle_{\mathcal{X}} \\ \mbox{s.t.} & & c_x(x_k) t = 0 \\ & & \|t+n_k\|_{\mathcal{X}} \le \Delta_k. \end{eqnarray*}

It is solved using the Steihaug variant of the conjugate gradient (CG) method. We use a projected CG method in the full space, rather than explicitly computing the null space of $c_x(x_k)$. This method relies on applications of the orthogonal projection operator $P_k:\mathcal{X} \rightarrow \mathcal{X}$, defined by

\begin{equation*} P_k = I - c_x(x_k)^* \left( c_x(x_k) c_x(x_k)^* \right) ^{-1} c_x(x_k), \end{equation*}

where $I:\mathcal{X} \rightarrow \mathcal{X}$ is the identity operator.

INEXACT VERSION:

As in the case of the quasi-normal step computation, the tangential step only needs to satisfy a fraction of the Cauchy decrease condition. The model function is in this case related to a quadratic model of the Lagrangian functional (we leave out the details).

Convergence in presence of inexactness is guaranteed through our choice of the stopping criteria for an iterative linear system solver that applies the projection operator $P_k$.

References Aristos::Vector::createVector(), Aristos::Vector::innerProd(), Aristos::Vector::linComb(), and Aristos::Vector::Set().

bool Aristos::SQPAlgo::runTangentialStepInx ( const Vector x,
const Vector g,
const Vector v,
const Vector l,
Vector s,
double  delta,
int  cgitmax,
double  cgtol,
double  fixedtol,
bool  istolfixed,
bool  fullortho,
bool  orthocheck,
bool  fcdcheck,
int &  cgiter,
int &  iflag 
)

Computes the tangential step using the inexact CG algorithm with full orthogonalization and inexactness control.

Parameters
x[in] - Current SQP iterate.
g[in] - Current gradient of the objective vector.
v[in] - Current quasi-normal step.
l[in] - Current Lagrange multiplier estimate.
s[out] - The computed tangential step vector.
delta[in] - Trust-region radius.
cgitmax[in] - Maximum number of conjugate gradient iterations (internal use only).
cgtol[in] - Stopping tolerance for the conjugate gradient method (internal use only).
fixedtol[in] - Tolerance for inexact computations.
istolfixed[in] - true, if fixed inexact solver tolerances are desired false, if the SQP algorithm is managing inexactness
fullortho[in] - true, if full orthogonalization of search directions is desired false, if standard CG 2-vector orthogonalization is sufficient
orthocheck[in] - true, if orthoganality of projected residuals is to be verified false, otherwise
fcdcheck[in] - true, if the fraction of Cauchy decrease condition is to be verified false, otherwise
cgiter[out] - Returns the required number of conjugate gradient iterations.
iflag[out] - Return flag: 0 - normal convergence 1 - negative curvature detected 2 - trust-region radius active 3 - maximum number of conjugate gradient iterations exceeded
Returns
true upon successful termination (not the same as convergence!) of runTangentialStep.
Detailed Description:

The standard requirements on the tangential step are that it lie in the tangent space of the constraints and that it move the quadratic model of the Lagrangian toward optimality. Thus we ask that the tangential step $t_k$ approximately solve a problem of the type

\begin{eqnarray*} \min & & \frac{1}{2} \langle H_k (t+n_k), t+n_k \rangle_{\mathcal{X}} + \langle \nabla_x \mathcal{L} (x_k, \lambda_k), t+n_k \rangle_{\mathcal{X}} \\ \mbox{s.t.} & & c_x(x_k) t = 0 \\ & & \|t+n_k\|_{\mathcal{X}} \le \Delta_k. \end{eqnarray*}

It is solved using the Steihaug variant of the conjugate gradient (CG) method. We use a projected CG method in the full space, rather than explicitly computing the null space of $c_x(x_k)$. This method relies on applications of the orthogonal projection operator $P_k:\mathcal{X} \rightarrow \mathcal{X}$, defined by

\begin{equation*} P_k = I - c_x(x_k)^* \left( c_x(x_k) c_x(x_k)^* \right) ^{-1} c_x(x_k), \end{equation*}

where $I:\mathcal{X} \rightarrow \mathcal{X}$ is the identity operator.

The tangential step only needs to satisfy a fraction of the Cauchy decrease condition. The model function is in this case related to a quadratic model of the Lagrangian functional (we leave out the details).

Convergence in presence of inexactness is guaranteed through our choice of the stopping criteria for an iterative linear system solver that applies the projection operator $P_k$.

References Aristos::Vector::createVector(), Aristos::Vector::innerProd(), Aristos::Vector::linComb(), and Aristos::Vector::Set().

Referenced by run().


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