SAMRAI::solv::PETScAbstractVectorReal< TYPE > Class Template Reference

#include <source/solvers/packages/petsc/PETScAbstractVectorReal.h>

Inheritance diagram for SAMRAI::solv::PETScAbstractVectorReal< TYPE >:

Inheritance graph
[legend]
List of all members.

Protected Member Functions

 PETScAbstractVectorReal (bool vector_created_via_duplicate, MPI_Comm comm)
virtual ~PETScAbstractVectorReal ()
Vec getPETScVector ()
virtual PETScAbstractVectorReal<
TYPE > * 
makeNewVector ()=0
virtual void freeVector ()=0
virtual void viewVector () const=0
virtual TYPE dotWith (const PETScAbstractVectorReal< TYPE > *y, bool local_only=false) const =0
virtual TYPE TdotWith (const PETScAbstractVectorReal< TYPE > *y, bool local_only=false) const =0
virtual TYPE L1Norm (bool local_only=false) const =0
virtual TYPE L2Norm (bool local_only=false) const =0
virtual TYPE maxNorm (bool local_only=false) const =0
virtual void scaleVector (const TYPE alpha)=0
virtual void copyVector (const PETScAbstractVectorReal< TYPE > *v_src)=0
virtual void setToScalar (const TYPE alpha)=0
virtual void swapWith (PETScAbstractVectorReal< TYPE > *v_other)=0
virtual void setAXPY (const TYPE alpha, const PETScAbstractVectorReal< TYPE > *x)=0
virtual void setAXPBY (const TYPE alpha, const PETScAbstractVectorReal< TYPE > *x, const TYPE beta)=0
virtual void setWAXPY (const TYPE alpha, const PETScAbstractVectorReal< TYPE > *x, const PETScAbstractVectorReal< TYPE > *y)=0
virtual void pointwiseMultiply (const PETScAbstractVectorReal< TYPE > *x, const PETScAbstractVectorReal< TYPE > *y)=0
virtual void pointwiseDivide (const PETScAbstractVectorReal< TYPE > *x, const PETScAbstractVectorReal< TYPE > *y)=0
virtual double maxPointwiseDivide (const PETScAbstractVectorReal< TYPE > *y)=0
virtual void vecMax (int &i, TYPE &max) const =0
virtual void vecMin (int &i, TYPE &min) const =0
virtual void setRandomValues (const TYPE width, const TYPE low)=0
virtual void getDataArray (TYPE **array)=0
virtual int getDataSize () const=0
virtual int getLocalDataSize () const=0
virtual void restoreDataArray (TYPE **array)=0

Detailed Description

template<class TYPE>
class SAMRAI::solv::PETScAbstractVectorReal< TYPE >

Class PETScAbstractVectorReal serves as an abstract base class for a C++ vector class that can be used with the PETSc solver framework. Specifically, this class provides an interface for real-valued PETSc vectors (i.e., where the data is either float or double). PETSc allows the use of user-defined vectors. Thus, the intent of this base class is that one may provide his/her own vector implementation in a subclass of this base class that provides the the necessary vector data structures and implements the pure virtual functions. This class declares private static member functions for linkage with the vector object function calls understood by PETSc. Each of the static member functions calls an associated function in a subclass via the virtual function mechanism. Note that the virtual members of this class are all protected. They should not be used outside of a subclass of this class. The data member of the PETSc vector object is set to an instantiation of the user-supplied vector class, when an object of this class is constructed.

PETSc was developed in the Mathematics and Computer Science Division at Argonne National Laboratory (ANL). For more information about PETSc, see http://www-fp.mcs.anl.gov/petsc/.

Important notes:


Constructor & Destructor Documentation

template<class TYPE>
SAMRAI::solv::PETScAbstractVectorReal< TYPE >::PETScAbstractVectorReal ( bool  vector_created_via_duplicate,
MPI_Comm  comm 
) [protected]

Constructor PETScAbstractVectorReal class that provides a wrapper for a SAMRAI vector so that it can be manipulated within PETSc. The constructor allocates the PETSc vector and sets its data structures and member functions so that it can operate on the SAMRAI vector.

template<class TYPE>
SAMRAI::solv::PETScAbstractVectorReal< TYPE >::~PETScAbstractVectorReal (  )  [protected, virtual]

Destructor for PETScAbstractVectorReal class destroys the PETSc vector created by the constructor.


Member Function Documentation

template<class TYPE>
Vec SAMRAI::solv::PETScAbstractVectorReal< TYPE >::getPETScVector (  )  [protected]

Return PETSc "Vec" object for this PETScAbstractVectorReal object.

template<class TYPE>
virtual PETScAbstractVectorReal<TYPE>* SAMRAI::solv::PETScAbstractVectorReal< TYPE >::makeNewVector (  )  [protected, pure virtual]

Clone the vector structure and allocate storage for the vector data. Then, return a pointer to the new vector instance. This function is distinct from the vector constructor since it is called from within PETSc (via the duplicateVec(), duplicateVecs() functions) to allocate new vector objects.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::freeVector (  )  [protected, pure virtual]

Destroy vector structure and its storage. This function is distinct from the destructor since it will be called from within PETSc to deallocate vectors (via the destroyVec(), destroyVecs() functions).

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::viewVector (  )  const [protected, pure virtual]

View all vector data. Note that the user-supplied vector must choose how to view the vector and its data; e.g., print to file, dump to standard out, etc.

