ROL
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ROL::PinTVector< Real > Class Template Reference

Defines a parallel in time vector. More...

#include <ROL_PinTVector.hpp>

+ Inheritance diagram for ROL::PinTVector< Real >:

Public Types

enum  ESendRecv { SEND_AND_RECV, SEND_ONLY, RECV_ONLY }
 

Public Member Functions

 PinTVector ()
 
 PinTVector (const PinTVector &v)
 
 PinTVector (const Ptr< const PinTCommunicators > &comm, const Ptr< Vector< Real >> &localVector, int steps, const std::vector< int > &stencil)
 
virtual ~PinTVector ()
 
void initialize (const Ptr< const PinTCommunicators > &comm, const Ptr< Vector< Real >> &localVector, int steps, const std::vector< int > &stencil)
 
int numOwnedVectors () const
 
int numOwnedSteps () const
 
const std::vector< int > & stencil () const
 
const PinTCommunicatorscommunicators () const
 
bool isValidIndex (int i) const
 Determine if an index is valid including the stencil. More...
 
Ptr< Vector< Real > > getVectorPtr (int i) const
 
Ptr< Vector< Real > > getRemoteBufferPtr (int i) const
 
void boundaryExchange (ESendRecv send_recv=SEND_AND_RECV) const
 Exchange unknowns with neighboring processors. More...
 
void boundaryExchangeSumInto ()
 Exchange unknowns with neighboring processors using summation. More...
 
virtual void plus (const Vector< Real > &x)
 Compute \(y \leftarrow y + x\), where \(y = \mathtt{*this}\). More...
 
virtual 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...
 
virtual Real norm () const
 Returns \( \| y \| \) where \(y = \mathtt{*this}\). More...
 
virtual ROL::Ptr< Vector< Real > > clone () const
 Clone to make a new (uninitialized) vector. More...
 
virtual void print (std::ostream &outStream) const override
 
virtual void applyUnary (const Elementwise::UnaryFunction< Real > &f) override
 
virtual void set (const Vector< Real > &x) override
 Set \(y \leftarrow x\) where \(y = \mathtt{*this}\). More...
 
virtual void axpy (const Real alpha, const Vector< Real > &x) override
 Compute \(y \leftarrow \alpha x + y\) where \(y = \mathtt{*this}\). More...
 
virtual void zeroAll ()
 Zeros all entries of the vector, including boundary exchange fields. More...
 
- Public Member Functions inherited from ROL::Vector< Real >
virtual ~Vector ()
 
virtual void plus (const Vector &x)=0
 Compute \(y \leftarrow y + x\), where \(y = \mathtt{*this}\). More...
 
virtual Real dot (const Vector &x) const =0
 Compute \( \langle y,x \rangle \) where \(y = \mathtt{*this}\). More...
 
virtual void axpy (const Real alpha, const Vector &x)
 Compute \(y \leftarrow \alpha x + y\) where \(y = \mathtt{*this}\). More...
 
virtual void zero ()
 Set to zero vector. More...
 
virtual ROL::Ptr< Vectorbasis (const int i) const
 Return i-th basis vector. More...
 
virtual int dimension () const
 Return dimension of the vector space. More...
 
virtual void set (const Vector &x)
 Set \(y \leftarrow x\) where \(y = \mathtt{*this}\). 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 void applyBinary (const Elementwise::BinaryFunction< Real > &f, const Vector &x)
 
virtual Real reduce (const Elementwise::ReductionOp< Real > &r) const
 
virtual void setScalar (const Real C)
 Set \(y \leftarrow C\) where \(C\in\mathbb{R}\). More...
 
virtual std::vector< RealcheckVector (const Vector< Real > &x, const Vector< Real > &y, const bool printToStream=true, std::ostream &outStream=std::cout) const
 Verify vector-space methods. More...
 

Protected Member Functions

void computeStepStartEnd (int steps)
 
void allocateBoundaryExchangeVectors ()
 

Protected Attributes

bool isInitialized_
 
Ptr< const PinTCommunicatorscommunicators_
 
Ptr< Vector< Real > > localVector_
 
int steps_
 
std::vector< int > stencil_
 
int stepStart_
 
int stepEnd_
 
Ptr< PinTVector< Real > > dual_
 
Ptr< PartitionedVector< Real > > stepVectors_
 
std::vector< Ptr< Vector< Real > > > leftVectors_
 
std::vector< Ptr< Vector< Real > > > rightVectors_
 
Ptr< PinTVectorCommunication
< Real > > 
vectorComm_
 

Detailed Description

