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

Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraints. More...

#include <ROL_StdVector.hpp>

+ Inheritance diagram for ROL::StdVector< Real, Element >:

Public Member Functions

 StdVector (const Ptr< std::vector< Element >> &std_vec)
 
 StdVector (const int dim, const Element val=0.0)
 
 StdVector (std::initializer_list< Element > ilist)
 
Real & operator[] (int i)
 
const Real & operator[] (int i) const
 
void set (const Vector< Real > &x)
 Set \(y \leftarrow x\) where \(y = \mathtt{*this}\). More...
 
void plus (const Vector< Real > &x)
 Compute \(y \leftarrow y + x\), where \(y = \mathtt{*this}\). More...
 
void axpy (const Real alpha, const Vector< Real > &x)
 Compute \(y \leftarrow \alpha x + y\) where \(y = \mathtt{*this}\). More...
 
void scale (const Real alpha)
 Compute \(y \leftarrow \alpha y\) where \(y = \mathtt{*this}\). More...
 
virtual Real dot (const Vector< Real > &x) const
 Compute \( \langle y,x \rangle \) where \(y = \mathtt{*this}\). More...
 
Real norm () const
 Returns \( \| y \| \) where \(y = \mathtt{*this}\). More...
 
virtual Ptr< Vector< Real > > clone () const
 Clone to make a new (uninitialized) vector. More...
 
Ptr< const std::vector< Element > > getVector () const
 
Ptr< std::vector< Element > > getVector ()
 
Ptr< Vector< Real > > basis (const int i) const
 Return i-th basis vector. More...
 
int dimension () const
 Return dimension of the vector space. More...
 
void applyUnary (const Elementwise::UnaryFunction< Real > &f)
 
void applyBinary (const Elementwise::BinaryFunction< Real > &f, const Vector< Real > &x)
 
Real reduce (const Elementwise::ReductionOp< Real > &r) const
 
void setScalar (const Real C)
 Set \(y \leftarrow C\) where \(C\in\mathbb{R}\). More...
 
void randomize (const Real l=0.0, const Real u=1.0)
 Set vector to be uniform random between [l,u]. More...
 
virtual void print (std::ostream &outStream) const
 
- Public Member Functions inherited from ROL::Vector< Real >
virtual ~Vector ()
 
virtual void zero ()
 Set to zero vector. More...
 
virtual const Vectordual () const
 Return dual representation of \(\mathtt{*this}\), for example, the result of applying a Riesz map, or change of basis, or change of memory layout. More...
 
virtual Real apply (const Vector< Real > &x) const
 Apply \(\mathtt{*this}\) to a dual vector. This is equivalent to the call \(\mathtt{this->dot(x.dual())}\). More...
 
virtual std::vector< Real > checkVector (const Vector< Real > &x, const Vector< Real > &y, const bool printToStream=true, std::ostream &outStream=std::cout) const
 Verify vector-space methods. More...
 

Private Types

using size_type = typename std::vector< Real >::size_type
 

Private Attributes

Ptr< std::vector< Element > > std_vec_
 

Detailed Description

template<class Real, class Element = Real>
class ROL::StdVector< Real, Element >

Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraints.

Provides the std::vector implementation of the ROL::Vector interface.

Definition at line 61 of file ROL_StdVector.hpp.

Member Typedef Documentation

template<class Real, class Element = Real>
using ROL::StdVector< Real, Element >::size_type = typename std::vector<Real>::size_type
private

Definition at line 63 of file ROL_StdVector.hpp.

Constructor & Destructor Documentation

template<class Real, class Element = Real>
ROL::StdVector< Real, Element >::StdVector ( const Ptr< std::vector< Element >> &  std_vec)
inline

Definition at line 67 of file ROL_StdVector.hpp.

template<class Real, class Element = Real>
ROL::StdVector< Real, Element >::StdVector ( const int  dim,
const Element  val = 0.0 
)
inline

Definition at line 69 of file ROL_StdVector.hpp.

template<class Real, class Element = Real>
ROL::StdVector< Real, Element >::StdVector ( std::initializer_list< Element >  ilist)
inline

Definition at line 73 of file ROL_StdVector.hpp.

Member Function Documentation

template<class Real, class Element = Real>
Real& ROL::StdVector< Real, Element >::operator[] ( int  i)
inline

Definition at line 76 of file ROL_StdVector.hpp.

template<class Real, class Element = Real>
const Real& ROL::StdVector< Real, Element >::operator[] ( int  i) const
inline

Definition at line 77 of file ROL_StdVector.hpp.

template<class Real, class Element = Real>
void ROL::StdVector< Real, Element >::set ( const Vector< Real > &  x)
inlinevirtual

Set \(y \leftarrow x\) where \(y = \mathtt{*this}\).

