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

#include <ROL_Constraint_PinTSimOpt.hpp>

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

Public Member Functions

 Constraint_PinTSimOpt ()
 Default constructor. More...
 
 Constraint_PinTSimOpt (const Ptr< Constraint_TimeSimOpt< Real >> &stepConstraint)
 Constructor. More...
 
void initialize (const Ptr< Constraint_TimeSimOpt< Real >> &stepConstraint)
 Initialize this class, setting up parallel distribution. More...
 
virtual void update_1 (const Vector< Real > &u, bool flag=true, int iter=-1)
 Update constraint functions with respect to Opt variable. u is the state variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void update_2 (const Vector< Real > &z, bool flag=true, int iter=-1)
 Update constraint functions with respect to Opt variable. z is the control variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void value (Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Evaluate the constraint operator \(c:\mathcal{U}\times\mathcal{Z} \rightarrow \mathcal{C}\) at \((u,z)\). More...
 
virtual void solve (Vector< Real > &c, Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Given \(z\), solve \(c(u,z)=0\) for \(u\). More...
 
virtual void applyJacobian_1 (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the partial constraint Jacobian at \((u,z)\), \(c_u(u,z) \in L(\mathcal{U}, \mathcal{C})\), to the vector \(v\). More...
 
virtual void applyJacobian_2 (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the partial constraint Jacobian at \((u,z)\), \(c_z(u,z) \in L(\mathcal{Z}, \mathcal{C})\), to the vector \(v\). More...
 
virtual void applyInverseJacobian_1 (Vector< Real > &ijv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) overridefinal
 Apply the inverse partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^{-1} \in L(\mathcal{C}, \mathcal{U})\), to the vector \(v\). More...
 
virtual void applyAdjointJacobian_1 (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^* \in L(\mathcal{C}^*, \mathcal{U}^*)\), to the vector \(v\). This is the primary interface. More...
 
virtual void applyAdjointJacobian_2 (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_z(u,z)^* \in L(\mathcal{C}^*, \mathcal{Z}^*)\), to vector \(v\). This is the primary interface. More...
 
virtual void applyInverseAdjointJacobian_1 (Vector< Real > &iajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) overridefinal
 Apply the inverse of the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^{-*} \in L(\mathcal{U}^*, \mathcal{C}^*)\), to the vector \(v\). More...
 
virtual void applyAdjointHessian_11 (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{uu}(u,z)(v,\cdot)^*w\). More...
 
virtual void applyAdjointHessian_12 (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{uz}(u,z)(v,\cdot)^*w\). More...
 
virtual void applyAdjointHessian_21 (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{zu}(u,z)(v,\cdot)^*w\). More...
 
virtual void applyAdjointHessian_22 (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{zz}(u,z)(v,\cdot)^*w\). More...
 
- Public Member Functions inherited from ROL::ROL::Constraint_SimOpt< Real >
 Constraint_SimOpt ()
 
virtual void update (const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1)
 Update constraint functions. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void setSolveParameters (Teuchos::ParameterList &parlist)
 Set solve parameters. More...
 
virtual void applyAdjointJacobian_1 (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualv, Real &tol)
 Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^* \in L(\mathcal{C}^*, \mathcal{U}^*)\), to the vector \(v\). This is the secondary interface, for use with dual spaces where the user does not define the dual() operation. More...
 
virtual void applyAdjointJacobian_2 (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualv, Real &tol)
 Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_z(u,z)^* \in L(\mathcal{C}^*, \mathcal{Z}^*)\), to vector \(v\). This is the secondary interface, for use with dual spaces where the user does not define the dual() operation. More...
 
virtual std::vector< RealsolveAugmentedSystem (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

\[ \begin{pmatrix} I & c'(x)^* \\ c'(x) & 0 \end{pmatrix} \begin{pmatrix} v_{1} \\ v_{2} \end{pmatrix} = \begin{pmatrix} b_{1} \\ b_{2} \end{pmatrix} \]

where \(v_{1} \in \mathcal{X}\), \(v_{2} \in \mathcal{C}^*\), \(b_{1} \in \mathcal{X}^*\), \(b_{2} \in \mathcal{C}\), \(I : \mathcal{X} \rightarrow \mathcal{X}^*\) is an identity operator, and \(0 : \mathcal{C}^* \rightarrow \mathcal{C}\) is a zero operator. More...

 
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\). In general, this preconditioner satisfies the following relationship:

