Reference documentation for deal.II version GIT relicensing-245-g36f19064f7 2024-03-29 07:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Private Types | Private Attributes | List of all members
Rol::VectorAdaptor< VectorType > Class Template Reference

#include <deal.II/optimization/rol/vector_adaptor.h>

Inheritance diagram for Rol::VectorAdaptor< VectorType >:
Inheritance graph
[legend]

Public Member Functions

 VectorAdaptor (const Teuchos::RCP< VectorType > &vector_ptr)
 
Teuchos::RCP< VectorType > getVector ()
 
Teuchos::RCP< const VectorType > getVector () const
 
int dimension () const
 
void set (const ROL::Vector< value_type > &rol_vector)
 
void plus (const ROL::Vector< value_type > &rol_vector)
 
void axpy (const value_type alpha, const ROL::Vector< value_type > &rol_vector)
 
void scale (const value_type alpha)
 
value_type dot (const ROL::Vector< value_type > &rol_vector) const
 
value_type norm () const
 
Teuchos::RCP< ROL::Vector< value_type > > clone () const
 
Teuchos::RCP< ROL::Vector< value_type > > basis (const int i) const
 
void applyUnary (const ROL::Elementwise::UnaryFunction< value_type > &f)
 
void applyBinary (const ROL::Elementwise::BinaryFunction< value_type > &f, const ROL::Vector< value_type > &rol_vector)
 
value_type reduce (const ROL::Elementwise::ReductionOp< value_type > &r) const
 
void print (std::ostream &outStream) const
 

Private Types

using size_type = typename VectorType::size_type
 
using value_type = typename VectorType::value_type
 
using real_type = typename VectorType::real_type
 

Private Attributes

Teuchos::RCP< VectorType > vector_ptr
 

Detailed Description

template<typename VectorType>
class Rol::VectorAdaptor< VectorType >

An adaptor that provides the implementation of the ROL::Vector interface for vectors of type VectorType.

This class supports vectors that satisfy the following requirements:

The VectorType should contain the following types.

VectorType::size_type; // The type for size of the vector.
VectorType::value_type; // The type for elements stored in the vector.
VectorType::real_type; // The type for real-valued numbers.

However, ROL doesn't distinguish VectorAdaptor::value_type from VectorAdaptor::real_type. This is due to ROL's assumption that the VectorAdaptor::value_type itself is a type for real-valued numbers. Therefore, VectorAdaptor supports vectors whose real_type is convertible to value_type in the sense that std::is_convertible_v<real_type, value_type> yields true.

The VectorType should contain the following methods.

// Reinitialize the current vector using a given vector's
// size (and the parallel distribution) without copying
// the elements.
VectorType::reinit(const VectorType &, ...);
// Globally add a given vector to the current.
VectorType::operator+=(const VectorType &);
// Scale all elements by a given scalar.
VectorType::operator*=(const VectorType::value_type &);
// Perform dot product with a given vector.
VectorType::operator*=(const VectorType &);
// Scale all elements of the current vector and globally
// add a given vector to it.
VectorType::add(const VectorType::value_type, const VectorType &);
// Copies the data of a given vector to the current.
// Resize the current vector if necessary (MPI safe).
VectorType::operation=(const VectorType &);
// Return the global size of the current vector.
VectorType::size();
// Return L^2 norm of the current vector
VectorType::l2_norm();
// Iterator to the start of the (locally owned) element
// of the current vector.
VectorType::begin();
// Iterator to the one past the last (locally owned)
// element of the current vector.
VectorType::end();
// Compress the vector i.e., flush the buffers of the
// vector object if it has any.
VectorType::compress(VectorOperation::insert);
Note
The current implementation in ROL doesn't support vector sizes above the largest value of int type. Some of the vectors in deal.II (see Vector) may not satisfy the above requirements.

Definition at line 110 of file vector_adaptor.h.

Member Typedef Documentation

◆ size_type

template<typename VectorType >
using Rol::VectorAdaptor< VectorType >::size_type = typename VectorType::size_type
private

An alias for size type of VectorType.

Definition at line 115 of file vector_adaptor.h.

◆ value_type

template<typename VectorType >
using Rol::VectorAdaptor< VectorType >::value_type = typename VectorType::value_type
private

An alias for element type stored in the VectorType.

Definition at line 120 of file vector_adaptor.h.

◆ real_type

template<typename VectorType >
using Rol::VectorAdaptor< VectorType >::real_type = typename VectorType::real_type
private

