ROL
|
#include <ROL_Reduced_Constraint_SimOpt.hpp>
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 > ¶m) |
![]() | |
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 | |
![]() | |
const std::vector< Real > | getParameter (void) const |
Definition at line 20 of file ROL_Reduced_Constraint_SimOpt.hpp.
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.
[in] | conVal | is a pointer to a SimOpt constraint, to be evaluated. |
[in] | conRed | is a pointer to a SimOpt constraint, to be reduced. |
[in] | stateStore | is a pointer to a VectorController object. |
[in] | state | is a pointer to a state space vector, \(\mathcal{U}\). |
[in] | control | is a pointer to a optimization space vector, \(\mathcal{Z}\). |
[in] | adjoint | is a pointer to a dual constraint space vector, \(\mathcal{C}_{\text{red}}^*\). |
[in] | residual | is a pointer to a primal constraint space vector, \(\mathcal{C}_{\text{val}}\). |
[in] | storage | is a flag whether or not to store computed states and adjoints. |
[in] | useFDhessVec | is a flag whether or not to use a finite-difference Hessian approximation. |
Definition at line 122 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.
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.
[in] | conVal | is a pointer to a SimOpt constraint, to be evaluated. |
[in] | conRed | is a pointer to a SimOpt constraint, to be reduced. |
[in] | stateStore | is a pointer to a VectorController object. |
[in] | state | is a pointer to a state space vector, \(\mathcal{U}\). |
[in] | control | is a pointer to a optimization space vector, \(\mathcal{Z}\). |
[in] | adjoint | is a pointer to a dual constraint space vector, \(\mathcal{C}_{\text{red}}^*\). |
[in] | residual | is a pointer to a primal constraint space vector, \(\mathcal{C}_{\text{val}}\). |
[in] | dualstate | is a pointer to a dual state space vector, \(\mathcal{U}^*\). |
[in] | dualadjoint | is a pointer to a constraint space vector, \(\mathcal{C}_{\text{red}}\). |
[in] | dualresidual | is a pointer to a dual constraint space vector, \(\mathcal{C}_{\text{val}}^*\). |
[in] | storage | is a flag whether or not to store computed states and adjoints. |
[in] | useFDhessVec | is a flag whether or not to use a finite-difference Hessian approximation. |
Definition at line 153 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.
|
private |
Definition at line 51 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.
|
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.
|
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.
|
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.
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.
void ROL::Reduced_Constraint_SimOpt< Real >::reset | ( | ) |
Definition at line 234 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.
|
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.
|
virtual |
Update constraint function.
This function updates the constraint function at new iterations.
[in] | x | is the new iterate. |
[in] | type | is the type of update requested. |
[in] | iter | is the outer algorithm iterations count. |
Reimplemented from ROL::Constraint< Real >.
Definition at line 248 of file ROL_Reduced_Constraint_SimOpt_Def.hpp.
|
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.
|
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().
|
virtual |
Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\).
[out] | ajv | is the result of applying the adjoint of the constraint Jacobian to v at x; a dual optimization-space vector |
[in] | v | is a dual constraint-space vector |
[in] | x | is the constraint argument; an optimization-space vector |
[in,out] | tol | is 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().
|
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().
|
inlinevirtual |
Reimplemented from ROL::Constraint< Real >.
Definition at line 158 of file ROL_Reduced_Constraint_SimOpt.hpp.
References ROL::Reduced_Constraint_SimOpt< Real >::conRed_, ROL::Reduced_Constraint_SimOpt< Real >::conVal_, and ROL::Constraint< Real >::setParameter().
|
private |
Definition at line 22 of file ROL_Reduced_Constraint_SimOpt.hpp.
Referenced by ROL::Reduced_Constraint_SimOpt< Real >::setParameter().
|
private |
Definition at line 22 of file ROL_Reduced_Constraint_SimOpt.hpp.
Referenced by ROL::Reduced_Constraint_SimOpt< Real >::setParameter().
|
private |
Definition at line 23 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 23 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 26 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 26 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 26 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 27 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 27 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 30 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 30 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 30 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 31 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 31 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 33 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 34 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 36 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 36 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 36 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 36 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 36 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 37 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 37 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 37 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 37 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 39 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 40 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 41 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 42 of file ROL_Reduced_Constraint_SimOpt.hpp.
|
private |
Definition at line 43 of file ROL_Reduced_Constraint_SimOpt.hpp.