\[ c'(x) c'(x)^* P(x) v \approx v \,. \]

It is used by the solveAugmentedSystem method. More...

 
virtual void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update constraint functions. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void value (Vector< Real > &c, const Vector< Real > &x, Real &tol)
 Evaluate the constraint operator \(c:\mathcal{X} \rightarrow \mathcal{C}\) at \(x\). More...
 
virtual void applyJacobian (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply the constraint Jacobian at \(x\), \(c'(x) \in L(\mathcal{X}, \mathcal{C})\), to vector \(v\). More...
 
virtual void applyAdjointJacobian (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\). More...
 
virtual void applyAdjointHessian (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply the derivative of the adjoint of the constraint Jacobian at \(x\) to vector \(u\) in direction \(v\), according to \( v \mapsto c''(x)(v,\cdot)^*u \). More...
 
virtual Real checkSolve (const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, const ROL::Vector< Real > &c, const bool printToStream=true, std::ostream &outStream=std::cout)
 
virtual Real checkAdjointConsistencyJacobian_1 (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
 Check the consistency of the Jacobian and its adjoint. This is the primary interface. More...
 
virtual Real checkAdjointConsistencyJacobian_1 (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout)
 Check the consistency of the Jacobian and its adjoint. This is the secondary interface, for use with dual spaces where the user does not define the dual() operation. More...
 
virtual Real checkAdjointConsistencyJacobian_2 (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
 Check the consistency of the Jacobian and its adjoint. This is the primary interface. More...
 
virtual Real checkAdjointConsistencyJacobian_2 (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout)
 Check the consistency of the Jacobian and its adjoint. This is the secondary interface, for use with dual spaces where the user does not define the dual() operation. More...
 
virtual Real checkInverseJacobian_1 (const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
 
virtual Real checkInverseAdjointJacobian_1 (const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
 
std::vector< std::vector< Real > > checkApplyJacobian_1 (const Vector< Real > &u, const Vector< Real > &z, 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)
 
std::vector< std::vector< Real > > checkApplyJacobian_1 (const Vector< Real > &u, const Vector< Real > &z, 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)
 
std::vector< std::vector< Real > > checkApplyJacobian_2 (const Vector< Real > &u, const Vector< Real > &z, 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)
 
std::vector< std::vector< Real > > checkApplyJacobian_2 (const Vector< Real > &u, const Vector< Real > &z, 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)
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_11 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_11 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_21 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 \( u\in U \), \( z\in Z \), \( p\in C^\ast \), \( v \in U \), \( hv \in U^\ast \) More...
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_21 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 \( u\in U \), \( z\in Z \), \( p\in C^\ast \), \( v \in U \), \( hv \in U^\ast \) More...
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_12 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 \( u\in U \), \( z\in Z \), \( p\in C^\ast \), \( v \in U \), \( hv \in U^\ast \) More...
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_12 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_22 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_22 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
- 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\). More...
 
void activate (void)
 Turn on constraints. More...
 
void deactivate (void)
 Turn off constraints. More...
 
bool isActivated (void)
 Check if constraints are on. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
virtual void setParameter (const std::vector< Real > &param)
 

Private Types

enum  ETimeAccessor { PAST, CURRENT, FUTURE, ALL }
 Enumeration to define which components of a a vector are required. More...
 

Private Member Functions

Ptr< const Vector< Real > > getVector (const PinTVector< double > &src, int step, ETimeAccessor accessor)
 Build a vector composed of all vectors required by accessor. More...
 
Ptr< Vector< Real > > getNonconstVector (const PinTVector< double > &src, int step, ETimeAccessor accessor)
 Build a vector composed of all vectors required by accessor. More...
 
Ptr< Vector< Real > > getStateVector (const PinTVector< double > &src, int step)
 

Private Attributes

Ptr< Constraint_TimeSimOpt
< Real > > 
stepConstraint_
 
bool isInitialized_
 

Additional Inherited Members

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

Detailed Description

template<class Real>
class ROL::Constraint_PinTSimOpt< Real >

Definition at line 121 of file ROL_Constraint_PinTSimOpt.hpp.

Member Enumeration Documentation

template<class Real >
enum ROL::Constraint_PinTSimOpt::ETimeAccessor
private

Enumeration to define which components of a a vector are required.

Enumerator
PAST 
CURRENT 
FUTURE 
ALL 

Definition at line 130 of file ROL_Constraint_PinTSimOpt.hpp.

Constructor & Destructor Documentation

template<class Real >
ROL::Constraint_PinTSimOpt< Real >::Constraint_PinTSimOpt ( )
inline

Default constructor.

Definition at line 199 of file ROL_Constraint_PinTSimOpt.hpp.

template<class Real >
ROL::Constraint_PinTSimOpt< Real >::Constraint_PinTSimOpt ( const Ptr< Constraint_TimeSimOpt< Real >> &  stepConstraint)
inline

Constructor.

Build a parallel-in-time constraint with a specified step constraint. This specifies any communication you might need between processors and steps using "stencils".

Parameters
[in]stepConstraintConstraint for a single step.
[in]stepsThe number of steps to take.
[in]stepStateState vector step value, also the initial condition
[in]stateStencilStencil for parallel stepping.
[in]controlStateControl vector step value, also the initial condition
[in]controlStencilStencil for control parallel stepping.

Definition at line 217 of file ROL_Constraint_PinTSimOpt.hpp.

References ROL::Constraint_PinTSimOpt< Real >::initialize().

Member Function Documentation

template<class Real >
Ptr<const Vector<Real> > ROL::Constraint_PinTSimOpt< Real >::getVector ( const PinTVector< double > &  src,
int  step,
ETimeAccessor  accessor 
)
inlineprivate

Build a vector composed of all vectors required by accessor.

Build a vector composed of all vectors require by accessor. If this is only one vector, that vector is returned directly. If there is more than one vector than a partioned vector is returned with each entry.

Definition at line 138 of file ROL_Constraint_PinTSimOpt.hpp.

References ROL::Constraint_PinTSimOpt< Real >::getNonconstVector().

Referenced by ROL::Constraint_PinTSimOpt< Real >::applyAdjointHessian_11(), ROL::Constraint_PinTSimOpt< Real >::applyAdjointJacobian_1(), ROL::Constraint_PinTSimOpt< Real >::applyAdjointJacobian_2(), ROL::Constraint_PinTSimOpt< Real >::applyJacobian_1(), ROL::Constraint_PinTSimOpt< Real >::applyJacobian_2(), ROL::Constraint_PinTSimOpt< Real >::solve(), and ROL::Constraint_PinTSimOpt< Real >::value().

template<class Real >
Ptr<Vector<Real> > ROL::Constraint_PinTSimOpt< Real >::getNonconstVector ( const PinTVector< double > &  src,
int  step,
ETimeAccessor  accessor 
)
inlineprivate
template<class Real >
Ptr<Vector<Real> > ROL::Constraint_PinTSimOpt< Real >::getStateVector ( const PinTVector< double > &  src,
int  step 
)
inlineprivate
template<class Real >
void ROL::Constraint_PinTSimOpt< Real >::initialize ( const Ptr< Constraint_TimeSimOpt< Real >> &  stepConstraint)
inline
template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::update_1 ( const Vector< Real > &  u,
bool  flag = true,
int  iter = -1 
)
inlinevirtual

Update constraint functions with respect to Opt variable. u is the state variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.

Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 240 of file ROL_Constraint_PinTSimOpt.hpp.

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

Update constraint functions with respect to Opt variable. z is the control variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.

Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 247 of file ROL_Constraint_PinTSimOpt.hpp.

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

Evaluate the constraint operator \(c:\mathcal{U}\times\mathcal{Z} \rightarrow \mathcal{C}\) at \((u,z)\).

Parameters
[out]cis the result of evaluating the constraint operator at \((u,z)\); a constraint-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{c} = c(u,z)\), where \(\mathsf{c} \in \mathcal{C}\), \(\mathsf{u} \in \mathcal{U}\), and $ \(\mathsf{z} \in\mathcal{Z}\).


Implements ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 249 of file ROL_Constraint_PinTSimOpt.hpp.

References ROL::Constraint_PinTSimOpt< Real >::ALL, ROL::PinTVector< Real >::boundaryExchange(), ROL::Constraint_PinTSimOpt< Real >::CURRENT, ROL::Constraint_PinTSimOpt< Real >::getNonconstVector(), ROL::PinTVector< Real >::getRemoteBufferPtr(), ROL::Constraint_PinTSimOpt< Real >::getStateVector(), ROL::Constraint_PinTSimOpt< Real >::getVector(), ROL::PinTVector< Real >::getVectorPtr(), ROL::PinTVector< Real >::numOwnedSteps(), ROL::PinTVector< Real >::stencil(), and ROL::Constraint_PinTSimOpt< Real >::stepConstraint_.

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::solve ( Vector< Real > &  c,
Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlineoverridevirtual

Given \(z\), solve \(c(u,z)=0\) for \(u\).

Parameters
[out]cis the result of evaluating the constraint operator at \((u,z)\); a constraint-space vector
[in,out]uis the solution vector; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

The defualt implementation is Newton's method globalized with a backtracking line search.


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 286 of file ROL_Constraint_PinTSimOpt.hpp.

References ROL::Constraint_PinTSimOpt< Real >::ALL, ROL::PinTVector< Real >::boundaryExchange(), ROL::Constraint_PinTSimOpt< Real >::CURRENT, ROL::Constraint_PinTSimOpt< Real >::getNonconstVector(), ROL::PinTVector< Real >::getRemoteBufferPtr(), ROL::Constraint_PinTSimOpt< Real >::getVector(), ROL::PinTVector< Real >::getVectorPtr(), ROL::PinTVector< Real >::numOwnedSteps(), ROL::Constraint_PinTSimOpt< Real >::PAST, ROL::PinTVector< Real >::stencil(), and ROL::Constraint_PinTSimOpt< Real >::stepConstraint_.

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::applyJacobian_1 ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlineoverridevirtual

Apply the partial constraint Jacobian at \((u,z)\), \(c_u(u,z) \in L(\mathcal{U}, \mathcal{C})\), to the vector \(v\).

Parameters
[out]jvis the result of applying the constraint Jacobian to v at \((u,z)\); a constraint-space vector
[in]vis a simulation-space vector
[in]uis the constraint argument; an simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{jv} = c_u(u,z)v\), where \(v \in \mathcal{U}\), \(\mathsf{jv} \in \mathcal{C}\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 334 of file ROL_Constraint_PinTSimOpt.hpp.

References ROL::Constraint_PinTSimOpt< Real >::ALL, ROL::PinTVector< Real >::boundaryExchange(), ROL::PinTVector< Real >::communicators(), ROL::Constraint_PinTSimOpt< Real >::CURRENT, ROL::Constraint_PinTSimOpt< Real >::getNonconstVector(), ROL::PinTVector< Real >::getRemoteBufferPtr(), ROL::Constraint_PinTSimOpt< Real >::getStateVector(), ROL::Constraint_PinTSimOpt< Real >::getVector(), ROL::PinTVector< Real >::getVectorPtr(), ROL::PinTVector< Real >::numOwnedSteps(), ROL::PinTVector< Real >::stencil(), and ROL::Constraint_PinTSimOpt< Real >::stepConstraint_.

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::applyJacobian_2 ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlineoverridevirtual

Apply the partial constraint Jacobian at \((u,z)\), \(c_z(u,z) \in L(\mathcal{Z}, \mathcal{C})\), to the vector \(v\).

Parameters
[out]jvis the result of applying the constraint Jacobian to v at \((u,z)\); a constraint-space vector
[in]vis an optimization-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{jv} = c_z(u,z)v\), where \(v \in \mathcal{Z}\), \(\mathsf{jv} \in \mathcal{C}\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 380 of file ROL_Constraint_PinTSimOpt.hpp.

References ROL::Constraint_PinTSimOpt< Real >::ALL, ROL::PinTVector< Real >::boundaryExchange(), ROL::Constraint_PinTSimOpt< Real >::CURRENT, ROL::Constraint_PinTSimOpt< Real >::getNonconstVector(), ROL::Constraint_PinTSimOpt< Real >::getStateVector(), ROL::Constraint_PinTSimOpt< Real >::getVector(), ROL::PinTVector< Real >::getVectorPtr(), ROL::PinTVector< Real >::numOwnedSteps(), ROL::PinTVector< Real >::stencil(), and ROL::Constraint_PinTSimOpt< Real >::stepConstraint_.

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::applyInverseJacobian_1 ( Vector< Real > &  ijv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlinefinaloverridevirtual

Apply the inverse partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^{-1} \in L(\mathcal{C}, \mathcal{U})\), to the vector \(v\).

Parameters
[out]ijvis the result of applying the inverse constraint Jacobian to v at \((u,z)\); a simulation-space vector
[in]vis a constraint-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{ijv} = c_u(u,z)^{-1}v\), where \(v \in \mathcal{C}\), \(\mathsf{ijv} \in \mathcal{U}\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 419 of file ROL_Constraint_PinTSimOpt.hpp.

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::applyAdjointJacobian_1 ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlineoverridevirtual

Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^* \in L(\mathcal{C}^*, \mathcal{U}^*)\), to the vector \(v\). This is the primary interface.

Parameters
[out]ajvis the result of applying the adjoint of the constraint Jacobian to v at (u,z); a dual simulation-space vector
[in]vis a dual constraint-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{ajv} = c_u(u,z)^*v\), where \(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{U}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 428 of file ROL_Constraint_PinTSimOpt.hpp.

References ROL::Constraint_PinTSimOpt< Real >::ALL, ROL::PinTVector< Real >::boundaryExchange(), ROL::PinTVector< Real >::boundaryExchangeSumInto(), ROL::Constraint_PinTSimOpt< Real >::CURRENT, ROL::PinTVector< Real >::getRemoteBufferPtr(), ROL::Constraint_PinTSimOpt< Real >::getStateVector(), ROL::Constraint_PinTSimOpt< Real >::getVector(), ROL::PinTVector< Real >::getVectorPtr(), ROL::PinTVector< Real >::numOwnedSteps(), ROL::PinTVector< Real >::stencil(), ROL::Constraint_PinTSimOpt< Real >::stepConstraint_, and ROL::PinTVector< Real >::zeroAll().

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::applyAdjointJacobian_2 ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlineoverridevirtual

Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_z(u,z)^* \in L(\mathcal{C}^*, \mathcal{Z}^*)\), to vector \(v\). This is the primary interface.

Parameters
[out]ajvis the result of applying the adjoint of the constraint Jacobian to v at \((u,z)\); a dual optimization-space vector
[in]vis a dual constraint-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{ajv} = c_z(u,z)^*v\), where \(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{Z}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 485 of file ROL_Constraint_PinTSimOpt.hpp.

References ROL::Constraint_PinTSimOpt< Real >::ALL, ROL::PinTVector< Real >::boundaryExchange(), ROL::Constraint_PinTSimOpt< Real >::CURRENT, ROL::Constraint_PinTSimOpt< Real >::getNonconstVector(), ROL::Constraint_PinTSimOpt< Real >::getStateVector(), ROL::Constraint_PinTSimOpt< Real >::getVector(), ROL::PinTVector< Real >::numOwnedSteps(), ROL::Constraint_PinTSimOpt< Real >::stepConstraint_, and ROL::PinTVector< Real >::zeroAll().

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::applyInverseAdjointJacobian_1 ( Vector< Real > &  iajv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlinefinaloverridevirtual

Apply the inverse of the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^{-*} \in L(\mathcal{U}^*, \mathcal{C}^*)\), to the vector \(v\).

Parameters
[out]iajvis the result of applying the inverse adjoint of the constraint Jacobian to v at (u,z); a dual constraint-space vector
[in]vis a dual simulation-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{iajv} = c_u(u,z)^{-*}v\), where \(v \in \mathcal{U}^*\), \(\mathsf{iajv} \in \mathcal{C}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 523 of file ROL_Constraint_PinTSimOpt.hpp.

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::applyAdjointHessian_11 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlineoverridevirtual

Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{uu}(u,z)(v,\cdot)^*w\).

Parameters
[out]ahwvis the result of applying the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at \((u,z)\) to the vector \(w\) in direction \(w\); a dual simulation-space vector
[in]wis the direction vector; a dual constraint-space vector
[in]vis a simulation-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{ahwv} = c_{uu}(u,z)(v,\cdot)^*w\), where \(w \in \mathcal{C}^*\), \(v \in \mathcal{U}\), and \(\mathsf{ahwv} \in \mathcal{U}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 549 of file ROL_Constraint_PinTSimOpt.hpp.

References ROL::Constraint_PinTSimOpt< Real >::ALL, ROL::PinTVector< Real >::boundaryExchange(), ROL::Constraint_PinTSimOpt< Real >::CURRENT, ROL::Constraint_PinTSimOpt< Real >::getStateVector(), ROL::Constraint_PinTSimOpt< Real >::getVector(), ROL::PinTVector< Real >::numOwnedSteps(), ROL::Constraint_PinTSimOpt< Real >::stepConstraint_, and ROL::PinTVector< Real >::zeroAll().

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::applyAdjointHessian_12 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlineoverridevirtual

Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{uz}(u,z)(v,\cdot)^*w\).

Parameters
[out]ahwvis the result of applying the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at \((u,z)\) to the vector \(w\) in direction \(w\); a dual optimization-space vector
[in]wis the direction vector; a dual constraint-space vector
[in]vis a simulation-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{ahwv} = c_{uz}(u,z)(v,\cdot)^*w\), where \(w \in \mathcal{C}^*\), \(v \in \mathcal{U}\), and \(\mathsf{ahwv} \in \mathcal{Z}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 616 of file ROL_Constraint_PinTSimOpt.hpp.

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

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::applyAdjointHessian_21 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlineoverridevirtual

Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{zu}(u,z)(v,\cdot)^*w\).

Parameters
[out]ahwvis the result of applying the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at \((u,z)\) to the vector \(w\) in direction \(w\); a dual simulation-space vector
[in]wis the direction vector; a dual constraint-space vector
[in]vis a optimization-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{ahwv} = c_{zu}(u,z)(v,\cdot)^*w\), where \(w \in \mathcal{C}^*\), \(v \in \mathcal{Z}\), and \(\mathsf{ahwv} \in \mathcal{U}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 643 of file ROL_Constraint_PinTSimOpt.hpp.

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

template<class Real >
virtual void ROL::Constraint_PinTSimOpt< Real >::applyAdjointHessian_22 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real tol 
)
inlineoverridevirtual

Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{zz}(u,z)(v,\cdot)^*w\).

Parameters
[out]ahwvis the result of applying the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian at \((u,z)\) to the vector \(w\) in direction \(w\); a dual optimization-space vector
[in]wis the direction vector; a dual constraint-space vector
[in]vis a optimization-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{ahwv} = c_{zz}(u,z)(v,\cdot)^*w\), where \(w \in \mathcal{C}^*\), \(v \in \mathcal{Z}\), and \(\mathsf{ahwv} \in \mathcal{Z}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 669 of file ROL_Constraint_PinTSimOpt.hpp.

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

Member Data Documentation

template<class Real >
Ptr<Constraint_TimeSimOpt<Real> > ROL::Constraint_PinTSimOpt< Real >::stepConstraint_
private
template<class Real >
bool ROL::Constraint_PinTSimOpt< Real >::isInitialized_
private

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