An alias for real-valued numbers.

Definition at line 125 of file vector_adaptor.h.

Constructor & Destructor Documentation

◆ VectorAdaptor()

template<typename VectorType >
Rol::VectorAdaptor< VectorType >::VectorAdaptor ( const Teuchos::RCP< VectorType > &  vector_ptr)

Constructor.

Member Function Documentation

◆ getVector() [1/2]

template<typename VectorType >
Teuchos::RCP< VectorType > Rol::VectorAdaptor< VectorType >::getVector ( )

Return the Teuchos smart reference counting pointer to the wrapper vector, vector_ptr.

◆ getVector() [2/2]

template<typename VectorType >
Teuchos::RCP< const VectorType > Rol::VectorAdaptor< VectorType >::getVector ( ) const

Return the Teuchos smart reference counting pointer to const vector.

◆ dimension()

template<typename VectorType >
int Rol::VectorAdaptor< VectorType >::dimension ( ) const

Return the dimension (global vector size) of the wrapped vector.

◆ set()

template<typename VectorType >
void Rol::VectorAdaptor< VectorType >::set ( const ROL::Vector< value_type > &  rol_vector)

Set the wrapper vector to a given ROL::Vector rol_vector by overwriting its contents.

If the current wrapper vector has ghost elements, then VectorType::operator=(const VectorType&) should still be allowed on it.

◆ plus()

template<typename VectorType >
void Rol::VectorAdaptor< VectorType >::plus ( const ROL::Vector< value_type > &  rol_vector)

Perform addition.

◆ axpy()

template<typename VectorType >
void Rol::VectorAdaptor< VectorType >::axpy ( const value_type  alpha,
const ROL::Vector< value_type > &  rol_vector 
)

Scale the wrapper vector by alpha and add ROL::Vector rol_vector to it.

◆ scale()

template<typename VectorType >
void Rol::VectorAdaptor< VectorType >::scale ( const value_type  alpha)

Scale the wrapper vector.

◆ dot()

template<typename VectorType >
value_type Rol::VectorAdaptor< VectorType >::dot ( const ROL::Vector< value_type > &  rol_vector) const

Return the dot product with a given ROL::Vector rol_vector.

◆ norm()

template<typename VectorType >
value_type Rol::VectorAdaptor< VectorType >::norm ( ) const

Return the \(L^{2}\) norm of the wrapped vector.

The returned type is of VectorAdaptor::value_type so as to maintain consistency with ROL::Vector<VectorAdaptor::value_type> and more importantly to not to create an overloaded version namely, VectorAdaptor::real_type norm() const; if real_type and value_type are not of the same type.

◆ clone()

template<typename VectorType >
Teuchos::RCP< ROL::Vector< value_type > > Rol::VectorAdaptor< VectorType >::clone ( ) const

Return a clone of the wrapped vector.

◆ basis()

template<typename VectorType >
Teuchos::RCP< ROL::Vector< value_type > > Rol::VectorAdaptor< VectorType >::basis ( const int  i) const

Create and return a Teuchos smart reference counting pointer to the basis vector corresponding to the i \({}^{th}\) element of the wrapper vector.

◆ applyUnary()

template<typename VectorType >
void Rol::VectorAdaptor< VectorType >::applyUnary ( const ROL::Elementwise::UnaryFunction< value_type > &  f)

Apply unary function f to all the elements of the wrapped vector.

◆ applyBinary()

template<typename VectorType >
void Rol::VectorAdaptor< VectorType >::applyBinary ( const ROL::Elementwise::BinaryFunction< value_type > &  f,
const ROL::Vector< value_type > &  rol_vector 
)

Apply binary function f along with ROL::Vector rol_vector to all the elements of the wrapped vector.

◆ reduce()

template<typename VectorType >
value_type Rol::VectorAdaptor< VectorType >::reduce ( const ROL::Elementwise::ReductionOp< value_type > &  r) const

Return the accumulated value on applying reduction operation r on all the elements of the wrapped vector.

◆ print()

template<typename VectorType >
void Rol::VectorAdaptor< VectorType >::print ( std::ostream &  outStream) const

Print the wrapped vector to the output stream outStream.

Member Data Documentation

◆ vector_ptr

template<typename VectorType >
Teuchos::RCP<VectorType> Rol::VectorAdaptor< VectorType >::vector_ptr
private

Teuchos smart reference counting pointer to the underlying vector of type VectorType.

Definition at line 136 of file vector_adaptor.h.


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