Parameters
[in]xis a vector.

On return \(\mathtt{*this} = x\). Uses zero and plus methods for the computation. Please overload if a more efficient implementation is needed.


Reimplemented from ROL::Vector< Real >.

Definition at line 79 of file ROL_StdVector.hpp.

Referenced by main(), and testRandomInputs().

template<class Real, class Element = Real>
void ROL::StdVector< Real, Element >::plus ( const Vector< Real > &  x)
inlinevirtual

Compute \(y \leftarrow y + x\), where \(y = \mathtt{*this}\).

Parameters
[in]xis the vector to be added to \(\mathtt{*this}\).

On return \(\mathtt{*this} = \mathtt{*this} + x\).


Implements ROL::Vector< Real >.

Definition at line 90 of file ROL_StdVector.hpp.

template<class Real, class Element = Real>
void ROL::StdVector< Real, Element >::axpy ( const Real  alpha,
const Vector< Real > &  x 
)
inlinevirtual

Compute \(y \leftarrow \alpha x + y\) where \(y = \mathtt{*this}\).

Parameters
[in]alphais the scaling of x.
[in]xis a vector.

On return \(\mathtt{*this} = \mathtt{*this} + \alpha x \). Uses clone, set, scale and plus for the computation. Please overload if a more efficient implementation is needed.


Reimplemented from ROL::Vector< Real >.

Definition at line 104 of file ROL_StdVector.hpp.

Referenced by main(), ROL::ZOO::Objective_PoissonInversion< Real >::solve_adjoint_sensitivity_equation(), and testRandomInputs().

template<class Real, class Element = Real>
void ROL::StdVector< Real, Element >::scale ( const Real  alpha)
inlinevirtual

Compute \(y \leftarrow \alpha y\) where \(y = \mathtt{*this}\).

Parameters
[in]alphais the scaling of \(\mathtt{*this}\).

On return \(\mathtt{*this} = \alpha (\mathtt{*this}) \).


Implements ROL::Vector< Real >.

Definition at line 118 of file ROL_StdVector.hpp.

Referenced by main(), and testRandomInputs().

template<class Real, class Element = Real>
virtual Real ROL::StdVector< Real, Element >::dot ( const Vector< Real > &  x) const
inlinevirtual
template<class Real, class Element = Real>
Real ROL::StdVector< Real, Element >::norm ( ) const
inlinevirtual

Returns \( \| y \| \) where \(y = \mathtt{*this}\).

Returns
A nonnegative number equal to the norm of \(\mathtt{*this}\).

Implements ROL::Vector< Real >.

Definition at line 143 of file ROL_StdVector.hpp.

Referenced by main().

template<class Real, class Element = Real>
virtual Ptr<Vector<Real> > ROL::StdVector< Real, Element >::clone ( ) const
inlinevirtual
template<class Real, class Element = Real>
Ptr<const std::vector<Element> > ROL::StdVector< Real, Element >::getVector ( void  ) const
inline

Definition at line 153 of file ROL_StdVector.hpp.

Referenced by TridiagonalToeplitzOperator< Real >::apply(), ROL::StdConstraint< Real >::applyAdjointHessian(), ROL::StdConstraint< Real >::applyAdjointJacobian(), ROL::StdVector< Real >::applyBinary(), TridiagonalToeplitzOperator< Real >::applyInverse(), ROL::StdConstraint< Real >::applyJacobian(), ROL::StdConstraint< Real >::applyPreconditioner(), ROL::StdVector< Real >::axpy(), ROL::StdVector< Real >::basis(), ROL::BatchStdVector< Real >::clone(), ROL::ProbabilityVector< Real >::clone(), ROL::PrimalScaledStdVector< Real, Element >::clone(), ROL::PrimalProbabilityVector< Real >::clone(), ROL::DualScaledStdVector< Real, Element >::clone(), ROL::DualProbabilityVector< Real >::clone(), L2VectorPrimal< Real >::clone(), L2VectorDual< Real >::clone(), H1VectorPrimal< Real >::clone(), H1VectorDual< Real >::clone(), ROL::StdObjective< Real >::dirDeriv(), ROL::BatchStdVector< Real >::dot(), ROL::PrimalScaledStdVector< Real, Element >::dot(), ROL::PrimalProbabilityVector< Real >::dot(), ROL::PrimalAtomVector< Real >::dot(), ROL::StdVector< Real >::dot(), ROL::DualScaledStdVector< Real, Element >::dot(), ROL::DualProbabilityVector< Real >::dot(), ROL::DualAtomVector< Real >::dot(), L2VectorPrimal< Real >::dot(), L2VectorDual< Real >::dot(), H1VectorPrimal< Real >::dot(), H1VectorDual< Real >::dot(), ROL::PrimalScaledStdVector< Real, Element >::dual(), ROL::PrimalProbabilityVector< Real >::dual(), ROL::PrimalAtomVector< Real >::dual(), ROL::DualScaledStdVector< Real, Element >::dual(), ROL::DualProbabilityVector< Real >::dual(), ROL::DualAtomVector< Real >::dual(), L2VectorPrimal< Real >::dual(), L2VectorDual< Real >::dual(), H1VectorPrimal< Real >::dual(), H1VectorDual< Real >::dual(), ROL::AtomVector< Real >::getAtom(), ROL::ProbabilityVector< Real >::getNumMyAtoms(), ROL::ProbabilityVector< Real >::getProbability(), ROL::StdObjective< Real >::gradient(), ROL::StdObjective< Real >::hessVec(), ROL::StdObjective< Real >::invHessVec(), ROL::StdVector< Real >::plus(), ROL::StdObjective< Real >::precond(), ROL::BatchStdVector< Real >::reduce(), ROL::StdVector< Real >::set(), ROL::AtomVector< Real >::setAtom(), ROL::ProbabilityVector< Real >::setProbability(), ROL::StdConstraint< Real >::solveAugmentedSystem(), ROL::StdConstraint< Real >::update(), ROL::StdObjective< Real >::update(), ROL::StdConstraint< Real >::value(), and ROL::StdObjective< Real >::value().

