Epetra_MultiVector: A class for constructing and using dense multi-vectors, vectors and matrices in parallel. More...
#include <Epetra_MultiVector.h>
Public Member Functions | |
int | ReplaceMap (const Epetra_BlockMap &map) |
int | Reduce () |
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. | |
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. | |
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_MultiVector & | operator= (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. | |
int | MyLength () const |
Returns the local vector length on the calling processor of vectors in the multi-vector. | |
int | GlobalLength () const |
Returns the global vector length of vectors in the multi-vector. | |
long long | GlobalLength64 () const |
Returns the 64-bit global vector length of vectors in the multi-vector. | |
int | Stride () const |
Returns the stride between vectors in the multi-vector (only meaningful if ConstantStride() is true). | |
bool | ConstantStride () const |
Returns true if this multi-vector has constant stride between vectors. | |
I/O methods | |
virtual void | Print (std::ostream &os) const |
Print method. | |
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. | |
double ** | Pointers () const |
Get pointer to individual vector pointers. | |
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. | |
virtual | ~Epetra_DistObject () |
Epetra_DistObject destructor. | |
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_BlockMap & | Map () const |
Returns the address of the Epetra_BlockMap for this multi-vector. | |
const Epetra_Comm & | Comm () const |
Returns the address of the Epetra_Comm for this multi-vector. | |
bool | DistributedGlobal () const |
Returns true if this multi-vector is distributed global, i.e., not local replicated. | |
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 int | ReportError (const std::string Message, int ErrorCode) const |
Error reporting method. | |
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... | |
Public Member Functions inherited from Epetra_SrcDistObject | |
virtual | ~Epetra_SrcDistObject () |
Epetra_SrcDistObject destructor. | |
Public Member Functions inherited from Epetra_CompObject | |
Epetra_CompObject & | operator= (const Epetra_CompObject &src) |
Epetra_CompObject () | |
Basic Epetra_CompObject constuctor. | |
Epetra_CompObject (const Epetra_CompObject &Source) | |
Epetra_CompObject copy constructor. | |
virtual | ~Epetra_CompObject () |
Epetra_CompObject destructor. | |
void | SetFlopCounter (const Epetra_Flops &FlopCounter_in) |
Set the internal Epetra_Flops() pointer. | |
void | SetFlopCounter (const Epetra_CompObject &CompObject) |
Set the internal Epetra_Flops() pointer to the flop counter of another Epetra_CompObject. | |
void | UnsetFlopCounter () |
Set the internal Epetra_Flops() pointer to 0 (no flops counted). | |
Epetra_Flops * | GetFlopCounter () const |
Get the pointer to the Epetra_Flops() object associated with this object, returns 0 if none. | |
void | ResetFlops () const |
Resets the number of floating point operations to zero for this multi-vector. | |
double | Flops () const |
Returns the number of floating point operations with this multi-vector. | |
void | UpdateFlops (int Flops_in) const |
Increment Flop count for this object. | |
void | UpdateFlops (long int Flops_in) const |
Increment Flop count for this object. | |
void | UpdateFlops (long long Flops_in) const |
Increment Flop count for this object. | |
void | UpdateFlops (double Flops_in) const |
Increment Flop count for this object. | |
void | UpdateFlops (float Flops_in) const |
Increment Flop count for this object. | |
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. | |
float | ASUM (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS one norm function (SASUM). | |
double | ASUM (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS one norm function (DASUM). | |
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). | |
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). | |
float | NRM2 (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS norm function (SNRM2). | |
double | NRM2 (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS norm function (DNRM2). | |
void | SCAL (const int N, const float ALPHA, float *X, const int INCX=1) const |
Epetra_BLAS vector scale function (SSCAL) | |
void | SCAL (const int N, const double ALPHA, double *X, const int INCX=1) const |
Epetra_BLAS vector scale function (DSCAL) | |
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) | |
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) | |
int | IAMAX (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS arg maximum of absolute value function (ISAMAX) | |
int | IAMAX (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS arg maximum of absolute value function (IDAMAX) | |
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) | |
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) | |
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) | |
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) | |
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) | |
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) | |
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) | |
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) | |
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) | |
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) | |
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) | |
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) | |
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. | |
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_Comm * | Comm_ |
char * | Exports_ |
char * | Imports_ |
int | LenExports_ |
int | LenImports_ |
int * | Sizes_ |
Protected Attributes inherited from Epetra_CompObject | |
Epetra_Flops * | FlopCounter_ |
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. | |
static std::ostream & | GetTracebackStream () |
Get the output stream for error reporting. | |
Static Public Attributes inherited from Epetra_Object | |
static int | TracebackMode |
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:
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.
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.
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.
In | Map - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap. |
In | NumVectors - Number of vectors in multi-vector. |
In | zeroOut - 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. |
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.
In | Epetra_DataAccess - Enumerated type set to Copy or View. |
In | Map - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap. |
In | A - 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. |
In | MyLDA - The "Leading Dimension", or stride between vectors in memory. |
In | NumVectors - Number of vectors in multi-vector. |
See Detailed Description section for further discussion.
Epetra_MultiVector::Epetra_MultiVector | ( | Epetra_DataAccess | CV, |
const Epetra_BlockMap & | Map, | ||
double ** | ArrayOfPointers, | ||
int | NumVectors | ||
) |
Set multi-vector values from array of pointers.
In | Epetra_DataAccess - Enumerated type set to Copy or View. |
In | Map - A Epetra_LocalMap, Epetra_Map or Epetra_BlockMap. |
In | ArrayOfPointers - An array of pointers such that ArrayOfPointers[i] points to the memory location containing ith vector to be copied. |
In | NumVectors - Number of vectors in multi-vector. |
See Detailed Description section for further discussion.
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.
In | Epetra_DataAccess - Enumerated type set to Copy or View. |
In | Source - An existing fully constructed Epetra_MultiVector. |
In | Indices - Integer list of the vectors to copy. |
In | NumVectors - Number of vectors in multi-vector. |
See Detailed Description section for further discussion.
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.
In | Epetra_DataAccess - Enumerated type set to Copy or View. |
In | Source - An existing fully constructed Epetra_MultiVector. |
In | StartIndex - First of the vectors to copy. |
In | NumVectors - Number of vectors in multi-vector. |
See Detailed Description section for further discussion.
int Epetra_MultiVector::Abs | ( | const Epetra_MultiVector & | A | ) |
Puts element-wise absolute values of input Multi-vector in target.
In | A - Input Multi-vector. |
Out | this will contain the absolute values of the entries of A. |
Note: It is possible to use the same argument for A and this.
int Epetra_MultiVector::Dot | ( | const Epetra_MultiVector & | A, |
double * | Result | ||
) | const |
Computes dot product of each corresponding pair of vectors.
In | A - Multi-vector to be used with the "\e this" multivector. |
Out | Result - Result[i] will contain the ith dot product result. |
int Epetra_MultiVector::ExtractCopy | ( | double * | A, |
int | MyLDA | ||
) | const |
Put multi-vector values into user-provided two-dimensional array.
Out | A - 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. |
In | MyLDA - The "Leading Dimension", or stride between vectors in memory. |
See Detailed Description section for further discussion.
int Epetra_MultiVector::ExtractCopy | ( | double ** | ArrayOfPointers | ) | const |
Put multi-vector values into user-provided array of pointers.
Out | ArrayOfPointers - 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. |
See Detailed Description section for further discussion.
int Epetra_MultiVector::ExtractView | ( | double ** | A, |
int * | MyLDA | ||
) | const |
Set user-provided addresses of A and MyLDA.
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. |
See Detailed Description section for further discussion.
int Epetra_MultiVector::ExtractView | ( | double *** | ArrayOfPointers | ) | const |
Set user-provided addresses of ArrayOfPointers.
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. |
See Detailed Description section for further discussion.
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.
Out | Result - Result[i] contains maximum value of ith vector. |
int Epetra_MultiVector::MeanValue | ( | double * | Result | ) | const |
Compute mean (average) value of each vector in multi-vector.
Out | Result - Result[i] contains mean value of ith vector. |
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.
Out | Result - Result[i] contains minimum value of ith vector. |
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
In | TransA - Operate with the transpose of A if = 'T', else no transpose if = 'N'. |
In | TransB - Operate with the transpose of B if = 'T', else no transpose if = 'N'. |
In | ScalarAB - Scalar to multiply with A*B. |
In | A - Multi-vector. |
In | B - Multi-vector. |
In | ScalarThis - Scalar to multiply with this. |
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.
int Epetra_MultiVector::Norm1 | ( | double * | Result | ) | const |
Compute 1-norm of each vector in multi-vector.
Out | Result - Result[i] contains 1-norm of ith vector. |
int Epetra_MultiVector::Norm2 | ( | double * | Result | ) | const |
Compute 2-norm of each vector in multi-vector.
Out | Result - Result[i] contains 2-norm of ith vector. |
int Epetra_MultiVector::NormInf | ( | double * | Result | ) | const |
Compute Inf-norm of each vector in multi-vector.
Out | Result - Result[i] contains Inf-norm of ith vector. |
int Epetra_MultiVector::NormWeighted | ( | const Epetra_MultiVector & | Weights, |
double * | Result | ||
) | const |
Compute Weighted 2-norm (RMS Norm) of each vector in multi-vector.
In | Weights - 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. |
Out | Result - Result[i] contains the weighted 2-norm of ith vector. Specifically if we denote the ith vector in the multivector by , and the ith weight vector by and let j represent the jth entry of each vector, on return Result[i] will contain the following result: , where is the global length of the vectors. |
Epetra_Vector* & Epetra_MultiVector::operator() | ( | int | i | ) |
Vector access function.
const Epetra_Vector* & Epetra_MultiVector::operator() | ( | int | i | ) | const |
Vector access function.
Epetra_MultiVector& Epetra_MultiVector::operator= | ( | const Epetra_MultiVector & | Source | ) |
|
inline |
Vector access function.
|
inline |
Vector access function.
int Epetra_MultiVector::PutScalar | ( | double | ScalarConstant | ) |
Initialize all values in a multi-vector with constant value.
In | ScalarConstant - Value to use. |
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).
int Epetra_MultiVector::Reciprocal | ( | const Epetra_MultiVector & | A | ) |
Puts element-wise reciprocal values of input Multi-vector in target.
In | A - Input Multi-vector. |
Out | this will contain the element-wise reciprocal values of the entries of A. |
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.
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.
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
In | GlobalRow - Row of Multivector to modify in global index space. |
In | VectorIndex - Vector within MultiVector that should to modify. |
In | ScalarValue - Value to add to existing value. |
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.
In | GlobalBlockRow - BlockRow of Multivector to modify in global index space. |
In | BlockRowOffset - Offset into BlockRow of Multivector to modify in global index space. |
In | VectorIndex - Vector within MultiVector that should to modify. |
In | ScalarValue - Value to add to existing value. |
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.
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 )
In | MyRow - Row of Multivector to modify in local index space. |
In | VectorIndex - Vector within MultiVector that should to modify. |
In | ScalarValue - Value to add to existing value. |
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.
In | MyBlockRow - BlockRow of Multivector to modify in local index space. |
In | BlockRowOffset - Offset into BlockRow of Multivector to modify in local index space. |
In | VectorIndex - Vector within MultiVector that should to modify. |
In | ScalarValue - Value to add to existing value. |
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.
ArrayOfPointers | Contains the array of pointers containing the multivector data. |
Referenced by Epetra_Vector::ResetView().
int Epetra_MultiVector::Scale | ( | double | ScalarValue | ) |
Scale the current values of a multi-vector, this = ScalarValue*this.
In | ScalarValue - Scale value. |
Out | This - Multi-vector with scaled values. |
int Epetra_MultiVector::Scale | ( | double | ScalarA, |
const Epetra_MultiVector & | A | ||
) |
Replace multi-vector values with scaled values of A, this = ScalarA*A.
In | ScalarA - Scale value. |
In | A - Multi-vector to copy. |
Out | This - Multi-vector with values overwritten by scaled values of A. |
|
inline |
Get seed from Random function.
|
inline |
Set seed for Random function.
In | Seed - Should be an integer on the interval (0, 2^31-1). |
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
In | GlobalRow - Row of Multivector to modify in global index space. |
In | VectorIndex - Vector within MultiVector that should to modify. |
In | ScalarValue - Value to add to existing value. |
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.
In | GlobalBlockRow - BlockRow of Multivector to modify in global index space. |
In | BlockRowOffset - Offset into BlockRow of Multivector to modify in global index space. |
In | VectorIndex - Vector within MultiVector that should to modify. |
In | ScalarValue - Value to add to existing value. |
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
In | MyRow - Row of Multivector to modify in local index space. |
In | VectorIndex - Vector within MultiVector that should to modify. |
In | ScalarValue - Value to add to existing value. |
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.
In | MyBlockRow - BlockRow of Multivector to modify in local index space. |
In | BlockRowOffset - Offset into BlockRow of Multivector to modify in local index space. |
In | VectorIndex - Vector within MultiVector that should to modify. |
In | ScalarValue - Value to add to existing value. |
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.
In | ScalarA - Scale value for A. |
In | A - Multi-vector to add. |
In | ScalarThis - Scale value for this. |
Out | This - Multi-vector with updatede values. |
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.
In | ScalarA - Scale value for A. |
In | A - Multi-vector to add. |
In | ScalarB - Scale value for B. |
In | B - Multi-vector to add. |
In | ScalarThis - Scale value for this. |
Out | This - Multi-vector with updatede values. |