template<class Real>
class ROL::PinTVector< Real >

Defines a parallel in time vector.

Definition at line 169 of file ROL_PinTVector.hpp.

Member Enumeration Documentation

template<class Real>
enum ROL::PinTVector::ESendRecv
Enumerator
SEND_AND_RECV 
SEND_ONLY 
RECV_ONLY 

Definition at line 258 of file ROL_PinTVector.hpp.

Constructor & Destructor Documentation

template<class Real>
ROL::PinTVector< Real >::PinTVector ( )
inline
template<class Real>
ROL::PinTVector< Real >::PinTVector ( const PinTVector< Real > &  v)
inline
template<class Real>
ROL::PinTVector< Real >::PinTVector ( const Ptr< const PinTCommunicators > &  comm,
const Ptr< Vector< Real >> &  localVector,
int  steps,
const std::vector< int > &  stencil 
)
inline

Definition at line 277 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::initialize().

template<class Real>
virtual ROL::PinTVector< Real >::~PinTVector ( )
inlinevirtual

Definition at line 286 of file ROL_PinTVector.hpp.

Member Function Documentation

template<class Real>
void ROL::PinTVector< Real >::computeStepStartEnd ( int  steps)
inlineprotected
template<class Real>
void ROL::PinTVector< Real >::allocateBoundaryExchangeVectors ( )
inlineprotected
template<class Real>
void ROL::PinTVector< Real >::initialize ( const Ptr< const PinTCommunicators > &  comm,
const Ptr< Vector< Real >> &  localVector,
int  steps,
const std::vector< int > &  stencil 
)
inline
template<class Real>
int ROL::PinTVector< Real >::numOwnedVectors ( ) const
inline

How many vectors are "owned" by this processor. This also includes any duplicate vectors.

It may be the case that this processor belongs to a sub group of processors. All processors in that sub group will return the same number of owned vectors.

Definition at line 323 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::stepVectors_.

template<class Real>
int ROL::PinTVector< Real >::numOwnedSteps ( ) const
inline
template<class Real>
const std::vector<int>& ROL::PinTVector< Real >::stencil ( ) const
inline

What is the stencil used to build this vector?

This accessor is directly based on what the user intiializes the vector with. Its used by ROL algorithms and constraints to ensure the communication patterns are understood.

Definition at line 345 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::stencil_.

Referenced by ROL::Constraint_PinTSimOpt< Real >::applyAdjointJacobian_1(), ROL::Constraint_PinTSimOpt< Real >::applyJacobian_1(), ROL::Constraint_PinTSimOpt< Real >::applyJacobian_2(), ROL::Constraint_PinTSimOpt< Real >::getNonconstVector(), ROL::PinTVector< Real >::initialize(), ROL::Constraint_PinTSimOpt< Real >::solve(), and ROL::Constraint_PinTSimOpt< Real >::value().

template<class Real>
const PinTCommunicators& ROL::PinTVector< Real >::communicators ( ) const
inline

What is the communicators object used to build this vector?

Definition at line 349 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::communicators_.

Referenced by ROL::Constraint_PinTSimOpt< Real >::applyJacobian_1().

template<class Real>
bool ROL::PinTVector< Real >::isValidIndex ( int  i) const
inline

Determine if an index is valid including the stencil.

An index is valid if is in [min(stencil),max(stencil)+numOwnedSteps()-1] Note that this treats owned vectors that live in the "shared" region as negative numbers, or numbers greater than numOwnedSteps()-1

Definition at line 358 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::leftVectors_, ROL::PinTVector< Real >::numOwnedSteps(), and ROL::PinTVector< Real >::rightVectors_.

Referenced by ROL::PinTVector< Real >::getRemoteBufferPtr(), and ROL::PinTVector< Real >::getVectorPtr().

template<class Real>
Ptr<Vector<Real> > ROL::PinTVector< Real >::getVectorPtr ( int  i) const
inline
template<class Real>
Ptr<Vector<Real> > ROL::PinTVector< Real >::getRemoteBufferPtr ( int  i) const
inline
template<class Real>
void ROL::PinTVector< Real >::boundaryExchange ( ESendRecv  send_recv = SEND_AND_RECV) const
inline
template<class Real>
void ROL::PinTVector< Real >::boundaryExchangeSumInto ( )
inline

Exchange unknowns with neighboring processors using summation.

Using the stencil information, do global communication with time neighbor processors. This sums in the the reverse direction.

Note
This method is not const because it does change the behavior of the vector.