template<class Real, class Element = Real>
Ptr<std::vector<Element> > ROL::StdVector< Real, Element >::getVector ( void  )
inline

Definition at line 157 of file ROL_StdVector.hpp.

template<class Real, class Element = Real>
Ptr<Vector<Real> > ROL::StdVector< Real, Element >::basis ( const int  i) const
inlinevirtual

Return i-th basis vector.

Parameters
[in]iis the index of the basis function.
Returns
A reference-counted pointer to the basis vector with index i.

Overloading the basis is only required if the default gradient implementation is used, which computes a finite-difference approximation.


Reimplemented from ROL::Vector< Real >.

Definition at line 161 of file ROL_StdVector.hpp.

Referenced by main().

template<class Real, class Element = Real>
int ROL::StdVector< Real, Element >::dimension ( void  ) const
inlinevirtual
template<class Real, class Element = Real>
void ROL::StdVector< Real, Element >::applyUnary ( const Elementwise::UnaryFunction< Real > &  f)
inlinevirtual

Reimplemented from ROL::Vector< Real >.

Definition at line 175 of file ROL_StdVector.hpp.

Referenced by testRandomInputs().

template<class Real, class Element = Real>
void ROL::StdVector< Real, Element >::applyBinary ( const Elementwise::BinaryFunction< Real > &  f,
const Vector< Real > &  x 
)
inlinevirtual

Reimplemented from ROL::Vector< Real >.

Definition at line 183 of file ROL_StdVector.hpp.

Referenced by testCases().

template<class Real, class Element = Real>
Real ROL::StdVector< Real, Element >::reduce ( const Elementwise::ReductionOp< Real > &  r) const
inlinevirtual

Reimplemented from ROL::Vector< Real >.

Definition at line 198 of file ROL_StdVector.hpp.

template<class Real, class Element = Real>
void ROL::StdVector< Real, Element >::setScalar ( const Real  C)
inlinevirtual

Set \(y \leftarrow C\) where \(C\in\mathbb{R}\).

Parameters
[in]Cis a scalar.

On return \(\mathtt{*this} = C\). Uses applyUnary methods for the computation. Please overload if a more efficient implementation is needed.


Reimplemented from ROL::Vector< Real >.

Definition at line 207 of file ROL_StdVector.hpp.

template<class Real, class Element = Real>
void ROL::StdVector< Real, Element >::randomize ( const Real  l = 0.0,
const Real  u = 1.0 
)
inlinevirtual

Set vector to be uniform random between [l,u].

Parameters
[in]lis a the lower bound.
[in]uis a the upper bound.

On return the components of \(\mathtt{*this}\) are uniform random numbers on the interval \([l,u]\). The default implementation uses applyUnary methods for the computation. Please overload if a more efficient implementation is needed.


Reimplemented from ROL::Vector< Real >.

Definition at line 212 of file ROL_StdVector.hpp.

Referenced by main(), and testRandomInputs().

template<class Real, class Element = Real>
virtual void ROL::StdVector< Real, Element >::print ( std::ostream &  outStream) const
inlinevirtual

Reimplemented from ROL::Vector< Real >.

Definition at line 228 of file ROL_StdVector.hpp.

Member Data Documentation

template<class Real, class Element = Real>
Ptr<std::vector<Element> > ROL::StdVector< Real, Element >::std_vec_
private

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