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

#include <example_01.hpp>

+ Inheritance diagram for Objective_GrossPitaevskii< Real >:

Public Member Functions

 Objective_GrossPitaevskii (const Real &g, const Vector< Real > &V)
 
Real value (const Vector< Real > &psi, Real &tol)
 Evaluate \(J[\psi]\). More...
 
void gradient (Vector< Real > &g, const Vector< Real > &psi, Real &tol)
 Evaluate \(\nabla J[\psi]\). More...
 
void hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
 Evaluate \(\nabla^2 J[\psi] v\). More...
 
 Objective_GrossPitaevskii (const Real &g, const Vector< Real > &V, ROL::Ptr< FiniteDifference< Real > > fd)
 
Real value (const Vector< Real > &psi, Real &tol)
 Evaluate \(J[\psi]\). More...
 
void gradient (Vector< Real > &g, const Vector< Real > &psi, Real &tol)
 Evaluate \(\nabla J[\psi]\). More...
 
void hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
 Evaluate \(\nabla^2 J[\psi] v\). More...
 
- Public Member Functions inherited from ROL::Objective< Real >
virtual ~Objective ()
 
virtual void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update objective function. More...
 
virtual Real dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol)
 Compute directional derivative. More...
 
virtual void invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply inverse Hessian approximation to vector. More...
 
virtual void precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply preconditioner to vector. More...
 
virtual std::vector
< std::vector< Real > > 
checkGradient (const Vector< Real > &x, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference gradient check. More...
 
virtual std::vector
< std::vector< Real > > 
checkGradient (const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference gradient check. More...
 
virtual std::vector
< std::vector< Real > > 
checkGradient (const Vector< Real > &x, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference gradient check with specified step sizes. More...
 
virtual std::vector
< std::vector< Real > > 
checkGradient (const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference gradient check with specified step sizes. More...
 
virtual std::vector
< std::vector< Real > > 
checkHessVec (const Vector< Real > &x, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference Hessian-applied-to-vector check. More...
 
virtual std::vector
< std::vector< Real > > 
checkHessVec (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference Hessian-applied-to-vector check. More...
 
virtual std::vector
< std::vector< Real > > 
checkHessVec (const Vector< Real > &x, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference Hessian-applied-to-vector check with specified step sizes. More...
 
virtual std::vector
< std::vector< Real > > 
checkHessVec (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference Hessian-applied-to-vector check with specified step sizes. More...
 
virtual std::vector< Real > checkHessSym (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout)
 Hessian symmetry check. More...
 
virtual std::vector< Real > checkHessSym (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout)
 Hessian symmetry check. 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)
 
void applyK (const Vector< Real > &v, Vector< Real > &kv)
 Apply finite difference operator. More...
 
void applyK (const Vector< Real > &v, Vector< Real > &kv)
 Apply finite difference operator. More...
 

Private Attributes

Real g_
 
uint nx_
 
Real dx_
 
ROL::Ptr< const vectorVp_
 
ROL::Ptr< const std::vector
< Real > > 
Vp_
 
ROL::Ptr< FiniteDifference
< Real > > 
fd_
 

Additional Inherited Members

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

Detailed Description

template<class Real>
class Objective_GrossPitaevskii< Real >

Objective Function Class

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

Member Typedef Documentation

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

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

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

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

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

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

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

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

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

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

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

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

Constructor & Destructor Documentation

template<class Real>
Objective_GrossPitaevskii< Real >::Objective_GrossPitaevskii ( const Real &  g,
const Vector< Real > &  V 
)
inline

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

template<class Real>
Objective_GrossPitaevskii< Real >::Objective_GrossPitaevskii ( const Real &  g,
const Vector< Real > &  V,
ROL::Ptr< FiniteDifference< Real > >  fd 
)
inline

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

Member Function Documentation

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

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

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

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

template<class Real>
void Objective_GrossPitaevskii< Real >::applyK ( const Vector< Real > &  v,
Vector< Real > &  kv 
)
inlineprivate

Apply finite difference operator.

Compute \(K\psi\), where \(K\) is the finite difference approximation of \(-D_x^2\)

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

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

Evaluate \(J[\psi]\).

\[ J[\psi]=\frac{1}{2} \int\limits_0^1 |\psi'|^2 + V(x)|\psi|^2+g|\psi|^4\,\mathrm{d}x \]

where the integral is approximated with the trapezoidal rule and the derivative is approximated using finite differences

Implements ROL::Objective< Real >.

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

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

template<class Real>
void Objective_GrossPitaevskii< Real >::gradient ( Vector< Real > &  g,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \(\nabla J[\psi]\).

\[ \nabla J[\psi] = -\psi'' + V(x)\psi+2g|\psi|^3 \]

Reimplemented from ROL::Objective< Real >.

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

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

template<class Real>
void Objective_GrossPitaevskii< Real >::hessVec ( Vector< Real > &  hv,
const Vector< Real > &  v,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \(\nabla^2 J[\psi] v\).

\[ \nabla^2 J[\psi]v = -v'' + V(x)v+6g|\psi|^2 v \]

Reimplemented from ROL::Objective< Real >.

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

template<class Real>
void Objective_GrossPitaevskii< Real >::applyK ( const Vector< Real > &  v,
Vector< Real > &  kv 
)
inlineprivate

Apply finite difference operator.

Compute \(K\psi\), where \(K\) is the finite difference approximation of \(-D_x^2\)

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

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

Evaluate \(J[\psi]\).

\[ J[\psi]=\frac{1}{2} \int\limits_0^1 |\psi'|^2 + V(x)|\psi|^2+g|\psi|^4\,\mathrm{d}x \]

where the integral is approximated with the trapezoidal rule and the derivative is approximated using finite differences

Implements ROL::Objective< Real >.

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

template<class Real>
void Objective_GrossPitaevskii< Real >::gradient ( Vector< Real > &  g,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \(\nabla J[\psi]\).

\[ \nabla J[\psi] = -\psi'' + V(x)\psi+2g|\psi|^3 \]

Reimplemented from ROL::Objective< Real >.

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

template<class Real>
void Objective_GrossPitaevskii< Real >::hessVec ( Vector< Real > &  hv,
const Vector< Real > &  v,
const Vector< Real > &  psi,
Real &  tol 
)
inlinevirtual

Evaluate \(\nabla^2 J[\psi] v\).

\[ \nabla^2 J[\psi]v = -v'' + V(x)v+6g|\psi|^2 v \]

Reimplemented from ROL::Objective< Real >.

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

Member Data Documentation

template<class Real>
Real Objective_GrossPitaevskii< Real >::g_
private

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

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

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

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

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

template<class Real>
ROL::Ptr<const vector> Objective_GrossPitaevskii< Real >::Vp_
private

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

template<class Real>
ROL::Ptr<const std::vector<Real> > Objective_GrossPitaevskii< Real >::Vp_
private

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

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

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


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