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

#include <example_02.hpp>

+ Inheritance diagram for DiffusionConstraint< Real >:

Public Member Functions

 DiffusionConstraint (const ROL::Ptr< FEM< Real > > &FEM)
 
int getNumSolves (void) const
 
void value (ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
 Evaluate the constraint operator \(c:\mathcal{U}\times\mathcal{Z} \rightarrow \mathcal{C}\) at \((u,z)\). More...
 
void solve (ROL::Vector< Real > &c, ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
 Given \(z\), solve \(c(u,z)=0\) for \(u\). More...
 
void applyJacobian_1 (ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
 Apply the partial constraint Jacobian at \((u,z)\), \(c_u(u,z) \in L(\mathcal{U}, \mathcal{C})\), to the vector \(v\). More...
 
void applyJacobian_2 (ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
 Apply the partial constraint Jacobian at \((u,z)\), \(c_z(u,z) \in L(\mathcal{Z}, \mathcal{C})\), to the vector \(v\). More...
 
void applyInverseJacobian_1 (ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
 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...
 
void applyAdjointJacobian_1 (ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, 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 primary interface. More...
 
void applyAdjointJacobian_2 (ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, 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 primary interface. More...
 
void applyInverseAdjointJacobian_1 (ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
 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...
 
void applyAdjointHessian_11 (ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
 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...
 
void applyAdjointHessian_12 (ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
 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...
 
void applyAdjointHessian_21 (ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
 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...
 
void applyAdjointHessian_22 (ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
 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::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 update_1 (const Vector< Real > &u, bool flag=true, int iter=-1)
 Update constraint functions with respect to Sim variable. x is the optimization 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. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void setSolveParameters (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< 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

\[ \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 Vector< Real > &u, const Vector< Real > &z, const 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 Member Functions

void plus (std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.0)
 
void scale (std::vector< Real > &u, const Real alpha=0.0)
 

Private Attributes

const ROL::Ptr< FEM< Real > > FEM_
 
int num_solves_
 

Additional Inherited Members

- Protected Member Functions inherited from ROL::Constraint< Real >
const std::vector< Real > getParameter (void) const
 
- Protected Attributes inherited from ROL::Constraint_SimOpt< Real >
Real atol_
 
Real rtol_
 
Real stol_
 
Real factor_
 
Real decr_
 
int maxit_
 
bool print_
 
bool zero_
 
int solverType_
 
bool firstSolve_
 

Detailed Description

template<class Real>
class DiffusionConstraint< Real >

Definition at line 403 of file poisson-control/example_02.hpp.

Constructor & Destructor Documentation

template<class Real >
DiffusionConstraint< Real >::DiffusionConstraint ( const ROL::Ptr< FEM< Real > > &  FEM)
inline

Definition at line 432 of file poisson-control/example_02.hpp.

Member Function Documentation

template<class Real >
void DiffusionConstraint< Real >::plus ( std::vector< Real > &  u,
const std::vector< Real > &  s,
const Real  alpha = 1.0 
)
inlineprivate
template<class Real >
void DiffusionConstraint< Real >::scale ( std::vector< Real > &  u,
const Real  alpha = 0.0 
)
inlineprivate
template<class Real >
int DiffusionConstraint< Real >::getNumSolves ( void  ) const
inline
template<class Real >
void DiffusionConstraint< Real >::value ( ROL::Vector< Real > &  c,
const ROL::Vector< Real > &  u,
const ROL::Vector< Real > &  z,
Real &  tol 
)
inlinevirtual

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::Constraint_SimOpt< Real >.

Definition at line 438 of file poisson-control/example_02.hpp.

References DiffusionConstraint< Real >::FEM_, ROL::Constraint< Real >::getParameter(), and DiffusionConstraint< Real >::plus().

Referenced by DiffusionConstraint< Real >::solve().

template<class Real >
void DiffusionConstraint< Real >::solve ( ROL::Vector< Real > &  c,
ROL::Vector< Real > &  u,
const ROL::Vector< Real > &  z,
Real &  tol 
)
inlinevirtual

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::Constraint_SimOpt< Real >.

Definition at line 454 of file poisson-control/example_02.hpp.

References DiffusionConstraint< Real >::FEM_, ROL::Constraint< Real >::getParameter(), DiffusionConstraint< Real >::num_solves_, DiffusionConstraint< Real >::plus(), and DiffusionConstraint< Real >::value().

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

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

  @param[out]      jv  is the result of applying the constraint Jacobian to @b v at @b  \form#196; a constraint-space vector
  @param[in]       v   is a simulation-space vector
  @param[in]       u   is the constraint argument; an simulation-space vector
  @param[in]       z   is the constraint argument; an optimization-space vector
  @param[in,out]   tol is a tolerance for inexact evaluations; currently unused

  On return, \form#201, where

\(v \in \mathcal{U}\), \(\mathsf{jv} \in \mathcal{C}\).


Reimplemented from ROL::Constraint_SimOpt< Real >.

Definition at line 476 of file poisson-control/example_02.hpp.

References DiffusionConstraint< Real >::FEM_, and ROL::Constraint< Real >::getParameter().

Referenced by DiffusionConstraint< Real >::applyAdjointJacobian_1().

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

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

  @param[out]      jv  is the result of applying the constraint Jacobian to @b v at @b  \form#196; a constraint-space vector
  @param[in]       v   is an optimization-space vector
  @param[in]       u   is the constraint argument; a simulation-space vector
  @param[in]       z   is the constraint argument; an optimization-space vector
  @param[in,out]   tol is a tolerance for inexact evaluations; currently unused

  On return, \form#204, where

\(v \in \mathcal{Z}\), \(\mathsf{jv} \in \mathcal{C}\).


Reimplemented from ROL::Constraint_SimOpt< Real >.

Definition at line 483 of file poisson-control/example_02.hpp.

References DiffusionConstraint< Real >::FEM_, ROL::Constraint< Real >::getParameter(), and DiffusionConstraint< Real >::scale().

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

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

  @param[out]      ijv is the result of applying the inverse constraint Jacobian to @b v at @b  \form#196; a simulation-space vector
  @param[in]       v   is a constraint-space vector
  @param[in]       u   is the constraint argument; a simulation-space vector
  @param[in]       z   is the constraint argument; an optimization-space vector
  @param[in,out]   tol is a tolerance for inexact evaluations; currently unused

  On return, \form#207, where

\(v \in \mathcal{C}\), \(\mathsf{ijv} \in \mathcal{U}\).


Reimplemented from ROL::Constraint_SimOpt< Real >.

Definition at line 491 of file poisson-control/example_02.hpp.

References DiffusionConstraint< Real >::FEM_, ROL::Constraint< Real >::getParameter(), and DiffusionConstraint< Real >::num_solves_.

Referenced by DiffusionConstraint< Real >::applyInverseAdjointJacobian_1().

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

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.

  @param[out]      ajv    is the result of applying the adjoint of the constraint Jacobian to @b v at @b (u,z); a dual simulation-space vector
  @param[in]       v      is a dual constraint-space vector
  @param[in]       u      is the constraint argument; a simulation-space vector
  @param[in]       z      is the constraint argument; an optimization-space vector
  @param[in,out]   tol    is a tolerance for inexact evaluations; currently unused

  On return, \form#210, where

\(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{U}^*\).


Reimplemented from ROL::Constraint_SimOpt< Real >.

Definition at line 505 of file poisson-control/example_02.hpp.

References DiffusionConstraint< Real >::applyJacobian_1().

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

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.

  @param[out]      ajv    is the result of applying the adjoint of the constraint Jacobian to @b v at @b  \form#196; a dual optimization-space vector
  @param[in]       v      is a dual constraint-space vector
  @param[in]       u      is the constraint argument; a simulation-space vector
  @param[in]       z      is the constraint argument; an optimization-space vector
  @param[in,out]   tol    is a tolerance for inexact evaluations; currently unused

  On return, \form#213, where

\(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{Z}^*\).


Reimplemented from ROL::Constraint_SimOpt< Real >.

Definition at line 510 of file poisson-control/example_02.hpp.

References DiffusionConstraint< Real >::FEM_, ROL::Constraint< Real >::getParameter(), and DiffusionConstraint< Real >::scale().

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

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\).

  @param[out]      iajv is the result of applying the inverse adjoint of the constraint Jacobian to @b v at @b (u,z); a dual constraint-space vector
  @param[in]       v   is a dual simulation-space vector
  @param[in]       u   is the constraint argument; a simulation-space vector
  @param[in]       z   is the constraint argument; an optimization-space vector
  @param[in,out]   tol is a tolerance for inexact evaluations; currently unused

  On return, \form#216, where

\(v \in \mathcal{U}^*\), \(\mathsf{iajv} \in \mathcal{C}^*\).


Reimplemented from ROL::Constraint_SimOpt< Real >.

Definition at line 518 of file poisson-control/example_02.hpp.

References DiffusionConstraint< Real >::applyInverseJacobian_1().

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

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\).

  @param[out]      ahwv is the result of applying the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at @b  \form#196 to the vector @b  \form#219 in direction @b  \form#219; a dual simulation-space vector
  @param[in]       w    is the direction vector; a dual constraint-space vector
  @param[in]       v    is a simulation-space vector
  @param[in]       u    is the constraint argument; a simulation-space vector
  @param[in]       z    is the constraint argument; an optimization-space vector
  @param[in,out]   tol  is a tolerance for inexact evaluations; currently unused

  On return, \form#221, where

\(w \in \mathcal{C}^*\), \(v \in \mathcal{U}\), and \(\mathsf{ahwv} \in \mathcal{U}^*\).


Reimplemented from ROL::Constraint_SimOpt< Real >.

Definition at line 523 of file poisson-control/example_02.hpp.

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

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

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\).

  @param[out]      ahwv is the result of applying the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at @b  \form#196 to the vector @b  \form#219 in direction @b  \form#219; a dual optimization-space vector
  @param[in]       w    is the direction vector; a dual constraint-space vector
  @param[in]       v    is a simulation-space vector
  @param[in]       u    is the constraint argument; a simulation-space vector
  @param[in]       z    is the constraint argument; an optimization-space vector
  @param[in,out]   tol  is a tolerance for inexact evaluations; currently unused

  On return, \form#225, where

\(w \in \mathcal{C}^*\), \(v \in \mathcal{U}\), and \(\mathsf{ahwv} \in \mathcal{Z}^*\).


Reimplemented from ROL::Constraint_SimOpt< Real >.

Definition at line 528 of file poisson-control/example_02.hpp.

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

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

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\).

  @param[out]      ahwv is the result of applying the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at @b  \form#196 to the vector @b  \form#219 in direction @b  \form#219; a dual simulation-space vector
  @param[in]       w    is the direction vector; a dual constraint-space vector
  @param[in]       v    is a optimization-space vector
  @param[in]       u    is the constraint argument; a simulation-space vector
  @param[in]       z    is the constraint argument; an optimization-space vector
  @param[in,out]   tol  is a tolerance for inexact evaluations; currently unused

  On return, \form#228, where

\(w \in \mathcal{C}^*\), \(v \in \mathcal{Z}\), and \(\mathsf{ahwv} \in \mathcal{U}^*\).


Reimplemented from ROL::Constraint_SimOpt< Real >.

Definition at line 533 of file poisson-control/example_02.hpp.

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

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

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\).

  @param[out]      ahwv is the result of applying the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian at @b  \form#196 to the vector @b  \form#219 in direction @b  \form#219; a dual optimization-space vector
  @param[in]       w    is the direction vector; a dual constraint-space vector
  @param[in]       v    is a optimization-space vector
  @param[in]       u    is the constraint argument; a simulation-space vector
  @param[in]       z    is the constraint argument; an optimization-space vector
  @param[in,out]   tol  is a tolerance for inexact evaluations; currently unused

  On return, \form#230, where

\(w \in \mathcal{C}^*\), \(v \in \mathcal{Z}\), and \(\mathsf{ahwv} \in \mathcal{Z}^*\).


Reimplemented from ROL::Constraint_SimOpt< Real >.

Definition at line 538 of file poisson-control/example_02.hpp.

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

Member Data Documentation

template<class Real >
const ROL::Ptr<FEM<Real> > DiffusionConstraint< Real >::FEM_
private
template<class Real >
int DiffusionConstraint< Real >::num_solves_
private

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