ROL
Public Member Functions | Private Member Functions | Private Attributes | List of all members
ROL::Reduced_Constraint_SimOpt< Real > Class Template Reference

#include <ROL_Reduced_Constraint_SimOpt.hpp>

+ Inheritance diagram for ROL::Reduced_Constraint_SimOpt< Real >:

Public Member Functions

 Reduced_Constraint_SimOpt (const ROL::Ptr< Constraint_SimOpt< Real > > &conVal, const ROL::Ptr< Constraint_SimOpt< Real > > &conRed, const ROL::Ptr< VectorController< Real > > &stateStore, const ROL::Ptr< Vector< Real > > &state, const ROL::Ptr< Vector< Real > > &control, const ROL::Ptr< Vector< Real > > &adjoint, const ROL::Ptr< Vector< Real > > &residual, bool storage=true, bool useFDhessVec=false)
 Constructor.
 
 Reduced_Constraint_SimOpt (const ROL::Ptr< Constraint_SimOpt< Real > > &conVal, const ROL::Ptr< Constraint_SimOpt< Real > > &conRed, const ROL::Ptr< VectorController< Real > > &stateStore, const ROL::Ptr< Vector< Real > > &state, const ROL::Ptr< Vector< Real > > &control, const ROL::Ptr< Vector< Real > > &adjoint, const ROL::Ptr< Vector< Real > > &residual, const ROL::Ptr< Vector< Real > > &dualstate, const ROL::Ptr< Vector< Real > > &dualcontrol, const ROL::Ptr< Vector< Real > > &dualadjoint, const ROL::Ptr< Vector< Real > > &dualresidual, bool storage=true, bool useFDhessVec=false)
 Secondary, general constructor for use with dual optimization vector spaces where the user does not define the dual() method.
 
void summarize (std::ostream &stream, const Ptr< BatchManager< Real > > &bman=nullPtr) const
 
void reset ()
 
void update (const Vector< Real > &z, bool flag=true, int iter=-1)
 Update the SimOpt objective function and equality constraint.
 
void update (const Vector< Real > &z, UpdateType type, int iter=-1)
 Update constraint function.
 
void value (Vector< Real > &c, const Vector< Real > &z, Real &tol)
 Given \(z\in\mathcal{Z}\), evaluate the equality constraint \(\widehat{c}(z) = c(u(z),z)\) where \(u=u(z)\in\mathcal{U}\) solves \(e(u,z) = 0\).
 
void applyJacobian (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &z, Real &tol)
 Given \(z\in\mathcal{Z}\), apply the Jacobian to a vector \(\widehat{c}'(z)v = c_u(u,z)s + c_z(u,z)v\) where \(s=s(u,z,v)\in\mathcal{U}^*\) solves \(e_u(u,z)s+e_z(u,z)v = 0\).
 
void applyAdjointJacobian (Vector< Real > &ajw, const Vector< Real > &w, const Vector< Real > &z, Real &tol)
 Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\).
 
void applyAdjointHessian (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &z, Real &tol)
 Given \(z\in\mathcal{Z}\), evaluate the Hessian of the objective function \(\nabla^2\widehat{J}(z)\) in the direction \(v\in\mathcal{Z}\).
 
void setParameter (const std::vector< Real > &param)
 
- Public Member Functions inherited from ROL::Constraint< Real >
virtual ~Constraint (void)
 
 Constraint (void)
 
virtual void applyAdjointJacobian (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &dualv, Real &tol)
 Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\).
 
virtual std::vector< Real > solveAugmentedSystem (Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &x, Real &tol)
 Approximately solves the augmented system
 
virtual void applyPreconditioner (Vector< Real > &pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol)
 Apply a constraint preconditioner at \(x\), \(P(x) \in L(\mathcal{C}, \mathcal{C}^*)\), to vector \(v\). Ideally, this preconditioner satisfies the following relationship:
 
void activate (void)
 Turn on constraints.
 
void deactivate (void)
 Turn off constraints.
 
bool isActivated (void)
 Check if constraints are on.
 
