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

#include <example_01.hpp>

+ Inheritance diagram for Normalization_Constraint< Real >:

Public Member Functions

 Normalization_Constraint (int n, Real dx)
 
void value (Vector< Real > &c, const Vector< Real > &psi, Real &tol)
 Evaluate \(c[\psi]\). More...
 
void applyJacobian (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
 Evaluate \(c'[\psi]v\). More...
 
void applyAdjointJacobian (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
 Evaluate \((c'[\psi])^\ast v\). More...
 
void applyAdjointHessian (Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
 Evaluate \(((c''[\psi])^\ast v)u\). More...
 
 Normalization_Constraint (int n, Real dx, ROL::Ptr< FiniteDifference< Real > > fd, bool exactsolve)
 
void value (Vector< Real > &c, const Vector< Real > &psi, Real &tol)
 Evaluate \(c[\psi]\). More...
 
void applyJacobian (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
 Evaluate \(c'[\psi]v\). More...
 
void applyAdjointJacobian (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
 Evaluate \((c'[\psi])^\ast v\). More...
 
void applyAdjointHessian (Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
 Evaluate \(((c''[\psi])^\ast v)u\). More...
 
std::vector< Real > solveAugmentedSystem (Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &psi, Real &tol)
 
- Public Member Functions inherited from ROL::Constraint< Real >
virtual ~Constraint (void)
 
 Constraint (void)
 
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 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...
 
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:

\[ \left[c'(x) \circ R \circ c'(x)^* \circ P(x)\right] v = v \,, \]

where R is the appropriate Riesz map in \(L(\mathcal{X}^*, \mathcal{X})\). It is used by the solveAugmentedSystem method. 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

typedef std::vector< Real > vector
 
typedef Vector< Real > V
 
typedef StdVector< Real > SV
 
typedef vector::size_type uint
 
typedef std::vector< Real > vector
 
typedef vector::size_type uint
 

Private Member Functions

ROL::Ptr< const vectorgetVector (const V &x)
 
ROL::Ptr< vectorgetVector (V &x)
 

Private Attributes

uint nx_
 
Real dx_
 
ROL::Ptr< FiniteDifference
< Real > > 
fd_
 
bool exactsolve_
 

Additional Inherited Members

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

Detailed Description

template<class Real>
class Normalization_Constraint< Real >

Constraint class

Definition at line 237 of file gross-pitaevskii/example_01.hpp.

Member Typedef Documentation

template<class Real>
typedef std::vector<Real> Normalization_Constraint< Real >::vector
private

Definition at line 239 of file gross-pitaevskii/example_01.hpp.

template<class Real>
typedef Vector<Real> Normalization_Constraint< Real >::V
private

Definition at line 240 of file gross-pitaevskii/example_01.hpp.

template<class Real>
typedef StdVector<Real> Normalization_Constraint< Real >::SV
private

Definition at line 241 of file gross-pitaevskii/example_01.hpp.

template<class Real>
typedef vector::size_type Normalization_Constraint< Real >::uint
private

Definition at line 243 of file gross-pitaevskii/example_01.hpp.

template<class Real>
typedef std::vector<Real> Normalization_Constraint< Real >::vector
private

Definition at line 597 of file gross-pitaevskii/example_02.hpp.

template<class Real>
typedef vector::size_type Normalization_Constraint< Real >::uint
private

Definition at line 598 of file gross-pitaevskii/example_02.hpp.

Constructor & Destructor Documentation

template<class Real>
Normalization_Constraint< Real >::Normalization_Constraint ( int  n,
Real  dx 
)
inline

Definition at line 261 of file gross-pitaevskii/example_01.hpp.

template<class Real>
Normalization_Constraint< Real >::Normalization_Constraint ( int  n,
Real  dx,
ROL::Ptr< FiniteDifference< Real > >  fd,
bool  exactsolve 
)
inline

Definition at line 607 of file gross-pitaevskii/example_02.hpp.

Member Function Documentation

template<class Real>
ROL::Ptr<const vector> Normalization_Constraint< Real >::getVector ( const V x)
inlineprivate

Definition at line 250 of file gross-pitaevskii/example_01.hpp.

template<class Real>
ROL::Ptr<vector> Normalization_Constraint< Real >::getVector ( V x)
inlineprivate

Definition at line 255 of file gross-pitaevskii/example_01.hpp.

template<class Real>
void Normalization_Constraint< Real >::value ( Vector< Real > &  c,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \(c[\psi]\).

\[ c[\psi]= \int\limits_0^1 |\psi|^2\,\mathrm{d}x - 1 \]

where the integral is approximated with the trapezoidal rule and the derivative is approximated using finite differences. This constraint is a scalar

Implements ROL::Constraint< Real >.

Definition at line 268 of file gross-pitaevskii/example_01.hpp.

template<class Real>
void Normalization_Constraint< Real >::applyJacobian ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \(c'[\psi]v\).

\[ c'[\psi]v= 2 \int\limits_0^1 \psi v\,\mathrm{d}x \]

The action of the Jacobian on a vector produces a scalar

Reimplemented from ROL::Constraint< Real >.

Definition at line 287 of file gross-pitaevskii/example_01.hpp.

template<class Real>
void Normalization_Constraint< Real >::applyAdjointJacobian ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \((c'[\psi])^\ast v\).

\[ (c'[\psi])^\ast v = 2 \int\limits_0^1 \psi v\,\mathrm{d}x \]

The action of the Jacobian adjoint on a scalar produces a vector

Reimplemented from ROL::Constraint< Real >.

Definition at line 309 of file gross-pitaevskii/example_01.hpp.

template<class Real>
void Normalization_Constraint< Real >::applyAdjointHessian ( Vector< Real > &  ahuv,
const Vector< Real > &  u,
const Vector< Real > &  v,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \(((c''[\psi])^\ast v)u\).

\[ ((c''[\psi])^\ast v)u = 2 v u \]

The action of the Hessian adjoint on a on a vector v in a direction u produces a vector of the same size as \(\psi\)

Reimplemented from ROL::Constraint< Real >.

Definition at line 331 of file gross-pitaevskii/example_01.hpp.

template<class Real>
void Normalization_Constraint< Real >::value ( Vector< Real > &  c,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \(c[\psi]\).

\[ c[\psi]= \int\limits_0^1 |\psi|^2\,\mathrm{d}x - 1 \]

where the integral is approximated with the trapezoidal rule and the derivative is approximated using finite differences. This constraint is a scalar

Implements ROL::Constraint< Real >.

Definition at line 615 of file gross-pitaevskii/example_02.hpp.

template<class Real>
void Normalization_Constraint< Real >::applyJacobian ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \(c'[\psi]v\).

\[ c'[\psi]v= 2 \int\limits_0^1 \psi v\,\mathrm{d}x \]

The action of the Jacobian on a vector produces a scalar

Reimplemented from ROL::Constraint< Real >.

Definition at line 634 of file gross-pitaevskii/example_02.hpp.

template<class Real>
void Normalization_Constraint< Real >::applyAdjointJacobian ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \((c'[\psi])^\ast v\).

\[ (c'[\psi])^\ast v = 2 \int\limits_0^1 \psi v\,\mathrm{d}x \]

The action of the Jacobian adjoint on a scalar produces a vector

Reimplemented from ROL::Constraint< Real >.

Definition at line 656 of file gross-pitaevskii/example_02.hpp.

template<class Real>
void Normalization_Constraint< Real >::applyAdjointHessian ( Vector< Real > &  ahuv,
const Vector< Real > &  u,
const Vector< Real > &  v,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \(((c''[\psi])^\ast v)u\).

\[ ((c''[\psi])^\ast v)u = 2 v u \]

The action of the Hessian adjoint on a on a vector v in a direction u produces a vector of the same size as \(\psi\)

Reimplemented from ROL::Constraint< Real >.

Definition at line 678 of file gross-pitaevskii/example_02.hpp.

template<class Real>
std::vector<Real> Normalization_Constraint< Real >::solveAugmentedSystem ( Vector< Real > &  v1,
Vector< Real > &  v2,
const Vector< Real > &  b1,
const Vector< Real > &  b2,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Solve the system

\[ \begin{\pmatrix} K & c'^\ast(\psi)\\ c'(\psi) & 0 \end{pmatrix} \begin{pmatrix} v_1\\v_2 \end{pmatrix}=\begin{pmatrix} b_1\\b_2\end{pmatrix}\]

In this example, \(K\) is the finite difference Laplacian the constraint is a scalar and the Jacobian is a vector and the exact inverse can be computed using the Schur complement method

Reimplemented from ROL::Constraint< Real >.

Definition at line 706 of file gross-pitaevskii/example_02.hpp.

References ROL::Vector< Real >::plus(), ROL::Vector< Real >::scale(), ROL::Vector< Real >::set(), and ROL::Constraint< Real >::solveAugmentedSystem().

Member Data Documentation

template<class Real>
uint Normalization_Constraint< Real >::nx_
private

Definition at line 247 of file gross-pitaevskii/example_01.hpp.

template<class Real>
Real Normalization_Constraint< Real >::dx_
private

Definition at line 248 of file gross-pitaevskii/example_01.hpp.

template<class Real>
ROL::Ptr<FiniteDifference<Real> > Normalization_Constraint< Real >::fd_
private

Definition at line 603 of file gross-pitaevskii/example_02.hpp.

template<class Real>
bool Normalization_Constraint< Real >::exactsolve_
private

Definition at line 604 of file gross-pitaevskii/example_02.hpp.


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