Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
Epetra_MultiVector Class Reference

Epetra_MultiVector: A class for constructing and using dense multi-vectors, vectors and matrices in parallel. More...

#include <Epetra_MultiVector.h>

Inheritance diagram for Epetra_MultiVector:
Inheritance graph
[legend]

Public Member Functions

int ReplaceMap (const Epetra_BlockMap &map)
 Replace map, only if new map has same point-structure as current map. More...
 
int Reduce ()
 
- Public Member Functions inherited from Epetra_DistObject
 Epetra_DistObject (const Epetra_BlockMap &Map)
 Basic Epetra_DistObject constuctor. More...
 
 Epetra_DistObject (const Epetra_BlockMap &Map, const char *const Label)
 
 Epetra_DistObject (const Epetra_DistObject &Source)
 Epetra_DistObject copy constructor. More...
 
virtual ~Epetra_DistObject ()
 Epetra_DistObject destructor. More...
 
int Import (const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
 Imports an Epetra_DistObject using the Epetra_Import object. More...
 
int Import (const Epetra_SrcDistObject &A, const Epetra_Export &Exporter, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
 Imports an Epetra_DistObject using the Epetra_Export object. More...
 
int Export (const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
 Exports an Epetra_DistObject using the Epetra_Import object. More...
 
int Export (const Epetra_SrcDistObject &A, const Epetra_Export &Exporter, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
 Exports an Epetra_DistObject using the Epetra_Export object. More...
 
const Epetra_BlockMapMap () const
 Returns the address of the Epetra_BlockMap for this multi-vector. More...
 
const Epetra_CommComm () const
 Returns the address of the Epetra_Comm for this multi-vector. More...
 
bool DistributedGlobal () const
 Returns true if this multi-vector is distributed global, i.e., not local replicated. More...
 
- Public Member Functions inherited from Epetra_Object
 Epetra_Object (int TracebackModeIn=-1, bool set_label=true)
 Epetra_Object Constructor. More...
 
 Epetra_Object (const char *const Label, int TracebackModeIn=-1)
 Epetra_Object Constructor. More...
 
 Epetra_Object (const Epetra_Object &Object)
 Epetra_Object Copy Constructor. More...
 
virtual ~Epetra_Object ()
 Epetra_Object Destructor. More...
 
virtual void SetLabel (const char *const Label)
 Epetra_Object Label definition using char *. More...
 
virtual const char * Label () const
 Epetra_Object Label access funtion. More...
 
virtual int ReportError (const std::string Message, int ErrorCode) const
 Error reporting method. More...
 
- Public Member Functions inherited from Epetra_SrcDistObject
virtual ~Epetra_SrcDistObject ()
 Epetra_SrcDistObject destructor. More...
 
- Public Member Functions inherited from Epetra_CompObject
Epetra_CompObjectoperator= (const Epetra_CompObject &src)
 
 Epetra_CompObject ()
 Basic Epetra_CompObject constuctor. More...
 
 Epetra_CompObject (const Epetra_CompObject &Source)
 Epetra_CompObject copy constructor. More...
 
virtual ~Epetra_CompObject ()
 Epetra_CompObject destructor. More...
 
void SetFlopCounter (const Epetra_Flops &FlopCounter_in)
 Set the internal Epetra_Flops() pointer. More...
 
void SetFlopCounter (const Epetra_CompObject &CompObject)
 Set the internal Epetra_Flops() pointer to the flop counter of another Epetra_CompObject. More...
 
void UnsetFlopCounter ()
 Set the internal Epetra_Flops() pointer to 0 (no flops counted). More...
 
Epetra_FlopsGetFlopCounter () const
 Get the pointer to the Epetra_Flops() object associated with this object, returns 0 if none. More...
 
void ResetFlops () const
 Resets the number of floating point operations to zero for this multi-vector. More...
 
double Flops () const
 Returns the number of floating point operations with this multi-vector. More...
 
void UpdateFlops (int Flops_in) const
 Increment Flop count for this object. More...
 
void UpdateFlops (long int Flops_in) const
 Increment Flop count for this object. More...
 
void UpdateFlops (long long Flops_in) const
 Increment Flop count for this object. More...
 
void UpdateFlops (double Flops_in) const
 Increment Flop count for this object. More...
 
void UpdateFlops (float Flops_in) const
 Increment Flop count for this object. More...
 
- Public Member Functions inherited from Epetra_BLAS
 Epetra_BLAS (void)
 Epetra_BLAS Constructor. More...
 
 Epetra_BLAS (const Epetra_BLAS &BLAS)
 Epetra_BLAS Copy Constructor. More...
 
virtual ~Epetra_BLAS (void)
 Epetra_BLAS Destructor. More...
 
float ASUM (const int N, const float *X, const int INCX=1) const
 Epetra_BLAS one norm function (SASUM). More...
 
double ASUM (const int N, const double *X, const int INCX=1) const
 Epetra_BLAS one norm function (DASUM). More...
 
float DOT (const int N, const float *X, const float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS dot product function (SDOT). More...
 
double DOT (const int N, const double *X, const double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS dot product function (DDOT). More...
 
float NRM2 (const int N, const float *X, const int INCX=1) const
 Epetra_BLAS norm function (SNRM2). More...
 
double NRM2 (const int N, const double *X, const int INCX=1) const
 Epetra_BLAS norm function (DNRM2). More...
 
void SCAL (const int N, const float ALPHA, float *X, const int INCX=1) const
 Epetra_BLAS vector scale function (SSCAL) More...
 
void SCAL (const int N, const double ALPHA, double *X, const int INCX=1) const
 Epetra_BLAS vector scale function (DSCAL) More...
 
void COPY (const int N, const float *X, float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector copy function (SCOPY) More...
 
void COPY (const int N, const double *X, double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector scale function (DCOPY) More...
 
int IAMAX (const int N, const float *X, const int INCX=1) const
 Epetra_BLAS arg maximum of absolute value function (ISAMAX) More...
 
int IAMAX (const int N, const double *X, const int INCX=1) const
 Epetra_BLAS arg maximum of absolute value function (IDAMAX) More...
 
void AXPY (const int N, const float ALPHA, const float *X, float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector update function (SAXPY) More...
 
void AXPY (const int N, const double ALPHA, const double *X, double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector update function (DAXPY) More...
 
void GEMV (const char TRANS, const int M, const int N, const float ALPHA, const float *A, const int LDA, const float *X, const float BETA, float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS matrix-vector multiply function (SGEMV) More...
 
void GEMV (const char TRANS, const int M, const int N, const double ALPHA, const double *A, const int LDA, const double *X, const double BETA, double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS matrix-vector multiply function (DGEMV) More...
 
void GEMM (const char TRANSA, const char TRANSB, const int M, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const
 Epetra_BLAS matrix-matrix multiply function (SGEMM) More...
 
void GEMM (const char TRANSA, const char TRANSB, const int M, const int N, const int K, const double ALPHA, const double *A, const int LDA, const double *B, const int LDB, const double BETA, double *C, const int LDC) const
 Epetra_BLAS matrix-matrix multiply function (DGEMM) More...
 
void SYMM (const char SIDE, const char UPLO, const int M, const int N, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const
 Epetra_BLAS symmetric matrix-matrix multiply function (SSYMM) More...
 
void SYMM (const char SIDE, const char UPLO, const int M, const int N, const double ALPHA, const double *A, const int LDA, const double *B, const int LDB, const double BETA, double *C, const int LDC) const
 Epetra_BLAS matrix-matrix multiply function (DSYMM) More...
 
void TRMM (const char SIDE, const char UPLO, const char TRANSA, const char DIAG, const int M, const int N, const float ALPHA, const float *A, const int LDA, float *B, const int LDB) const
 Epetra_BLAS triangular matrix-matrix multiply function (STRMM) More...
 
void TRMM (const char SIDE, const char UPLO, const char TRANSA, const char DIAG, const int M, const int N, const double ALPHA, const double *A, const int LDA, double *B, const int LDB) const
 Epetra_BLAS triangular matrix-matrix multiply function (DTRMM) More...
 
void SYRK (const char UPLO, const char TRANS, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float BETA, float *C, const int LDC) const
 Eperta_BLAS symetric rank k funtion (ssyrk) More...
 
void SYRK (const char UPLO, const char TRANS, const int N, const int K, const double ALPHA, const double *A, const int LDA, const double BETA, double *C, const int LDC) const
 Eperta_BLAS symetric rank k funtion (dsyrk) More...
 

Protected Member Functions

void Assign (const Epetra_MultiVector &rhs)
 
int CheckInput ()
 
- Protected Member Functions inherited from Epetra_DistObject
virtual int DoTransfer (const Epetra_SrcDistObject &A, Epetra_CombineMode CombineMode, int NumSameIDs, int NumPermuteIDs, int NumRemoteIDs, int NumExportIDs, int *PermuteToLIDs, int *PermuteFromLIDs, int *RemoteLIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &LenImports, char *&Imports, Epetra_Distributor &Distor, bool DoReverse, const Epetra_OffsetIndex *Indexor)
 Perform actual transfer (redistribution) of data across memory images, using Epetra_Distributor object. More...
 
- Protected Member Functions inherited from Epetra_Object
std::string toString (const int &x) const
 
std::string toString (const long long &x) const
 
std::string toString (const double &x) const
 

Protected Attributes

double * Values_
 
- Protected Attributes inherited from Epetra_DistObject
Epetra_BlockMap Map_
 
const Epetra_CommComm_
 
char * Exports_
 
char * Imports_
 
int LenExports_
 
int LenImports_
 
int * Sizes_
 
- Protected Attributes inherited from Epetra_CompObject
Epetra_FlopsFlopCounter_
 

Private Member Functions

int AllocateForCopy (void)
 
int DoCopy (void)
 
void UpdateDoubleTemp () const
 
void UpdateVectors () const
 
int AllocateForView (void)
 
int DoView (void)
 
template<typename int_type >
int ChangeGlobalValue (int_type GlobalBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue, bool SumInto)
 
int ChangeMyValue (int MyBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue, bool SumInto)
 
int CheckSizes (const Epetra_SrcDistObject &A)
 Allows the source and target (this) objects to be compared for compatibility, return nonzero if not. More...
 
int CopyAndPermute (const Epetra_SrcDistObject &Source, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode=Zero)
 Perform ID copies and permutations that are on processor. More...
 
int PackAndPrepare (const Epetra_SrcDistObject &Source, int NumExportIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &SizeOfPacket, int *Sizes, bool &VarSizes, Epetra_Distributor &Distor)
 Perform any packing or preparation required for call to DoTransfer(). More...
 
int UnpackAndCombine (const Epetra_SrcDistObject &Source, int NumImportIDs, int *ImportLIDs, int LenImports, char *Imports, int &SizeOfPacket, Epetra_Distributor &Distor, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor)
 Perform any unpacking and combining after call to DoTransfer(). More...
 

Private Attributes

double ** Pointers_
 
int MyLength_
 
long long GlobalLength_
 
int NumVectors_
 
bool UserAllocated_
 
bool ConstantStride_
 
int Stride_
 
bool Allocated_
 
double * DoubleTemp_
 
Epetra_Vector ** Vectors_
 
Epetra_Util Util_
 

Constructors/destructors

 Epetra_MultiVector (const Epetra_BlockMap &Map, int NumVectors, bool zeroOut=true)
 Basic Epetra_MultiVector constuctor. More...
 
 Epetra_MultiVector (const Epetra_MultiVector &Source)
 Epetra_MultiVector copy constructor. More...
 
 Epetra_MultiVector (Epetra_DataAccess CV, const Epetra_BlockMap &Map, double *A, int MyLDA, int NumVectors)
 Set multi-vector values from two-dimensional array. More...
 
 Epetra_MultiVector (Epetra_DataAccess CV, const Epetra_BlockMap &Map, double **ArrayOfPointers, int NumVectors)
 Set multi-vector values from array of pointers. More...
 
 Epetra_MultiVector (Epetra_DataAccess CV, const Epetra_MultiVector &Source, int *Indices, int NumVectors)
 Set multi-vector values from list of vectors in an existing Epetra_MultiVector. More...
 
 Epetra_MultiVector (Epetra_DataAccess CV, const Epetra_MultiVector &Source, int StartIndex, int NumVectors)
 Set multi-vector values from range of vectors in an existing Epetra_MultiVector. More...
 
virtual ~Epetra_MultiVector ()
 Epetra_MultiVector destructor. More...
 

Post-construction modification routines

int ReplaceGlobalValue (int GlobalRow, int VectorIndex, double ScalarValue)
 Replace current value at the specified (GlobalRow, VectorIndex) location with ScalarValue. More...
 
int ReplaceGlobalValue (long long GlobalRow, int VectorIndex, double ScalarValue)
 
int ReplaceGlobalValue (int GlobalBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Replace current value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex) location with ScalarValue. More...
 
int ReplaceGlobalValue (long long GlobalBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 
int SumIntoGlobalValue (int GlobalRow, int VectorIndex, double ScalarValue)
 Adds ScalarValue to existing value at the specified (GlobalRow, VectorIndex) location. More...
 
int SumIntoGlobalValue (long long GlobalRow, int VectorIndex, double ScalarValue)
 
int SumIntoGlobalValue (int GlobalBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Adds ScalarValue to existing value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex) location. More...
 
int SumIntoGlobalValue (long long GlobalBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 
int ReplaceMyValue (int MyRow, int VectorIndex, double ScalarValue)
 Replace current value at the specified (MyRow, VectorIndex) location with ScalarValue. More...
 
int ReplaceMyValue (int MyBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Replace current value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location with ScalarValue. More...
 
int SumIntoMyValue (int MyRow, int VectorIndex, double ScalarValue)
 Adds ScalarValue to existing value at the specified (MyRow, VectorIndex) location. More...
 
int SumIntoMyValue (int MyBlockRow, int BlockRowOffset, int VectorIndex, double ScalarValue)
 Adds ScalarValue to existing value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location. More...
 
int PutScalar (double ScalarConstant)
 Initialize all values in a multi-vector with constant value. More...
 
int Random ()
 Set multi-vector values to random numbers. More...
 

Extraction methods

int ExtractCopy (double *A, int MyLDA) const
 Put multi-vector values into user-provided two-dimensional array. More...
 
int ExtractCopy (double **ArrayOfPointers) const
 Put multi-vector values into user-provided array of pointers. More...
 
int ExtractView (double **A, int *MyLDA) const
 Set user-provided addresses of A and MyLDA. More...
 
int ExtractView (double ***ArrayOfPointers) const
 Set user-provided addresses of ArrayOfPointers. More...
 

Mathematical methods

int Dot (const Epetra_MultiVector &A, double *Result) const
 Computes dot product of each corresponding pair of vectors. More...
 
int Abs (const Epetra_MultiVector &A)
 Puts element-wise absolute values of input Multi-vector in target. More...
 
int Reciprocal (const Epetra_MultiVector &A)
 Puts element-wise reciprocal values of input Multi-vector in target. More...
 
int Scale (double ScalarValue)
 Scale the current values of a multi-vector, this = ScalarValue*this. More...
 
int Scale (double ScalarA, const Epetra_MultiVector &A)
 Replace multi-vector values with scaled values of A, this = ScalarA*A. More...
 
int Update (double ScalarA, const Epetra_MultiVector &A, double ScalarThis)
 Update multi-vector values with scaled values of A, this = ScalarThis*this + ScalarA*A. More...
 
int Update (double ScalarA, const Epetra_MultiVector &A, double ScalarB, const Epetra_MultiVector &B, double ScalarThis)
 Update multi-vector with scaled values of A and B, this = ScalarThis*this + ScalarA*A + ScalarB*B. More...
 
int Norm1 (double *Result) const
 Compute 1-norm of each vector in multi-vector. More...
 
int Norm2 (double *Result) const
 Compute 2-norm of each vector in multi-vector. More...
 
int NormInf (double *Result) const
 Compute Inf-norm of each vector in multi-vector. More...
 
int NormWeighted (const Epetra_MultiVector &Weights, double *Result) const
 Compute Weighted 2-norm (RMS Norm) of each vector in multi-vector. More...
 
int MinValue (double *Result) const
 Compute minimum value of each vector in multi-vector. More...
 
int MaxValue (double *Result) const
 Compute maximum value of each vector in multi-vector. More...
 
int MeanValue (double *Result) const
 Compute mean (average) value of each vector in multi-vector. More...
 
int Multiply (char TransA, char TransB, double ScalarAB, const Epetra_MultiVector &A, const Epetra_MultiVector &B, double ScalarThis)
 Matrix-Matrix multiplication, this = ScalarThis*this + ScalarAB*A*B. More...
 
int Multiply (double ScalarAB, const Epetra_MultiVector &A, const Epetra_MultiVector &B, double ScalarThis)
 Multiply a Epetra_MultiVector with another, element-by-element. More...
 
int ReciprocalMultiply (double ScalarAB, const Epetra_MultiVector &A, const Epetra_MultiVector &B, double ScalarThis)
 Multiply a Epetra_MultiVector by the reciprocal of another, element-by-element. More...
 

Random number utilities

int SetSeed (unsigned int Seed_in)
 Set seed for Random function. More...
 
unsigned int Seed ()
 Get seed from Random function. More...
 

Overloaded operators

Epetra_MultiVectoroperator= (const Epetra_MultiVector &Source)
 = Operator. More...
 
double *& operator[] (int i)
 Vector access function. More...
 
double *const & operator[] (int i) const
 Vector access function. More...
 
Epetra_Vector *& operator() (int i)
 Vector access function. More...
 
const Epetra_Vector *& operator() (int i) const
 Vector access function. More...
 

Attribute access functions

int NumVectors () const
 Returns the number of vectors in the multi-vector. More...
 
int MyLength () const
 Returns the local vector length on the calling processor of vectors in the multi-vector. More...
 
int GlobalLength () const
 Returns the global vector length of vectors in the multi-vector. More...
 
long long GlobalLength64 () const
 Returns the 64-bit global vector length of vectors in the multi-vector. More...
 
int Stride () const
 Returns the stride between vectors in the multi-vector (only meaningful if ConstantStride() is true). More...
 
bool ConstantStride () const
 Returns true if this multi-vector has constant stride between vectors. More...
 

I/O methods

virtual void Print (std::ostream &os) const
 Print method. More...
 

Expert-only unsupported methods

int ResetView (double **ArrayOfPointers)
 Reset the view of an existing multivector to point to new user data. More...
 
double * Values () const
 Get pointer to MultiVector values. More...
 
double ** Pointers () const
 Get pointer to individual vector pointers. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from Epetra_Object
static void SetTracebackMode (int TracebackModeValue)
 Set the value of the Epetra_Object error traceback report mode. More...
 
static int GetTracebackMode ()
 Get the value of the Epetra_Object error report mode. More...
 
static std::ostream & GetTracebackStream ()
 Get the output stream for error reporting. More...
 
- Static Public Attributes inherited from Epetra_Object
static int TracebackMode
 

Detailed Description

Epetra_MultiVector: A class for constructing and using dense multi-vectors, vectors and matrices in parallel.

The Epetra_MultiVector class enables the construction and use of real-valued, double-precision dense vectors, multi-vectors, and matrices in a distributed memory environment. The dimensions and distribution of the dense multi-vectors is determined in part by a Epetra_Comm object, a Epetra_Map (or Epetra_LocalMap or Epetra_BlockMap) and the number of vectors passed to the constructors described below.

There are several concepts that important for understanding the Epetra_MultiVector class:

Constructing Epetra_MultiVectors

Except for the basic constructor and copy constructor, Epetra_MultiVector constructors have two data access modes:

  1. Copy mode - Allocates memory and makes a copy of the user-provided data. In this case, the user data is not needed after construction.
  2. View mode - Creates a "view" of the user data. In this case, the user data is required to remain intact for the life of the multi-vector.
Warning
View mode is extremely dangerous from a data hiding perspective. Therefore, we strongly encourage users to develop code using Copy mode first and only use the View mode in a secondary optimization phase.

All Epetra_MultiVector constructors require a map argument that describes the layout of elements on the parallel machine. Specifically, map is a Epetra_Map, Epetra_LocalMap or Epetra_BlockMap object describing the desired memory layout for the multi-vector.

There are six different Epetra_MultiVector constructors:

Extracting Data from Epetra_MultiVectors

Once a Epetra_MultiVector is constructed, it is possible to extract a copy of the values or create a view of them.

Warning
ExtractView functions are extremely dangerous from a data hiding perspective. For both ExtractView fuctions, there is a corresponding ExtractCopy function. We strongly encourage users to develop code using ExtractCopy functions first and only use the ExtractView functions in a secondary optimization phase.

There are four Extract functions:

Vector, Matrix and Utility Functions

Once a Epetra_MultiVector is constructed, a variety of mathematical functions can be applied to the individual vectors. Specifically:

In addition, a matrix-matrix multiply function supports a variety of operations on any viable combination of global distributed and local replicated multi-vectors using calls to DGEMM, a high performance kernel for matrix operations. In the near future we will add support for calls to other selected BLAS and LAPACK functions.

Counting Floating Point Operations

Each Epetra_MultiVector object keep track of the number of serial floating point operations performed using the specified object as the this argument to the function. The Flops() function returns this number as a double precision number. Using this information, in conjunction with the Epetra_Time class, one can get accurate parallel performance numbers. The ResetFlops() function resets the floating point counter.

Warning
A Epetra_Map, Epetra_LocalMap or Epetra_BlockMap object is required for all Epetra_MultiVector constructors.

Definition at line 184 of file Epetra_MultiVector.h.

Constructor & Destructor Documentation

Epetra_MultiVector::Epetra_MultiVector ( const Epetra_BlockMap Map,
int  NumVectors,
bool  zeroOut = true 
)

Basic Epetra_MultiVector constuctor.

Creates a Epetra_MultiVector object and, by default, fills with zero values.
Parameters
InMap - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
Warning
Note that, because Epetra_LocalMap derives from Epetra_Map and Epetra_Map derives from Epetra_BlockMap, this constructor works for all three types of Epetra map classes.
Parameters
InNumVectors - Number of vectors in multi-vector.
InzeroOut - If true then the allocated memory will be zeroed out initialy. If false then this memory will not be touched which can be significantly faster.
Returns
Pointer to a Epetra_MultiVector.

Definition at line 60 of file Epetra_MultiVector.cpp.

Epetra_MultiVector::Epetra_MultiVector ( const Epetra_MultiVector Source)

Epetra_MultiVector copy constructor.

Definition at line 85 of file Epetra_MultiVector.cpp.

Epetra_MultiVector::Epetra_MultiVector ( Epetra_DataAccess  CV,
const Epetra_BlockMap Map,
double *  A,
int  MyLDA,
int  NumVectors 
)

Set multi-vector values from two-dimensional array.

Parameters
InEpetra_DataAccess - Enumerated type set to Copy or View.
InMap - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
InA - Pointer to an array of double precision numbers. The first vector starts at A. The second vector starts at A+MyLDA, the third at A+2*MyLDA, and so on.
InMyLDA - The "Leading Dimension", or stride between vectors in memory.
Warning
This value refers to the stride on the calling processor. Thus it is a local quantity, not a global quantity.
Parameters
InNumVectors - Number of vectors in multi-vector.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Definition at line 111 of file Epetra_MultiVector.cpp.

Epetra_MultiVector::Epetra_MultiVector ( Epetra_DataAccess  CV,
const Epetra_BlockMap Map,
double **  ArrayOfPointers,
int  NumVectors 
)

Set multi-vector values from array of pointers.

Parameters
InEpetra_DataAccess - Enumerated type set to Copy or View.
InMap - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
InArrayOfPointers - An array of pointers such that ArrayOfPointers[i] points to the memory location containing ith vector to be copied.
InNumVectors - Number of vectors in multi-vector.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Definition at line 141 of file Epetra_MultiVector.cpp.

Epetra_MultiVector::Epetra_MultiVector ( Epetra_DataAccess  CV,
const Epetra_MultiVector Source,
int *  Indices,
int  NumVectors 
)

Set multi-vector values from list of vectors in an existing Epetra_MultiVector.

Parameters
InEpetra_DataAccess - Enumerated type set to Copy or View.
InSource - An existing fully constructed Epetra_MultiVector.
InIndices - Integer list of the vectors to copy.
InNumVectors - Number of vectors in multi-vector.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Definition at line 172 of file Epetra_MultiVector.cpp.

Epetra_MultiVector::Epetra_MultiVector ( Epetra_DataAccess  CV,
const Epetra_MultiVector Source,
int  StartIndex,
int  NumVectors 
)

Set multi-vector values from range of vectors in an existing Epetra_MultiVector.

Parameters
InEpetra_DataAccess - Enumerated type set to Copy or View.
InSource - An existing fully constructed Epetra_MultiVector.
InStartIndex - First of the vectors to copy.
InNumVectors - Number of vectors in multi-vector.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Definition at line 203 of file Epetra_MultiVector.cpp.

Epetra_MultiVector::~Epetra_MultiVector ( )
virtual

Epetra_MultiVector destructor.

Definition at line 230 of file Epetra_MultiVector.cpp.

Member Function Documentation

int Epetra_MultiVector::ReplaceGlobalValue ( int  GlobalRow,
int  VectorIndex,
double  ScalarValue 
)

Replace current value at the specified (GlobalRow, VectorIndex) location with ScalarValue.

Replaces the existing value for a single entry in the multivector. The specified global row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the global row will be modified. To modify a different point entry, use the other version of this method

Parameters
InGlobalRow - Row of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector that should to modify.
InScalarValue - Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().

Definition at line 363 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ReplaceGlobalValue ( long long  GlobalRow,
int  VectorIndex,
double  ScalarValue 
)

Definition at line 372 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ReplaceGlobalValue ( int  GlobalBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue 
)

Replace current value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex) location with ScalarValue.

Replaces the existing value for a single entry in the multivector. The specified global block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InGlobalBlockRow - BlockRow of Multivector to modify in global index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector that should to modify.
InScalarValue - Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.

Definition at line 381 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ReplaceGlobalValue ( long long  GlobalBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue 
)

Definition at line 390 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::SumIntoGlobalValue ( int  GlobalRow,
int  VectorIndex,
double  ScalarValue 
)

Adds ScalarValue to existing value at the specified (GlobalRow, VectorIndex) location.

Sums the given value into the existing value for a single entry in the multivector. The specified global row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the global row will be modified. To modify a different point entry, use the other version of this method

Parameters
InGlobalRow - Row of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector that should to modify.
InScalarValue - Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().

Definition at line 399 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::SumIntoGlobalValue ( long long  GlobalRow,
int  VectorIndex,
double  ScalarValue 
)

Definition at line 408 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::SumIntoGlobalValue ( int  GlobalBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue 
)

Adds ScalarValue to existing value at the specified (GlobalBlockRow, BlockRowOffset, VectorIndex) location.

Sums the given value into the existing value for a single entry in the multivector. The specified global block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InGlobalBlockRow - BlockRow of Multivector to modify in global index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in global index space.
InVectorIndex - Vector within MultiVector that should to modify.
InScalarValue - Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if GlobalRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.

Definition at line 417 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::SumIntoGlobalValue ( long long  GlobalBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue 
)

Definition at line 426 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ReplaceMyValue ( int  MyRow,
int  VectorIndex,
double  ScalarValue 
)

Replace current value at the specified (MyRow, VectorIndex) location with ScalarValue.

Replaces the existing value for a single entry in the multivector. The specified local row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

This method is intended for use with vectors based on an Epetra_Map. If used on a vector based on a non-trivial Epetra_BlockMap, this will update only block row 0, i.e.

Epetra_MultiVector::ReplaceMyValue ( MyRow, VectorIndex, ScalarValue ) is equivalent to: Epetra_MultiVector::ReplaceMyValue ( 0, MyRow, VectorIndex, ScalarValue )

Parameters
InMyRow - Row of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector that should to modify.
InScalarValue - Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().

Definition at line 434 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ReplaceMyValue ( int  MyBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue 
)

Replace current value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location with ScalarValue.

Replaces the existing value for a single entry in the multivector. The specified local block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InMyBlockRow - BlockRow of Multivector to modify in local index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector that should to modify.
InScalarValue - Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.

Definition at line 441 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::SumIntoMyValue ( int  MyRow,
int  VectorIndex,
double  ScalarValue 
)

Adds ScalarValue to existing value at the specified (MyRow, VectorIndex) location.

Sums the given value into the existing value for a single entry in the multivector. The specified local row must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

If the map associated with this multivector is an Epetra_BlockMap, only the first point entry associated with the local row will be modified. To modify a different point entry, use the other version of this method

Parameters
InMyRow - Row of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector that should to modify.
InScalarValue - Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors().

Definition at line 448 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::SumIntoMyValue ( int  MyBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue 
)

Adds ScalarValue to existing value at the specified (MyBlockRow, BlockRowOffset, VectorIndex) location.

Sums the given value into the existing value for a single entry in the multivector. The specified local block row and block row offset must correspond to a GID owned by the map of the multivector on the calling processor. In other words, this method does not perform cross-processor communication.

Parameters
InMyBlockRow - BlockRow of Multivector to modify in local index space.
InBlockRowOffset - Offset into BlockRow of Multivector to modify in local index space.
InVectorIndex - Vector within MultiVector that should to modify.
InScalarValue - Value to add to existing value.
Returns
Integer error code, set to 0 if successful, set to 1 if MyRow not associated with calling processor set to -1 if VectorIndex >= NumVectors(), set to -2 if BlockRowOffset is out-of-range.

Definition at line 454 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::PutScalar ( double  ScalarConstant)

Initialize all values in a multi-vector with constant value.

Parameters
InScalarConstant - Value to use.
Returns
Integer error code, set to 0 if successful.

Definition at line 595 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Random ( )

Set multi-vector values to random numbers.

MultiVector uses the random number generator provided by Epetra_Util. The multi-vector values will be set to random values on the interval (-1.0, 1.0).

Returns
Integer error code, set to 0 if successful.

Definition at line 490 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ExtractCopy ( double *  A,
int  MyLDA 
) const

Put multi-vector values into user-provided two-dimensional array.

Parameters
OutA - Pointer to memory space that will contain the multi-vector values. The first vector will be copied to the memory pointed to by A. The second vector starts at A+MyLDA, the third at A+2*MyLDA, and so on.
InMyLDA - The "Leading Dimension", or stride between vectors in memory.
Warning
This value refers to the stride on the calling processor. Thus it is a local quantity, not a global quantity.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Definition at line 521 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ExtractCopy ( double **  ArrayOfPointers) const

Put multi-vector values into user-provided array of pointers.

Parameters
OutArrayOfPointers - An array of pointers to memory space that will contain the multi-vector values, such that ArrayOfPointers[i] points to the memory location where the ith vector to be copied.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Definition at line 556 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ExtractView ( double **  A,
int *  MyLDA 
) const

Set user-provided addresses of A and MyLDA.

Parameters
A(Out) - Address of a pointer to that will be set to point to the values of the multi-vector. The first vector will be at the memory pointed to by A. The second vector starts at A+MyLDA, the third at A+2*MyLDA, and so on.
MyLDA(Out) - Address of the "Leading Dimension", or stride between vectors in memory.
Warning
This value refers to the stride on the calling processor. Thus it is a local quantity, not a global quantity.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Definition at line 574 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ExtractView ( double ***  ArrayOfPointers) const

Set user-provided addresses of ArrayOfPointers.

Parameters
ArrayOfPointers(Out) - Address of array of pointers to memory space that will set to the multi-vector array of pointers, such that ArrayOfPointers[i] points to the memory location where the ith vector is located.
Returns
Integer error code, set to 0 if successful.

See Detailed Description section for further discussion.

Definition at line 587 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Dot ( const Epetra_MultiVector A,
double *  Result 
) const

Computes dot product of each corresponding pair of vectors.

Parameters
InA - Multi-vector to be used with the "\e this" multivector.
OutResult - Result[i] will contain the ith dot product result.
Returns
Integer error code, set to 0 if successful.

Definition at line 1123 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Abs ( const Epetra_MultiVector A)

Puts element-wise absolute values of input Multi-vector in target.

Parameters
InA - Input Multi-vector.
Outthis will contain the absolute values of the entries of A.
Returns
Integer error code, set to 0 if successful.

Note: It is possible to use the same argument for A and this.

Definition at line 1158 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Reciprocal ( const Epetra_MultiVector A)

Puts element-wise reciprocal values of input Multi-vector in target.

Parameters
InA - Input Multi-vector.
Outthis will contain the element-wise reciprocal values of the entries of A.
Returns
Integer error code, set to 0 if successful. Returns 2 if some entry is too small, but not zero. Returns 1 if some entry is zero.

Note: It is possible to use the same argument for A and this. Also, if a given value of A is smaller than Epetra_DoubleMin (defined in Epetra_Epetra.h), but nonzero, then the return code is 2. If an entry is zero, the return code is 1. However, in all cases the reciprocal value is still used, even if a NaN is the result.

Definition at line 1180 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Scale ( double  ScalarValue)

Scale the current values of a multi-vector, this = ScalarValue*this.

Parameters
InScalarValue - Scale value.
OutThis - Multi-vector with scaled values.
Returns
Integer error code, set to 0 if successful.

Definition at line 1229 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Scale ( double  ScalarA,
const Epetra_MultiVector A 
)

Replace multi-vector values with scaled values of A, this = ScalarA*A.

Parameters
InScalarA - Scale value.
InA - Multi-vector to copy.
OutThis - Multi-vector with values overwritten by scaled values of A.
Returns
Integer error code, set to 0 if successful.

Definition at line 1251 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Update ( double  ScalarA,
const Epetra_MultiVector A,
double  ScalarThis 
)

Update multi-vector values with scaled values of A, this = ScalarThis*this + ScalarA*A.

Parameters
InScalarA - Scale value for A.
InA - Multi-vector to add.
InScalarThis - Scale value for this.
OutThis - Multi-vector with updatede values.
Returns
Integer error code, set to 0 if successful.

Definition at line 1276 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Update ( double  ScalarA,
const Epetra_MultiVector A,
double  ScalarB,
const Epetra_MultiVector B,
double  ScalarThis 
)

Update multi-vector with scaled values of A and B, this = ScalarThis*this + ScalarA*A + ScalarB*B.

Parameters
InScalarA - Scale value for A.
InA - Multi-vector to add.
InScalarB - Scale value for B.
InB - Multi-vector to add.
InScalarThis - Scale value for this.
OutThis - Multi-vector with updatede values.
Returns
Integer error code, set to 0 if successful.

Definition at line 1341 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Norm1 ( double *  Result) const

Compute 1-norm of each vector in multi-vector.

Parameters
OutResult - Result[i] contains 1-norm of ith vector.
Warning
Map of the this multivector must have unique GIDs (UniqueGIDs() must return true).
Returns
Integer error code, set to 0 if successful.

Definition at line 1508 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Norm2 ( double *  Result) const

Compute 2-norm of each vector in multi-vector.

Parameters
OutResult - Result[i] contains 2-norm of ith vector.
Warning
Map of the this multivector must have unique GIDs (UniqueGIDs() must return true).
Returns
Integer error code, set to 0 if successful.

Definition at line 1547 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::NormInf ( double *  Result) const

Compute Inf-norm of each vector in multi-vector.

Parameters
OutResult - Result[i] contains Inf-norm of ith vector.
Returns
Integer error code, set to 0 if successful.

Definition at line 1580 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::NormWeighted ( const Epetra_MultiVector Weights,
double *  Result 
) const

Compute Weighted 2-norm (RMS Norm) of each vector in multi-vector.

Parameters
InWeights - Multi-vector of weights. If Weights contains a single vector, that vector will be used as the weights for all vectors of this. Otherwise, Weights should have the same number of vectors as this.
OutResult - Result[i] contains the weighted 2-norm of ith vector. Specifically if we denote the ith vector in the multivector by $x$, and the ith weight vector by $w$ and let j represent the jth entry of each vector, on return Result[i] will contain the following result:

\[\sqrt{(1/n)\sum_{j=1}^n(x_j/w_j)^2}\]

, where $n$ is the global length of the vectors.
Returns
Integer error code, set to 0 if successful.

Definition at line 1624 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::MinValue ( double *  Result) const

Compute minimum value of each vector in multi-vector.

Note that the vector contents must be already initialized for this function to compute a well-defined result. The length of the vector need not be greater than zero on all processors. If length is greater than zero on any processor then a valid result will be computed.

Parameters
OutResult - Result[i] contains minimum value of ith vector.
Returns
Integer error code, set to 0 if successful.

Definition at line 1674 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::MaxValue ( double *  Result) const

Compute maximum value of each vector in multi-vector.

Note that the vector contents must be already initialized for this function to compute a well-defined result. The length of the vector need not be greater than zero on all processors. If length is greater than zero on any processor then a valid result will be computed.

Parameters
OutResult - Result[i] contains maximum value of ith vector.
Returns
Integer error code, set to 0 if successful.

Definition at line 1787 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::MeanValue ( double *  Result) const

Compute mean (average) value of each vector in multi-vector.

Parameters
OutResult - Result[i] contains mean value of ith vector.
Warning
Map of the this multivector must have unique GIDs (UniqueGIDs() must return true).
Returns
Integer error code, set to 0 if successful.

Definition at line 1900 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Multiply ( char  TransA,
char  TransB,
double  ScalarAB,
const Epetra_MultiVector A,
const Epetra_MultiVector B,
double  ScalarThis 
)

Matrix-Matrix multiplication, this = ScalarThis*this + ScalarAB*A*B.

This function performs a variety of matrix-matrix multiply operations, interpreting the Epetra_MultiVectors (this-aka C , A and B) as 2D matrices. Variations are due to the fact that A, B and C can be local replicated or global distributed Epetra_MultiVectors and that we may or may not operate with the transpose of A and B. Possible cases are:

Total of 32 case (2^5).
Num
OPERATIONS                        case  Notes
1) C(local) = A^X(local) * B^X(local)  4   (X=Transpose or Not, No comm needed)
2) C(local) = A^T(distr) * B  (distr)  1   (2D dot product, replicate C)
3) C(distr) = A  (distr) * B^X(local)  2   (2D vector update, no comm needed)

Note that the following operations are not meaningful for
1D distributions:

1) C(local) = A^T(distr) * B^T(distr)  1
2) C(local) = A  (distr) * B^X(distr)  2
3) C(distr) = A^X(local) * B^X(local)  4
4) C(distr) = A^X(local) * B^X(distr)  4
5) C(distr) = A^T(distr) * B^X(local)  2
6) C(local) = A^X(distr) * B^X(local)  4
7) C(distr) = A^X(distr) * B^X(local)  4
8) C(local) = A^X(local) * B^X(distr)  4
Parameters
InTransA - Operate with the transpose of A if = 'T', else no transpose if = 'N'.
InTransB - Operate with the transpose of B if = 'T', else no transpose if = 'N'.
InScalarAB - Scalar to multiply with A*B.
InA - Multi-vector.
InB - Multi-vector.
InScalarThis - Scalar to multiply with this.
Warning
Map of the distributed multivectors must have unique GIDs (UniqueGIDs() must return true).
Returns
Integer error code, set to 0 if successful.
Warning
{Each multi-vector A, B and this is checked if it has constant stride using the ConstantStride() query function. If it does not have constant stride, a temporary copy is made and used for the computation. This activity is transparent to the user, except that there is memory and computation overhead. All temporary space is deleted prior to exit.}

A.ConstantStride() || B.ConstantStride() ) EPETRA_CHK_ERR(-1); // Return error

Definition at line 1936 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::Multiply ( double  ScalarAB,
const Epetra_MultiVector A,
const Epetra_MultiVector B,
double  ScalarThis 
)

Multiply a Epetra_MultiVector with another, element-by-element.

This function supports diagonal matrix multiply. A is usually a single vector while B and this may have one or more columns. Note that B and this must have the same shape. A can be one vector or have the same shape as B. The actual computation is this = ScalarThis * this + ScalarAB * B @ A where @ denotes element-wise multiplication.

Definition at line 2087 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ReciprocalMultiply ( double  ScalarAB,
const Epetra_MultiVector A,
const Epetra_MultiVector B,
double  ScalarThis 
)

Multiply a Epetra_MultiVector by the reciprocal of another, element-by-element.

This function supports diagonal matrix scaling. A is usually a single vector while B and this may have one or more columns. Note that B and this must have the same shape. A can be one vector or have the same shape as B. The actual computation is this = ScalarThis * this + ScalarAB * B @ A where @ denotes element-wise division.

Definition at line 2210 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::SetSeed ( unsigned int  Seed_in)
inline

Set seed for Random function.

Parameters
InSeed - Should be an integer on the interval (0, 2^31-1).
Returns
Integer error code, set to 0 if successful.

Definition at line 870 of file Epetra_MultiVector.h.

unsigned int Epetra_MultiVector::Seed ( )
inline

Get seed from Random function.

Returns
Current random number seed.

Definition at line 876 of file Epetra_MultiVector.h.

Epetra_MultiVector & Epetra_MultiVector::operator= ( const Epetra_MultiVector Source)

= Operator.

Parameters
InA - Epetra_MultiVector to copy.
Returns
Epetra_MultiVector.

Definition at line 2368 of file Epetra_MultiVector.cpp.

double*& Epetra_MultiVector::operator[] ( int  i)
inline

Vector access function.

Returns
Pointer to the array of doubles containing the local values of the ith vector in the multi-vector.

Definition at line 900 of file Epetra_MultiVector.h.

double* const& Epetra_MultiVector::operator[] ( int  i) const
inline

Vector access function.

Returns
Pointer to the array of doubles containing the local values of the ith vector in the multi-vector.

Definition at line 906 of file Epetra_MultiVector.h.

Epetra_Vector *& Epetra_MultiVector::operator() ( int  i)

Vector access function.

Returns
An Epetra_Vector pointer to the ith vector in the multi-vector.

Definition at line 2335 of file Epetra_MultiVector.cpp.

const Epetra_Vector *& Epetra_MultiVector::operator() ( int  i) const

Vector access function.

Returns
An Epetra_Vector pointer to the ith vector in the multi-vector.

Definition at line 2351 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::NumVectors ( ) const
inline

Returns the number of vectors in the multi-vector.

Definition at line 925 of file Epetra_MultiVector.h.

int Epetra_MultiVector::MyLength ( ) const
inline

Returns the local vector length on the calling processor of vectors in the multi-vector.

Definition at line 928 of file Epetra_MultiVector.h.

int Epetra_MultiVector::GlobalLength ( ) const
inline

Returns the global vector length of vectors in the multi-vector.

Definition at line 931 of file Epetra_MultiVector.h.

long long Epetra_MultiVector::GlobalLength64 ( ) const
inline

Returns the 64-bit global vector length of vectors in the multi-vector.

Definition at line 934 of file Epetra_MultiVector.h.

int Epetra_MultiVector::Stride ( ) const
inline

Returns the stride between vectors in the multi-vector (only meaningful if ConstantStride() is true).

Definition at line 937 of file Epetra_MultiVector.h.

bool Epetra_MultiVector::ConstantStride ( ) const
inline

Returns true if this multi-vector has constant stride between vectors.

Definition at line 940 of file Epetra_MultiVector.h.

int Epetra_MultiVector::ReplaceMap ( const Epetra_BlockMap map)

Replace map, only if new map has same point-structure as current map.

return 0 if map is replaced, -1 if not.

Definition at line 536 of file Epetra_MultiVector.cpp.

void Epetra_MultiVector::Print ( std::ostream &  os) const
virtual

Print method.

Reimplemented from Epetra_DistObject.

Definition at line 2447 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ResetView ( double **  ArrayOfPointers)

Reset the view of an existing multivector to point to new user data.

Allows the (very) light-weight replacement of multivector values for an existing multivector that was constructed using an Epetra_DataAccess mode of View. No checking is performed to see if the array of values passed in contains valid data. It is assumed that the user has verified the integrity of data before calling this method. This method is useful for situations where a multivector is needed for use with an Epetra operator or matrix and the user is not passing in a multivector, or the multivector is being passed in with another map that is not exactly compatible with the operator, but has the correct number of entries.

This method is used by AztecOO and Ifpack in the matvec, and solve methods to improve performance and reduce repeated calls to constructors and destructors.

Parameters
ArrayOfPointersContains the array of pointers containing the multivector data.
Returns
Integer error code, set to 0 if successful, -1 if the multivector was not created as a View.
Warning
This method is extremely dangerous and should only be used by experts.

Definition at line 2435 of file Epetra_MultiVector.cpp.

double* Epetra_MultiVector::Values ( ) const
inline

Get pointer to MultiVector values.

Definition at line 981 of file Epetra_MultiVector.h.

double** Epetra_MultiVector::Pointers ( ) const
inline

Get pointer to individual vector pointers.

Definition at line 984 of file Epetra_MultiVector.h.

int Epetra_MultiVector::Reduce ( )

Definition at line 2400 of file Epetra_MultiVector.cpp.

void Epetra_MultiVector::Assign ( const Epetra_MultiVector rhs)
protected

Definition at line 2377 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::CheckInput ( )
protected
int Epetra_MultiVector::AllocateForCopy ( void  )
private

Definition at line 248 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::DoCopy ( void  )
private

Definition at line 278 of file Epetra_MultiVector.cpp.

void Epetra_MultiVector::UpdateDoubleTemp ( ) const
inlineprivate

Definition at line 1009 of file Epetra_MultiVector.h.

void Epetra_MultiVector::UpdateVectors ( ) const
inlineprivate

Definition at line 1012 of file Epetra_MultiVector.h.

int Epetra_MultiVector::AllocateForView ( void  )
private

Definition at line 306 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::DoView ( void  )
private

Definition at line 333 of file Epetra_MultiVector.cpp.

template<typename int_type >
int Epetra_MultiVector::ChangeGlobalValue ( int_type  GlobalBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue,
bool  SumInto 
)
private

Definition at line 462 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::ChangeMyValue ( int  MyBlockRow,
int  BlockRowOffset,
int  VectorIndex,
double  ScalarValue,
bool  SumInto 
)
private

Definition at line 473 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::CheckSizes ( const Epetra_SrcDistObject Source)
privatevirtual

Allows the source and target (this) objects to be compared for compatibility, return nonzero if not.

Implements Epetra_DistObject.

Definition at line 610 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::CopyAndPermute ( const Epetra_SrcDistObject Source,
int  NumSameIDs,
int  NumPermuteIDs,
int *  PermuteToLIDs,
int *  PermuteFromLIDs,
const Epetra_OffsetIndex Indexor,
Epetra_CombineMode  CombineMode = Zero 
)
privatevirtual

Perform ID copies and permutations that are on processor.

Implements Epetra_DistObject.

Definition at line 619 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::PackAndPrepare ( const Epetra_SrcDistObject Source,
int  NumExportIDs,
int *  ExportLIDs,
int &  LenExports,
char *&  Exports,
int &  SizeOfPacket,
int *  Sizes,
bool &  VarSizes,
Epetra_Distributor Distor 
)
privatevirtual

Perform any packing or preparation required for call to DoTransfer().

Implements Epetra_DistObject.

Definition at line 724 of file Epetra_MultiVector.cpp.

int Epetra_MultiVector::UnpackAndCombine ( const Epetra_SrcDistObject Source,
int  NumImportIDs,
int *  ImportLIDs,
int  LenImports,
char *  Imports,
int &  SizeOfPacket,
Epetra_Distributor Distor,
Epetra_CombineMode  CombineMode,
const Epetra_OffsetIndex Indexor 
)
privatevirtual

Perform any unpacking and combining after call to DoTransfer().

Implements Epetra_DistObject.

Definition at line 816 of file Epetra_MultiVector.cpp.

Member Data Documentation

double* Epetra_MultiVector::Values_
protected

Definition at line 999 of file Epetra_MultiVector.h.

double** Epetra_MultiVector::Pointers_
private

Definition at line 1061 of file Epetra_MultiVector.h.

int Epetra_MultiVector::MyLength_
private

Definition at line 1063 of file Epetra_MultiVector.h.

long long Epetra_MultiVector::GlobalLength_
private

Definition at line 1064 of file Epetra_MultiVector.h.

int Epetra_MultiVector::NumVectors_
private

Definition at line 1065 of file Epetra_MultiVector.h.

bool Epetra_MultiVector::UserAllocated_
private

Definition at line 1066 of file Epetra_MultiVector.h.

bool Epetra_MultiVector::ConstantStride_
private

Definition at line 1067 of file Epetra_MultiVector.h.

int Epetra_MultiVector::Stride_
private

Definition at line 1068 of file Epetra_MultiVector.h.

bool Epetra_MultiVector::Allocated_
private

Definition at line 1069 of file Epetra_MultiVector.h.

double* Epetra_MultiVector::DoubleTemp_
mutableprivate

Definition at line 1070 of file Epetra_MultiVector.h.

Epetra_Vector** Epetra_MultiVector::Vectors_
mutableprivate

Definition at line 1071 of file Epetra_MultiVector.h.

Epetra_Util Epetra_MultiVector::Util_
private

Definition at line 1072 of file Epetra_MultiVector.h.


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