virtual std::vector< std::vector< Real > > checkApplyJacobian (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference check for the constraint Jacobian application.
 
virtual std::vector< std::vector< Real > > checkApplyJacobian (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference check for the constraint Jacobian application.
 
virtual std::vector< std::vector< Real > > checkApplyAdjointJacobian (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &c, const Vector< Real > &ajv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS)
 Finite-difference check for the application of the adjoint of constraint Jacobian.
 
virtual Real checkAdjointConsistencyJacobian (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, const bool printToStream=true, std::ostream &outStream=std::cout)
 
virtual Real checkAdjointConsistencyJacobian (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout)
 
virtual std::vector< std::vector< Real > > checkApplyAdjointHessian (const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &step, const bool printToScreen=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference check for the application of the adjoint of constraint Hessian.
 
virtual std::vector< std::vector< Real > > checkApplyAdjointHessian (const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const bool printToScreen=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference check for the application of the adjoint of constraint Hessian.
 

Private Member Functions

void solve_state_equation (const Vector< Real > &z, Real &tol)
 
void solve_adjoint_equation (const Vector< Real > &w, const Vector< Real > &z, Real &tol)
 Given \((u,z)\in\mathcal{U}\times\mathcal{Z}\) which solves the state equation, solve the adjoint equation \(c_u(u,z)^*\lambda + c_u(u,z)^*w = 0\) for \(\lambda=\lambda(u,z)\in\mathcal{C}^*\).
 
void solve_state_sensitivity (const Vector< Real > &v, const Vector< Real > &z, Real &tol)
 Given \((u,z)\in\mathcal{U}\times\mathcal{Z}\) which solves the state equation and a direction \(v\in\mathcal{Z}\), solve the state senstivity equation \(c_u(u,z)s + c_z(u,z)v = 0\) for \(s=u_z(z)v\in\mathcal{U}\).
 
void solve_adjoint_sensitivity (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &z, Real &tol)
 Given \((u,z)\in\mathcal{U}\times\mathcal{Z}\), the adjoint variable \(\lambda\in\mathcal{C}^*\), and a direction \(v\in\mathcal{Z}\), solve the adjoint sensitvity equation \(c_u(u,z)^*p + J_{uu}(u,z)s + J_{uz}(u,z)v + c_{uu}(u,z)(\cdot,s)^*\lambda + c_{zu}(u,z)(\cdot,v)^*\lambda = 0\) for \(p = \lambda_z(u(z),z)v\in\mathcal{C}^*\).
 

Private Attributes

const ROL::Ptr< Constraint_SimOpt< Real > > conVal_
 
const ROL::Ptr< Constraint_SimOpt< Real > > conRed_
 
const ROL::Ptr< VectorController< Real > > stateStore_
 
const ROL::Ptr< VectorController< Real > > adjointStore_
 
const ROL::Ptr< Vector< Real > > state_
 
const ROL::Ptr< Vector< Real > > adjoint_
 
const ROL::Ptr< Vector< Real > > residual_
 
const ROL::Ptr< Vector< Real > > state_sens_
 
const ROL::Ptr< Vector< Real > > adjoint_sens_
 
const ROL::Ptr< Vector< Real > > dualstate_
 
const ROL::Ptr< Vector< Real > > dualstate1_
 
const ROL::Ptr< Vector< Real > > dualadjoint_
 
const ROL::Ptr< Vector< Real > > dualcontrol_
 
const ROL::Ptr< Vector< Real > > dualresidual_
 
const bool storage_
 
const bool useFDhessVec_
 
unsigned nupda_
 
unsigned nvalu_
 
unsigned njaco_
 
unsigned najac_
 
unsigned nhess_
 
unsigned nstat_
 
unsigned nadjo_
 
unsigned nssen_
 
unsigned nasen_
 
bool updateFlag_
 
int updateIter_
 
UpdateType updateType_
 
bool newUpdate_
 
bool isUpdated_
 

Additional Inherited Members

- Protected Member Functions inherited from ROL::Constraint< Real >
const std::vector< Real > getParameter (void) const
 

Detailed Description

template<class Real>
class ROL::Reduced_Constraint_SimOpt< Real >

Definition at line 20 of file ROL_Reduced_Constraint_SimOpt.hpp.

Constructor & Destructor Documentation

◆ Reduced_Constraint_SimOpt() [1/2]

template<class Real >
ROL::Reduced_Constraint_SimOpt< Real >::Reduced_Constraint_SimOpt ( const ROL::Ptr< Constraint_SimOpt< Real > > &  conVal,
const ROL::Ptr< Constraint_SimOpt< Real > > &  conRed,
const ROL::Ptr< VectorController< Real > > &  stateStore,
const ROL::Ptr< Vector< Real > > &  state,
const ROL::Ptr< Vector< Real > > &  control,
const ROL::Ptr< Vector< Real > > &  adjoint,
const ROL::Ptr< Vector< Real > > &  residual,
bool  storage = true,
bool  useFDhessVec = false 
)

Constructor.

Parameters
[in]conValis a pointer to a SimOpt constraint, to be evaluated.
[in]conRedis a pointer to a SimOpt constraint, to be reduced.
[in]stateStoreis a pointer to a VectorController object.
[in]stateis a pointer to a state space vector, \(\mathcal{U}\).
[in]controlis a pointer to a optimization space vector, \(\mathcal{Z}\).
[in]adjointis a pointer to a dual constraint space vector, \(\mathcal{C}_{\text{red}}^*\).
[in]residualis a pointer to a primal constraint space vector, \(\mathcal{C}_{\text{val}}\).
[in]storageis a flag whether or not to store computed states and adjoints.
[in]useFDhessVecis a flag whether or not to use a finite-difference Hessian approximation.

Definition at line 122 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

◆ Reduced_Constraint_SimOpt() [2/2]

template<class Real >
ROL::Reduced_Constraint_SimOpt< Real >::Reduced_Constraint_SimOpt ( const ROL::Ptr< Constraint_SimOpt< Real > > &  conVal,
const ROL::Ptr< Constraint_SimOpt< Real > > &  conRed,
const ROL::Ptr< VectorController< Real > > &  stateStore,
const ROL::Ptr< Vector< Real > > &  state,
const ROL::Ptr< Vector< Real > > &  control,
const ROL::Ptr< Vector< Real > > &  adjoint,
const ROL::Ptr< Vector< Real > > &  residual,
const ROL::Ptr< Vector< Real > > &  dualstate,
const ROL::Ptr< Vector< Real > > &  dualcontrol,
const ROL::Ptr< Vector< Real > > &  dualadjoint,
const ROL::Ptr< Vector< Real > > &  dualresidual,
bool  storage = true,
bool  useFDhessVec = false 
)

Secondary, general constructor for use with dual optimization vector spaces where the user does not define the dual() method.

Parameters
[in]conValis a pointer to a SimOpt constraint, to be evaluated.
[in]conRedis a pointer to a SimOpt constraint, to be reduced.
[in]stateStoreis a pointer to a VectorController object.
[in]stateis a pointer to a state space vector, \(\mathcal{U}\).
[in]controlis a pointer to a optimization space vector, \(\mathcal{Z}\).
[in]adjointis a pointer to a dual constraint space vector, \(\mathcal{C}_{\text{red}}^*\).
[in]residualis a pointer to a primal constraint space vector, \(\mathcal{C}_{\text{val}}\).
[in]dualstateis a pointer to a dual state space vector, \(\mathcal{U}^*\).
[in]dualadjointis a pointer to a constraint space vector, \(\mathcal{C}_{\text{red}}\).
[in]dualresidualis a pointer to a dual constraint space vector, \(\mathcal{C}_{\text{val}}^*\).
[in]storageis a flag whether or not to store computed states and adjoints.
[in]useFDhessVecis a flag whether or not to use a finite-difference Hessian approximation.

Definition at line 153 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

Member Function Documentation

◆ solve_state_equation()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::solve_state_equation ( const Vector< Real > &  z,
Real &  tol 
)
private

Definition at line 51 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

◆ solve_adjoint_equation()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::solve_adjoint_equation ( const Vector< Real > &  w,
const Vector< Real > &  z,
Real &  tol 
)
private

Given \((u,z)\in\mathcal{U}\times\mathcal{Z}\) which solves the state equation, solve the adjoint equation \(c_u(u,z)^*\lambda + c_u(u,z)^*w = 0\) for \(\lambda=\lambda(u,z)\in\mathcal{C}^*\).

Definition at line 79 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

◆ solve_state_sensitivity()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::solve_state_sensitivity ( const Vector< Real > &  v,
const Vector< Real > &  z,
Real &  tol 
)
private

Given \((u,z)\in\mathcal{U}\times\mathcal{Z}\) which solves the state equation and a direction \(v\in\mathcal{Z}\), solve the state senstivity equation \(c_u(u,z)s + c_z(u,z)v = 0\) for \(s=u_z(z)v\in\mathcal{U}\).

Definition at line 96 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

◆ solve_adjoint_sensitivity()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::solve_adjoint_sensitivity ( const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  z,
Real &  tol 
)
private

Given \((u,z)\in\mathcal{U}\times\mathcal{Z}\), the adjoint variable \(\lambda\in\mathcal{C}^*\), and a direction \(v\in\mathcal{Z}\), solve the adjoint sensitvity equation \(c_u(u,z)^*p + J_{uu}(u,z)s + J_{uz}(u,z)v + c_{uu}(u,z)(\cdot,s)^*\lambda + c_{zu}(u,z)(\cdot,v)^*\lambda = 0\) for \(p = \lambda_z(u(z),z)v\in\mathcal{C}^*\).

Definition at line 105 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

◆ summarize()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::summarize ( std::ostream &  stream,
const Ptr< BatchManager< Real > > &  bman = nullPtr 
) const

Definition at line 188 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

◆ reset()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::reset ( )

Definition at line 234 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

◆ update() [1/2]

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::update ( const Vector< Real > &  z,
bool  flag = true,
int  iter = -1 
)
virtual

Update the SimOpt objective function and equality constraint.

Reimplemented from ROL::Constraint< Real >.

Definition at line 240 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

◆ update() [2/2]

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::update ( const Vector< Real > &  x,
UpdateType  type,
int  iter = -1 
)
virtual

Update constraint function.

This function updates the constraint function at new iterations.

Parameters
[in]xis the new iterate.
[in]typeis the type of update requested.
[in]iteris the outer algorithm iterations count.

Reimplemented from ROL::Constraint< Real >.

Definition at line 248 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

◆ value()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::value ( Vector< Real > &  c,
const Vector< Real > &  z,
Real &  tol 
)
virtual

Given \(z\in\mathcal{Z}\), evaluate the equality constraint \(\widehat{c}(z) = c(u(z),z)\) where \(u=u(z)\in\mathcal{U}\) solves \(e(u,z) = 0\).

Implements ROL::Constraint< Real >.

Definition at line 259 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

◆ applyJacobian()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::applyJacobian ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  z,
Real &  tol 
)
virtual

Given \(z\in\mathcal{Z}\), apply the Jacobian to a vector \(\widehat{c}'(z)v = c_u(u,z)s + c_z(u,z)v\) where \(s=s(u,z,v)\in\mathcal{U}^*\) solves \(e_u(u,z)s+e_z(u,z)v = 0\).

Reimplemented from ROL::Constraint< Real >.

Definition at line 268 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

References ROL::Vector< Real >::plus().

◆ applyAdjointJacobian()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::applyAdjointJacobian ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  x,
Real &  tol 
)
virtual

Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\).

Parameters
[out]ajvis the result of applying the adjoint of the constraint Jacobian to v at x; a dual optimization-space vector
[in]vis a dual constraint-space vector
[in]xis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{ajv} = c'(x)^*v\), where \(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{X}^*\).

The default implementation is a finite-difference approximation.


Reimplemented from ROL::Constraint< Real >.

Definition at line 283 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

References ROL::Vector< Real >::plus().

◆ applyAdjointHessian()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::applyAdjointHessian ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  z,
Real &  tol 
)
virtual

Given \(z\in\mathcal{Z}\), evaluate the Hessian of the objective function \(\nabla^2\widehat{J}(z)\) in the direction \(v\in\mathcal{Z}\).

Reimplemented from ROL::Constraint< Real >.

Definition at line 298 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.

References ROL::Constraint< Real >::applyAdjointHessian(), and ROL::Vector< Real >::plus().

◆ setParameter()

template<class Real >
void ROL::Reduced_Constraint_SimOpt< Real >::setParameter ( const std::vector< Real > &  param)
inlinevirtual

Member Data Documentation

◆ conVal_

template<class Real >
const ROL::Ptr<Constraint_SimOpt<Real> > ROL::Reduced_Constraint_SimOpt< Real >::conVal_
private

◆ conRed_

template<class Real >
const ROL::Ptr<Constraint_SimOpt<Real> > ROL::Reduced_Constraint_SimOpt< Real >::conRed_
private

◆ stateStore_

template<class Real >
const ROL::Ptr<VectorController<Real> > ROL::Reduced_Constraint_SimOpt< Real >::stateStore_
private

Definition at line 23 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ adjointStore_

template<class Real >
const ROL::Ptr<VectorController<Real> > ROL::Reduced_Constraint_SimOpt< Real >::adjointStore_
private

Definition at line 23 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ state_

template<class Real >
const ROL::Ptr<Vector<Real> > ROL::Reduced_Constraint_SimOpt< Real >::state_
private

Definition at line 26 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ adjoint_

template<class Real >
const ROL::Ptr<Vector<Real> > ROL::Reduced_Constraint_SimOpt< Real >::adjoint_
private

Definition at line 26 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ residual_

template<class Real >
const ROL::Ptr<Vector<Real> > ROL::Reduced_Constraint_SimOpt< Real >::residual_
private

Definition at line 26 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ state_sens_

template<class Real >
const ROL::Ptr<Vector<Real> > ROL::Reduced_Constraint_SimOpt< Real >::state_sens_
private

Definition at line 27 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ adjoint_sens_

template<class Real >
const ROL::Ptr<Vector<Real> > ROL::Reduced_Constraint_SimOpt< Real >::adjoint_sens_
private

Definition at line 27 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ dualstate_

template<class Real >
const ROL::Ptr<Vector<Real> > ROL::Reduced_Constraint_SimOpt< Real >::dualstate_
private

Definition at line 30 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ dualstate1_

template<class Real >
const ROL::Ptr<Vector<Real> > ROL::Reduced_Constraint_SimOpt< Real >::dualstate1_
private

Definition at line 30 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ dualadjoint_

template<class Real >
const ROL::Ptr<Vector<Real> > ROL::Reduced_Constraint_SimOpt< Real >::dualadjoint_
private

Definition at line 30 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ dualcontrol_

template<class Real >
const ROL::Ptr<Vector<Real> > ROL::Reduced_Constraint_SimOpt< Real >::dualcontrol_
private

Definition at line 31 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ dualresidual_

template<class Real >
const ROL::Ptr<Vector<Real> > ROL::Reduced_Constraint_SimOpt< Real >::dualresidual_
private

Definition at line 31 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ storage_

template<class Real >
const bool ROL::Reduced_Constraint_SimOpt< Real >::storage_
private

Definition at line 33 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ useFDhessVec_

template<class Real >
const bool ROL::Reduced_Constraint_SimOpt< Real >::useFDhessVec_
private

Definition at line 34 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ nupda_

template<class Real >
unsigned ROL::Reduced_Constraint_SimOpt< Real >::nupda_
private

Definition at line 36 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ nvalu_

template<class Real >
unsigned ROL::Reduced_Constraint_SimOpt< Real >::nvalu_
private

Definition at line 36 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ njaco_

template<class Real >
unsigned ROL::Reduced_Constraint_SimOpt< Real >::njaco_
private

Definition at line 36 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ najac_

template<class Real >
unsigned ROL::Reduced_Constraint_SimOpt< Real >::najac_
private

Definition at line 36 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ nhess_

template<class Real >
unsigned ROL::Reduced_Constraint_SimOpt< Real >::nhess_
private

Definition at line 36 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ nstat_

template<class Real >
unsigned ROL::Reduced_Constraint_SimOpt< Real >::nstat_
private

Definition at line 37 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ nadjo_

template<class Real >
unsigned ROL::Reduced_Constraint_SimOpt< Real >::nadjo_
private

Definition at line 37 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ nssen_

template<class Real >
unsigned ROL::Reduced_Constraint_SimOpt< Real >::nssen_
private

Definition at line 37 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ nasen_

template<class Real >
unsigned ROL::Reduced_Constraint_SimOpt< Real >::nasen_
private

Definition at line 37 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ updateFlag_

template<class Real >
bool ROL::Reduced_Constraint_SimOpt< Real >::updateFlag_
private

Definition at line 39 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ updateIter_

template<class Real >
int ROL::Reduced_Constraint_SimOpt< Real >::updateIter_
private

Definition at line 40 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ updateType_

template<class Real >
UpdateType ROL::Reduced_Constraint_SimOpt< Real >::updateType_
private

Definition at line 41 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ newUpdate_

template<class Real >
bool ROL::Reduced_Constraint_SimOpt< Real >::newUpdate_
private

Definition at line 42 of file ROL_Reduced_Constraint_SimOpt.hpp.

◆ isUpdated_

template<class Real >
bool ROL::Reduced_Constraint_SimOpt< Real >::isUpdated_
private

Definition at line 43 of file ROL_Reduced_Constraint_SimOpt.hpp.


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