Definition at line 479 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::communicators_, ROL::PinTVector< Real >::getRemoteBufferPtr(), ROL::PinTVector< Real >::getVectorPtr(), ROL::PinTVector< Real >::numOwnedSteps(), ROL::PinTVector< Real >::stencil_, ROL::PinTVector< Real >::stepEnd_, ROL::PinTVector< Real >::steps_, ROL::PinTVector< Real >::stepStart_, and ROL::PinTVector< Real >::vectorComm_.

Referenced by ROL::Constraint_PinTSimOpt< Real >::applyAdjointJacobian_1().

template<class Real>
virtual void ROL::PinTVector< Real >::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\).


Definition at line 529 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::PinTVector(), and ROL::PinTVector< Real >::stepVectors_.

template<class Real>
virtual void ROL::PinTVector< Real >::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 546 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::stepVectors_.

template<class Real>
virtual Real ROL::PinTVector< Real >::dot ( const Vector< Real > &  x) const
inlinevirtual

Compute \( \langle y,x \rangle \) where \(y = \mathtt{*this}\).

Parameters
[in]xis the vector that forms the dot product with \(\mathtt{*this}\).
Returns
The number equal to \(\langle \mathtt{*this}, x \rangle\).

Definition at line 559 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::communicators_, ROL::PinTVector< Real >::PinTVector(), and ROL::PinTVector< Real >::stepVectors_.

Referenced by ROL::PinTVector< Real >::norm().

template<class Real>
virtual Real ROL::PinTVector< Real >::norm ( void  ) 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 584 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::dot().

template<class Real>
virtual ROL::Ptr<Vector<Real> > ROL::PinTVector< Real >::clone ( void  ) const
inlinevirtual

Clone to make a new (uninitialized) vector.

Returns
A reference-counted pointer to the cloned vector.

Provides the means of allocating temporary memory in ROL.


Implements ROL::Vector< Real >.

Definition at line 598 of file ROL_PinTVector.hpp.

template<class Real>
virtual void ROL::PinTVector< Real >::print ( std::ostream &  outStream) const
inlineoverridevirtual

Reimplemented from ROL::Vector< Real >.

Definition at line 603 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::stepVectors_.

template<class Real>
virtual void ROL::PinTVector< Real >::applyUnary ( const Elementwise::UnaryFunction< Real > &  f)
inlineoverridevirtual

Reimplemented from ROL::Vector< Real >.

Definition at line 608 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::stepVectors_.

template<class Real>
virtual void ROL::PinTVector< Real >::set ( const Vector< Real > &  x)
inlineoverridevirtual

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.


Definition at line 622 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::PinTVector(), and ROL::PinTVector< Real >::stepVectors_.

Referenced by ROL::PinTVector< Real >::PinTVector().

template<class Real>
virtual void ROL::PinTVector< Real >::axpy ( const Real  alpha,
const Vector< Real > &  x 
)
inlineoverridevirtual

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.


Definition at line 640 of file ROL_PinTVector.hpp.

References ROL::PinTVector< Real >::PinTVector(), and ROL::PinTVector< Real >::stepVectors_.

template<class Real>
virtual void ROL::PinTVector< Real >::zeroAll ( )
inlinevirtual

Member Data Documentation

template<class Real>
bool ROL::PinTVector< Real >::isInitialized_
protected

Definition at line 176 of file ROL_PinTVector.hpp.

Referenced by ROL::PinTVector< Real >::initialize().

template<class Real>
Ptr<const PinTCommunicators> ROL::PinTVector< Real >::communicators_
protected
template<class Real>
Ptr<Vector<Real> > ROL::PinTVector< Real >::localVector_
protected
template<class Real>
int ROL::PinTVector< Real >::steps_
protected
template<class Real>
std::vector<int> ROL::PinTVector< Real >::stencil_
protected
template<class Real>
int ROL::PinTVector< Real >::stepStart_
protected
template<class Real>
int ROL::PinTVector< Real >::stepEnd_
protected
template<class Real>
Ptr<PinTVector<Real> > ROL::PinTVector< Real >::dual_
mutableprotected

Definition at line 188 of file ROL_PinTVector.hpp.

template<class Real>
Ptr<PartitionedVector<Real> > ROL::PinTVector< Real >::stepVectors_
protected
template<class Real>
std::vector<Ptr<Vector<Real> > > ROL::PinTVector< Real >::leftVectors_
protected
template<class Real>
std::vector<Ptr<Vector<Real> > > ROL::PinTVector< Real >::rightVectors_
protected
template<class Real>
Ptr<PinTVectorCommunication<Real> > ROL::PinTVector< Real >::vectorComm_
protected

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