template<class TYPE>
virtual TYPE SAMRAI::solv::PETScAbstractVectorReal< TYPE >::dotWith ( const PETScAbstractVectorReal< TYPE > *  y,
bool  local_only = false 
) const [protected, pure virtual]

Return $ (x,y) = \sum_i ( x_i * conj(y_i) ) $ , where $ x $ is this vector. Note that for real vectors, this is the same as TdotWith(). If local_only is true, the operation is limited to parts owned by the local process.

template<class TYPE>
virtual TYPE SAMRAI::solv::PETScAbstractVectorReal< TYPE >::TdotWith ( const PETScAbstractVectorReal< TYPE > *  y,
bool  local_only = false 
) const [protected, pure virtual]

Limited to local data only, return $ (x,y) = \sum_i ( x_i * (y_i) ) $ , where $ x $ is this vector. Note that for real vectors, this is the same as dotWith(). If local_only is true, the operation is limited to parts owned by the local process.

template<class TYPE>
virtual TYPE SAMRAI::solv::PETScAbstractVectorReal< TYPE >::L1Norm ( bool  local_only = false  )  const [protected, pure virtual]

Return $ L_1 $ -norm of this vector.

Parameters:
local_only Flag to get result for local data only.

template<class TYPE>
virtual TYPE SAMRAI::solv::PETScAbstractVectorReal< TYPE >::L2Norm ( bool  local_only = false  )  const [protected, pure virtual]

Return $ L_2 $ -norm of this vector.

Parameters:
local_only Flag to get result for local data only.

template<class TYPE>
virtual TYPE SAMRAI::solv::PETScAbstractVectorReal< TYPE >::maxNorm ( bool  local_only = false  )  const [protected, pure virtual]

Return $ L_{\infty} $ -norm of this vector.

Parameters:
local_only Flag to get result for local data only.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::scaleVector ( const TYPE  alpha  )  [protected, pure virtual]

Multiply each entry of this vector by given scalar.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::copyVector ( const PETScAbstractVectorReal< TYPE > *  v_src  )  [protected, pure virtual]

Copy source vector data to this vector.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::setToScalar ( const TYPE  alpha  )  [protected, pure virtual]

Set each entry of this vector to given scalar.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::swapWith ( PETScAbstractVectorReal< TYPE > *  v_other  )  [protected, pure virtual]

Swap data between this vector and argument vector.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::setAXPY ( const TYPE  alpha,
const PETScAbstractVectorReal< TYPE > *  x 
) [protected, pure virtual]

Set $ y = \alpha x + y $ , where $ y $ is this vector.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::setAXPBY ( const TYPE  alpha,
const PETScAbstractVectorReal< TYPE > *  x,
const TYPE  beta 
) [protected, pure virtual]

Set $ y = \alpha x + @beta y $ , where $ y $ is this vector.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::setWAXPY ( const TYPE  alpha,
const PETScAbstractVectorReal< TYPE > *  x,
const PETScAbstractVectorReal< TYPE > *  y 
) [protected, pure virtual]

Set $ w = \alpha x + y $ , where $ w $ is this vector.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::pointwiseMultiply ( const PETScAbstractVectorReal< TYPE > *  x,
const PETScAbstractVectorReal< TYPE > *  y 
) [protected, pure virtual]

Set $ w_i = x_i y_i $ , where $ w_i $ is $ i $ -th entry of this vector.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::pointwiseDivide ( const PETScAbstractVectorReal< TYPE > *  x,
const PETScAbstractVectorReal< TYPE > *  y 
) [protected, pure virtual]

Set $ w_i = x_i / y_i $ , where $ w_i $ is $ i $ -th entry of this vector.

template<class TYPE>
virtual double SAMRAI::solv::PETScAbstractVectorReal< TYPE >::maxPointwiseDivide ( const PETScAbstractVectorReal< TYPE > *  y  )  [protected, pure virtual]

Compute $ max_i = abs(w_i / y_i) $ , where $ w_i $ is $ i $ -th entry of this vector.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::vecMax ( int &  i,
TYPE &  max 
) const [protected, pure virtual]

Find maximum vector entry and vector index at which maximum occurs.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::vecMin ( int &  i,
TYPE &  min 
) const [protected, pure virtual]

Find minimum vector entry and vector index at which minimum occurs.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::setRandomValues ( const TYPE  width,
const TYPE  low 
) [protected, pure virtual]

Set vector entries to random values. Note that PETSc uses the drand48() function and computes random value as width*drand48()+low.

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::getDataArray ( TYPE **  array  )  [protected, pure virtual]

Set argument to vector data in contiguous array (local to processor).

template<class TYPE>
virtual int SAMRAI::solv::PETScAbstractVectorReal< TYPE >::getDataSize (  )  const [protected, pure virtual]

Return total length of vector.

template<class TYPE>
virtual int SAMRAI::solv::PETScAbstractVectorReal< TYPE >::getLocalDataSize (  )  const [protected, pure virtual]

Return length of vector (local to processor).

template<class TYPE>
virtual void SAMRAI::solv::PETScAbstractVectorReal< TYPE >::restoreDataArray ( TYPE **  array  )  [protected, pure virtual]

Restore pointer to vector data in contiguous array (local to processor).


The documentation for this class was generated from the following files:
Generated on Thu Jun 18 11:28:56 2009 for SAMRAI by  doxygen 1.5.1