ROL
|
Implements the computation of optimization steps with the Newton primal-dual active set method. More...
#include <ROL_PrimalDualActiveSetStep.hpp>
Classes | |
class | HessianPD |
class | PrecondPD |
Public Member Functions | |
PrimalDualActiveSetStep (ROL::ParameterList &parlist) | |
Constructor. More... | |
void | initialize (Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state) |
Initialize step with bound constraint. More... | |
void | compute (Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state) |
Compute step. More... | |
void | update (Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state) |
Update step, if successful. More... | |
std::string | printHeader (void) const |
Print iterate header. More... | |
std::string | printName (void) const |
Print step name. More... | |
virtual std::string | print (AlgorithmState< Real > &algo_state, bool print_header=false) const |
Print iterate status. More... | |
Public Member Functions inherited from ROL::Step< Real > | |
virtual | ~Step () |
Step (void) | |
virtual void | initialize (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state) |
Initialize step with bound constraint. More... | |
virtual void | initialize (Vector< Real > &x, const Vector< Real > &g, Vector< Real > &l, const Vector< Real > &c, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state) |
Initialize step with equality constraint. More... | |
virtual void | initialize (Vector< Real > &x, const Vector< Real > &g, Vector< Real > &l, const Vector< Real > &c, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state) |
Initialize step with equality constraint. More... | |
virtual void | compute (Vector< Real > &s, const Vector< Real > &x, const Vector< Real > &l, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state) |
Compute step (equality constraints). More... | |
virtual void | update (Vector< Real > &x, Vector< Real > &l, const Vector< Real > &s, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state) |
Update step, if successful (equality constraints). More... | |
virtual void | compute (Vector< Real > &s, const Vector< Real > &x, const Vector< Real > &l, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state) |
Compute step (equality constraints). More... | |
virtual void | update (Vector< Real > &x, Vector< Real > &l, const Vector< Real > &s, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state) |
Update step, if successful (equality constraints). More... | |
const ROL::Ptr< const StepState< Real > > | getStepState (void) const |
Get state for step object. More... | |
void | reset (const Real searchSize=1.0) |
Get state for step object. More... | |
Private Member Functions | |
Real | computeCriticalityMeasure (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &con, Real tol) |
Compute the gradient-based criticality measure. More... | |
Private Attributes | |
ROL::Ptr< Krylov< Real > > | krylov_ |
int | iterCR_ |
CR iteration counter. More... | |
int | flagCR_ |
CR termination flag. More... | |
Real | itol_ |
Inexact CR tolerance. More... | |
int | maxit_ |
Maximum number of PDAS iterations. More... | |
int | iter_ |
PDAS iteration counter. More... | |
int | flag_ |
PDAS termination flag. More... | |
Real | stol_ |
PDAS minimum step size stopping tolerance. More... | |
Real | gtol_ |
PDAS gradient stopping tolerance. More... | |
Real | scale_ |
Scale for dual variables in the active set, \(c\). More... | |
Real | neps_ |
\(\epsilon\)-active set parameter More... | |
bool | feasible_ |
Flag whether the current iterate is feasible or not. More... | |
ROL::Ptr< Vector< Real > > | lambda_ |
Container for dual variables. More... | |
ROL::Ptr< Vector< Real > > | xlam_ |
Container for primal plus dual variables. More... | |
ROL::Ptr< Vector< Real > > | x0_ |
Container for initial priaml variables. More... | |
ROL::Ptr< Vector< Real > > | xbnd_ |
Container for primal variable bounds. More... | |
ROL::Ptr< Vector< Real > > | As_ |
Container for step projected onto active set. More... | |
ROL::Ptr< Vector< Real > > | xtmp_ |
Container for temporary primal storage. More... | |
ROL::Ptr< Vector< Real > > | res_ |
Container for optimality system residual for quadratic model. More... | |
ROL::Ptr< Vector< Real > > | Ag_ |
Container for gradient projected onto active set. More... | |
ROL::Ptr< Vector< Real > > | rtmp_ |
Container for temporary right hand side storage. More... | |
ROL::Ptr< Vector< Real > > | gtmp_ |
Container for temporary gradient storage. More... | |
ESecant | esec_ |
Enum for secant type. More... | |
ROL::Ptr< Secant< Real > > | secant_ |
Secant object. More... | |
bool | useSecantPrecond_ |
bool | useSecantHessVec_ |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::Step< Real > | |
ROL::Ptr< StepState< Real > > | getState (void) |
Implements the computation of optimization steps with the Newton primal-dual active set method.
To describe primal-dual active set (PDAS), we consider the following abstract setting. Suppose \(\mathcal{X}\) is a Hilbert space of functions mapping \(\Xi\) to \(\mathbb{R}\). For example, \(\Xi\subset\mathbb{R}^n\) and \(\mathcal{X}=L^2(\Xi)\) or \(\Xi = \{1,\ldots,n\}\) and \(\mathcal{X}=\mathbb{R}^n\). We assume \( f:\mathcal{X}\to\mathbb{R}\) is twice-continuously Fréchet differentiable and \(a,\,b\in\mathcal{X}\) with \(a\le b\) almost everywhere in \(\Xi\). Note that the PDAS algorithm will also work with secant approximations of the Hessian.
Traditionally, PDAS is an algorithm for the minimizing quadratic objective functions subject to bound constraints. ROL implements a Newton PDAS which extends PDAS to general bound-constrained nonlinear programs, i.e.,
\[ \min_x \quad f(x) \quad \text{s.t.} \quad a \le x \le b. \]
Given the \(k\)-th iterate \(x_k\), the Newton PDAS algorithm computes steps by applying PDAS to the quadratic subproblem
\[ \min_s \quad \langle \nabla^2 f(x_k)s + \nabla f(x_k),s \rangle_{\mathcal{X}} \quad \text{s.t.} \quad a \le x_k + s \le b. \]
For the \(k\)-th quadratic subproblem, PDAS builds an approximation of the active set \(\mathcal{A}_k\) using the dual variable \(\lambda_k\) as
\[ \mathcal{A}^+_k = \{\,\xi\in\Xi\,:\,(\lambda_k + c(x_k-b))(\xi) > 0\,\}, \quad \mathcal{A}^-_k = \{\,\xi\in\Xi\,:\,(\lambda_k + c(x_k-a))(\xi) < 0\,\}, \quad\text{and}\quad \mathcal{A}_k = \mathcal{A}^-_k\cup\mathcal{A}^+_k. \]
We define the inactive set \(\mathcal{I}_k=\Xi\setminus\mathcal{A}_k\). The solution to the quadratic subproblem is then computed iteratively by solving
\[ \nabla^2 f(x_k) s_k + \lambda_{k+1} = -\nabla f(x_k), \quad x_k+s_k = a \;\text{on}\;\mathcal{A}^-_k,\quad x_k+s_k = b\;\text{on}\;\mathcal{A}^+_k, \quad\text{and}\quad \lambda_{k+1} = 0\;\text{on}\;\mathcal{I}_k \]
and updating the active and inactive sets.
One can rewrite this system by consolidating active and inactive parts, i.e.,
\[ \begin{pmatrix} \nabla^2 f(x_k)_{\mathcal{A}_k,\mathcal{A}_k} & \nabla^2 f(x_k)_{\mathcal{A}_k,\mathcal{I}_k} \\ \nabla^2 f(x_k)_{\mathcal{I}_k,\mathcal{A}_k} & \nabla^2 f(x_k)_{\mathcal{I}_k,\mathcal{I}_k} \end{pmatrix} \begin{pmatrix} (s_k)_{\mathcal{A}_k} \\ (s_k)_{\mathcal{I}_k} \end{pmatrix} + \begin{pmatrix} (\lambda_{k+1})_{\mathcal{A}_k} \\ 0 \end{pmatrix} = - \begin{pmatrix} \nabla f(x_k)_{\mathcal{A}_k}\\ \nabla f(x_k)_{\mathcal{I}_k} \end{pmatrix}. \]
Here the subscripts \(\mathcal{A}_k\) and \(\mathcal{I}_k\) denote the active and inactive components, respectively. Moreover, the active components of \(s_k\) are \(s_k(\xi) = a(\xi)-x_k(\xi)\) if \(\xi\in\mathcal{A}^-_k\) and \(s_k(\xi) = b(\xi)-x_k(\xi)\) if \(\xi\in\mathcal{A}^+_k\). Since \((s_k)_{\mathcal{A}_k}\) is fixed, we only need to solve for the inactive components of \(s_k\) which we can do this using conjugate residuals (CR) (i.e., the Hessian operator corresponding to the inactive indices may not be positive definite). Once \((s_k)_{\mathcal{I}_k}\) is computed, it is straight forward to update the dual variables.
Definition at line 133 of file ROL_PrimalDualActiveSetStep.hpp.
|
inline |
Constructor.
[in] | parlist | is a parameter list containing relevent algorithmic information |
[in] | useSecant | is a bool which determines whether or not the algorithm uses a secant approximation of the Hessian |
Definition at line 280 of file ROL_PrimalDualActiveSetStep.hpp.
References ROL::PrimalDualActiveSetStep< Real >::esec_, ROL::PrimalDualActiveSetStep< Real >::gtol_, ROL::PrimalDualActiveSetStep< Real >::krylov_, ROL::PrimalDualActiveSetStep< Real >::maxit_, ROL::PrimalDualActiveSetStep< Real >::scale_, ROL::PrimalDualActiveSetStep< Real >::secant_, ROL::PrimalDualActiveSetStep< Real >::stol_, ROL::StringToESecant(), ROL::PrimalDualActiveSetStep< Real >::useSecantHessVec_, and ROL::PrimalDualActiveSetStep< Real >::useSecantPrecond_.
|
inlineprivate |
Compute the gradient-based criticality measure.
The criticality measure is
\(\|x_k - P_{[a,b]}(x_k-\nabla f(x_k))\|_{\mathcal{X}}\). Here, \(P_{[a,b]}\) denotes the projection onto the bound constraints.
[in] | x | is the current iteration |
[in] | obj | is the objective function |
[in] | con | are the bound constraints |
[in] | tol | is a tolerance for inexact evaluations of the objective function |
Definition at line 262 of file ROL_PrimalDualActiveSetStep.hpp.
References ROL::Step< Real >::getState(), ROL::Objective< Real >::gradient(), ROL::BoundConstraint< Real >::project(), and ROL::PrimalDualActiveSetStep< Real >::xtmp_.
Referenced by ROL::PrimalDualActiveSetStep< Real >::initialize(), and ROL::PrimalDualActiveSetStep< Real >::update().
|
inlinevirtual |
Initialize step with bound constraint.
Reimplemented from ROL::Step< Real >.
Definition at line 320 of file ROL_PrimalDualActiveSetStep.hpp.
References ROL::PrimalDualActiveSetStep< Real >::Ag_, ROL::PrimalDualActiveSetStep< Real >::As_, ROL::Vector< Real >::clone(), ROL::PrimalDualActiveSetStep< Real >::computeCriticalityMeasure(), ROL::Step< Real >::getState(), ROL::AlgorithmState< Real >::gnorm, ROL::PrimalDualActiveSetStep< Real >::gtmp_, ROL::AlgorithmState< Real >::iter, ROL::PrimalDualActiveSetStep< Real >::lambda_, ROL::AlgorithmState< Real >::nfval, ROL::AlgorithmState< Real >::ngrad, ROL::BoundConstraint< Real >::project(), ROL::PrimalDualActiveSetStep< Real >::res_, ROL::PrimalDualActiveSetStep< Real >::rtmp_, ROL::Objective< Real >::update(), ROL::Objective< Real >::value(), ROL::AlgorithmState< Real >::value, ROL::PrimalDualActiveSetStep< Real >::x0_, ROL::PrimalDualActiveSetStep< Real >::xbnd_, ROL::PrimalDualActiveSetStep< Real >::xlam_, ROL::PrimalDualActiveSetStep< Real >::xtmp_, and zero.
|
inlinevirtual |
Compute step.
Reimplemented from ROL::Step< Real >.
Definition at line 380 of file ROL_PrimalDualActiveSetStep.hpp.
References ROL::PrimalDualActiveSetStep< Real >::Ag_, ROL::PrimalDualActiveSetStep< Real >::As_, ROL::PrimalDualActiveSetStep< Real >::flag_, ROL::PrimalDualActiveSetStep< Real >::flagCR_, ROL::BoundConstraint< Real >::getLowerBound(), ROL::Step< Real >::getState(), ROL::BoundConstraint< Real >::getUpperBound(), ROL::AlgorithmState< Real >::gnorm, ROL::PrimalDualActiveSetStep< Real >::gtmp_, ROL::PrimalDualActiveSetStep< Real >::gtol_, ROL::Objective< Real >::hessVec(), ROL::PrimalDualActiveSetStep< Real >::iter_, ROL::AlgorithmState< Real >::iterateVec, ROL::PrimalDualActiveSetStep< Real >::iterCR_, ROL::PrimalDualActiveSetStep< Real >::itol_, ROL::PrimalDualActiveSetStep< Real >::krylov_, ROL::PrimalDualActiveSetStep< Real >::lambda_, ROL::PrimalDualActiveSetStep< Real >::maxit_, ROL::PrimalDualActiveSetStep< Real >::neps_, ROL::Vector< Real >::norm(), ROL::Vector< Real >::plus(), ROL::BoundConstraint< Real >::project(), ROL::BoundConstraint< Real >::pruneActive(), ROL::BoundConstraint< Real >::pruneLowerActive(), ROL::BoundConstraint< Real >::pruneUpperActive(), ROL::PrimalDualActiveSetStep< Real >::res_, ROL::PrimalDualActiveSetStep< Real >::rtmp_, ROL::PrimalDualActiveSetStep< Real >::scale_, ROL::PrimalDualActiveSetStep< Real >::secant_, ROL::PrimalDualActiveSetStep< Real >::stol_, ROL::PrimalDualActiveSetStep< Real >::useSecantHessVec_, ROL::PrimalDualActiveSetStep< Real >::useSecantPrecond_, ROL::PrimalDualActiveSetStep< Real >::x0_, ROL::PrimalDualActiveSetStep< Real >::xbnd_, ROL::PrimalDualActiveSetStep< Real >::xlam_, ROL::PrimalDualActiveSetStep< Real >::xtmp_, zero, and ROL::Vector< Real >::zero().
|
inlinevirtual |
Update step, if successful.
Reimplemented from ROL::Step< Real >.
Definition at line 513 of file ROL_PrimalDualActiveSetStep.hpp.
References ROL::PrimalDualActiveSetStep< Real >::computeCriticalityMeasure(), ROL::PrimalDualActiveSetStep< Real >::feasible_, ROL::PrimalDualActiveSetStep< Real >::flag_, ROL::PrimalDualActiveSetStep< Real >::flagCR_, ROL::Step< Real >::getState(), ROL::AlgorithmState< Real >::gnorm, ROL::PrimalDualActiveSetStep< Real >::gtmp_, ROL::BoundConstraint< Real >::isFeasible(), ROL::AlgorithmState< Real >::iter, ROL::PrimalDualActiveSetStep< Real >::iter_, ROL::AlgorithmState< Real >::iterateVec, ROL::PrimalDualActiveSetStep< Real >::iterCR_, ROL::PrimalDualActiveSetStep< Real >::maxit_, ROL::AlgorithmState< Real >::nfval, ROL::AlgorithmState< Real >::ngrad, ROL::Vector< Real >::norm(), ROL::Vector< Real >::plus(), ROL::PrimalDualActiveSetStep< Real >::secant_, ROL::AlgorithmState< Real >::snorm, ROL::Objective< Real >::update(), ROL::Objective< Real >::value(), and ROL::AlgorithmState< Real >::value.
|
inlinevirtual |
Print iterate header.
This function produces a string containing header information.
Reimplemented from ROL::Step< Real >.
Definition at line 545 of file ROL_PrimalDualActiveSetStep.hpp.
References ROL::PrimalDualActiveSetStep< Real >::maxit_.
Referenced by ROL::PrimalDualActiveSetStep< Real >::print().
|
inlinevirtual |
Print step name.
This function produces a string containing the algorithmic step information.
Reimplemented from ROL::Step< Real >.
Definition at line 572 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::print().
|
inlinevirtual |
Print iterate status.
This function prints the iteration status.
[in] | algo_state | is the current state of the algorithm |
[in] | printHeader | if set to true will print the header at each iteration |
Reimplemented from ROL::Step< Real >.
Definition at line 585 of file ROL_PrimalDualActiveSetStep.hpp.
References ROL::PrimalDualActiveSetStep< Real >::feasible_, ROL::PrimalDualActiveSetStep< Real >::flag_, ROL::PrimalDualActiveSetStep< Real >::flagCR_, ROL::AlgorithmState< Real >::gnorm, ROL::AlgorithmState< Real >::iter, ROL::PrimalDualActiveSetStep< Real >::iter_, ROL::PrimalDualActiveSetStep< Real >::iterCR_, ROL::PrimalDualActiveSetStep< Real >::maxit_, ROL::AlgorithmState< Real >::nfval, ROL::AlgorithmState< Real >::ngrad, ROL::PrimalDualActiveSetStep< Real >::printHeader(), ROL::PrimalDualActiveSetStep< Real >::printName(), ROL::AlgorithmState< Real >::snorm, and ROL::AlgorithmState< Real >::value.
|
private |
Definition at line 136 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::PrimalDualActiveSetStep().
|
private |
CR iteration counter.
Definition at line 139 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::PrimalDualActiveSetStep< Real >::print(), and ROL::PrimalDualActiveSetStep< Real >::update().
|
private |
CR termination flag.
Definition at line 140 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::PrimalDualActiveSetStep< Real >::print(), and ROL::PrimalDualActiveSetStep< Real >::update().
|
private |
Inexact CR tolerance.
Definition at line 141 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute().
|
private |
Maximum number of PDAS iterations.
Definition at line 144 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::PrimalDualActiveSetStep< Real >::PrimalDualActiveSetStep(), ROL::PrimalDualActiveSetStep< Real >::print(), ROL::PrimalDualActiveSetStep< Real >::printHeader(), and ROL::PrimalDualActiveSetStep< Real >::update().
|
private |
PDAS iteration counter.
Definition at line 145 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::PrimalDualActiveSetStep< Real >::print(), and ROL::PrimalDualActiveSetStep< Real >::update().
|
private |
PDAS termination flag.
Definition at line 146 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::PrimalDualActiveSetStep< Real >::print(), and ROL::PrimalDualActiveSetStep< Real >::update().
|
private |
PDAS minimum step size stopping tolerance.
Definition at line 147 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::PrimalDualActiveSetStep().
|
private |
PDAS gradient stopping tolerance.
Definition at line 148 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::PrimalDualActiveSetStep().
|
private |
Scale for dual variables in the active set, \(c\).
Definition at line 149 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::PrimalDualActiveSetStep().
|
private |
\(\epsilon\)-active set parameter
Definition at line 150 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute().
|
private |
Flag whether the current iterate is feasible or not.
Definition at line 151 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::print(), and ROL::PrimalDualActiveSetStep< Real >::update().
|
private |
Container for dual variables.
Definition at line 154 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::initialize().
|
private |
Container for primal plus dual variables.
Definition at line 155 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::initialize().
|
private |
Container for initial priaml variables.
Definition at line 156 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::initialize().
|
private |
Container for primal variable bounds.
Definition at line 157 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::initialize().
|
private |
Container for step projected onto active set.
Definition at line 158 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::initialize().
|
private |
Container for temporary primal storage.
Definition at line 159 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::PrimalDualActiveSetStep< Real >::computeCriticalityMeasure(), and ROL::PrimalDualActiveSetStep< Real >::initialize().
|
private |
Container for optimality system residual for quadratic model.
Definition at line 160 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::initialize().
|
private |
Container for gradient projected onto active set.
Definition at line 161 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::initialize().
|
private |
Container for temporary right hand side storage.
Definition at line 162 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::initialize().
|
private |
Container for temporary gradient storage.
Definition at line 163 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::PrimalDualActiveSetStep< Real >::initialize(), and ROL::PrimalDualActiveSetStep< Real >::update().
|
private |
Enum for secant type.
Definition at line 166 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::PrimalDualActiveSetStep().
|
private |
Secant object.
Definition at line 167 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::PrimalDualActiveSetStep< Real >::PrimalDualActiveSetStep(), and ROL::PrimalDualActiveSetStep< Real >::update().
|
private |
Definition at line 168 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::PrimalDualActiveSetStep().
|
private |
Definition at line 169 of file ROL_PrimalDualActiveSetStep.hpp.
Referenced by ROL::PrimalDualActiveSetStep< Real >::compute(), and ROL::PrimalDualActiveSetStep< Real >::PrimalDualActiveSetStep().