EpetraExt
Development
|
#include <EpetraExt_BlockVector.h>
Public Member Functions | |
int | ExtractBlockValues (Epetra_Vector &BaseVec, long long BlockRow) const |
Extract a single block from a Block Vector: block row is global, not a stencil value. More... | |
int | LoadBlockValues (const Epetra_Vector &BaseVec, long long BlockRow) |
Load a single block into a Block Vector: block row is global, not a stencil value. More... | |
int | BlockSumIntoGlobalValues (int NumIndices, double *Values, int *Indices, int BlockRow) |
Load entries into BlockVector with base vector indices offset by BlockRow. More... | |
int | BlockReplaceGlobalValues (int NumIndices, double *Values, int *Indices, int BlockRow) |
Load entries into BlockVector with base vector indices offset by BlockRow. More... | |
int | BlockSumIntoGlobalValues (int NumIndices, double *Values, long long *Indices, long long BlockRow) |
Load entries into BlockVector with base vector indices offset by BlockRow. More... | |
int | BlockReplaceGlobalValues (int NumIndices, double *Values, long long *Indices, long long BlockRow) |
Load entries into BlockVector with base vector indices offset by BlockRow. More... | |
Teuchos::RCP< const Epetra_Vector > | GetBlock (long long BlockRow) const |
Return Epetra_Vector for given block row. More... | |
Teuchos::RCP< Epetra_Vector > | GetBlock (long long BlockRow) |
Return Epetra_Vector for given block row. More... | |
const Epetra_BlockMap & | GetBaseMap () const |
Return base map. More... | |
Public Member Functions inherited from Epetra_Vector | |
Epetra_Vector (const Epetra_BlockMap &Map, bool zeroOut=true) | |
Epetra_Vector (const Epetra_Vector &Source) | |
Epetra_Vector (Epetra_DataAccess CV, const Epetra_BlockMap &Map, double *V) | |
Epetra_Vector (Epetra_DataAccess CV, const Epetra_MultiVector &Source, int Index) | |
virtual | ~Epetra_Vector () |
int | ReplaceGlobalValues (int NumEntries, const double *Values, const int *Indices) |
int | ReplaceMyValues (int NumEntries, const double *Values, const int *Indices) |
int | SumIntoGlobalValues (int NumEntries, const double *Values, const int *Indices) |
int | SumIntoMyValues (int NumEntries, const double *Values, const int *Indices) |
int | ReplaceGlobalValues (int NumEntries, int BlockOffset, const double *Values, const int *Indices) |
int | ReplaceMyValues (int NumEntries, int BlockOffset, const double *Values, const int *Indices) |
int | SumIntoGlobalValues (int NumEntries, int BlockOffset, const double *Values, const int *Indices) |
int | SumIntoMyValues (int NumEntries, int BlockOffset, const double *Values, const int *Indices) |
int | ExtractCopy (double *V) const |
int | ExtractView (double **V) const |
double & | operator[] (int index) |
const double & | operator[] (int index) const |
int | ResetView (double *Values_in) |
Public Member Functions inherited from Epetra_DistObject | |
virtual void | Print (std::ostream &os) const |
Public Member Functions inherited from Epetra_BLAS | |
Epetra_BLAS (void) | |
Epetra_BLAS (const Epetra_BLAS &BLAS) | |
virtual | ~Epetra_BLAS (void) |
float | ASUM (const int N, const float *X, const int INCX=1) const |
double | ASUM (const int N, const double *X, const int INCX=1) const |
float | DOT (const int N, const float *X, const float *Y, const int INCX=1, const int INCY=1) const |
double | DOT (const int N, const double *X, const double *Y, const int INCX=1, const int INCY=1) const |
float | NRM2 (const int N, const float *X, const int INCX=1) const |
double | NRM2 (const int N, const double *X, const int INCX=1) const |
void | SCAL (const int N, const float ALPHA, float *X, const int INCX=1) const |
void | SCAL (const int N, const double ALPHA, double *X, const int INCX=1) const |
void | COPY (const int N, const float *X, float *Y, const int INCX=1, const int INCY=1) const |
void | COPY (const int N, const double *X, double *Y, const int INCX=1, const int INCY=1) const |
int | IAMAX (const int N, const float *X, const int INCX=1) const |
int | IAMAX (const int N, const double *X, const int INCX=1) const |
void | AXPY (const int N, const float ALPHA, const float *X, float *Y, const int INCX=1, const int INCY=1) const |
void | AXPY (const int N, const double ALPHA, const double *X, double *Y, const int INCX=1, const int INCY=1) const |
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 |
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 |
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 |
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 |
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 |
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 |
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 |
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 |
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 |
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 |
Protected Attributes | |
Epetra_BlockMap | BaseMap_ |
long long | Offset_ |
BlockVector (const Epetra_BlockMap &BaseMap, const Epetra_BlockMap &GlobalMap) | |
BlockVector constuctor with one block row per processor. More... | |
BlockVector (Epetra_DataAccess CV, const Epetra_BlockMap &BaseMap, const Epetra_Vector &BlockVec) | |
BlockVector (const BlockVector &MV) | |
Copy constructor. More... | |
virtual | ~BlockVector () |
Destructor. More... | |
Additional Inherited Members | |
Protected Member Functions inherited from Epetra_DistObject | |
virtual int | CheckSizes (const Epetra_SrcDistObject &Source)=0 |
virtual int | CopyAndPermute (const Epetra_SrcDistObject &Source, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode=Zero)=0 |
virtual int | PackAndPrepare (const Epetra_SrcDistObject &Source, int NumExportIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &SizeOfPacket, int *Sizes, bool &VarSizes, Epetra_Distributor &Distor)=0 |
virtual 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)=0 |
Definition at line 60 of file EpetraExt_BlockVector.h.
EpetraExt::BlockVector::BlockVector | ( | const Epetra_BlockMap & | BaseMap, |
const Epetra_BlockMap & | GlobalMap | ||
) |
BlockVector constuctor with one block row per processor.
Creates a BlockVector object and allocates storage.
In | BaseMap - Map determining local structure, can be distrib. over subset of proc.'s |
In | GlobalMap - Full map describing the overall global structure, generally generated by the construction of a BlockCrsMatrix object |
In | NumBlocks - Number of local blocks |
Definition at line 51 of file EpetraExt_BlockVector.cpp.
EpetraExt::BlockVector::BlockVector | ( | Epetra_DataAccess | CV, |
const Epetra_BlockMap & | BaseMap, | ||
const Epetra_Vector & | BlockVec | ||
) |
Creates a BlockVector object from an existing Epetra_Vector.
In | Epetra_DataAccess - Enumerated type set to Copy or View. |
In | BaseMap - Map determining local structure, can be distrib. over subset of proc.'s |
In | BlockVec - Source Epetra vector whose map must be the full map for the block vector |
Definition at line 62 of file EpetraExt_BlockVector.cpp.
EpetraExt::BlockVector::BlockVector | ( | const BlockVector & | MV | ) |
Copy constructor.
Definition at line 74 of file EpetraExt_BlockVector.cpp.
|
virtual |
Destructor.
Definition at line 82 of file EpetraExt_BlockVector.cpp.
int EpetraExt::BlockVector::ExtractBlockValues | ( | Epetra_Vector & | BaseVec, |
long long | BlockRow | ||
) | const |
Extract a single block from a Block Vector: block row is global, not a stencil value.
Definition at line 87 of file EpetraExt_BlockVector.cpp.
int EpetraExt::BlockVector::LoadBlockValues | ( | const Epetra_Vector & | BaseVec, |
long long | BlockRow | ||
) |
Load a single block into a Block Vector: block row is global, not a stencil value.
Definition at line 108 of file EpetraExt_BlockVector.cpp.
int EpetraExt::BlockVector::BlockSumIntoGlobalValues | ( | int | NumIndices, |
double * | Values, | ||
int * | Indices, | ||
int | BlockRow | ||
) |
Load entries into BlockVector with base vector indices offset by BlockRow.
Definition at line 129 of file EpetraExt_BlockVector.cpp.
int EpetraExt::BlockVector::BlockReplaceGlobalValues | ( | int | NumIndices, |
double * | Values, | ||
int * | Indices, | ||
int | BlockRow | ||
) |
Load entries into BlockVector with base vector indices offset by BlockRow.
Definition at line 150 of file EpetraExt_BlockVector.cpp.
int EpetraExt::BlockVector::BlockSumIntoGlobalValues | ( | int | NumIndices, |
double * | Values, | ||
long long * | Indices, | ||
long long | BlockRow | ||
) |
Load entries into BlockVector with base vector indices offset by BlockRow.
Definition at line 173 of file EpetraExt_BlockVector.cpp.
int EpetraExt::BlockVector::BlockReplaceGlobalValues | ( | int | NumIndices, |
double * | Values, | ||
long long * | Indices, | ||
long long | BlockRow | ||
) |
Load entries into BlockVector with base vector indices offset by BlockRow.
Definition at line 194 of file EpetraExt_BlockVector.cpp.
Teuchos::RCP< const Epetra_Vector > EpetraExt::BlockVector::GetBlock | ( | long long | BlockRow | ) | const |
Return Epetra_Vector for given block row.
Definition at line 217 of file EpetraExt_BlockVector.cpp.
Teuchos::RCP< Epetra_Vector > EpetraExt::BlockVector::GetBlock | ( | long long | BlockRow | ) |
Return Epetra_Vector for given block row.
Definition at line 225 of file EpetraExt_BlockVector.cpp.
const Epetra_BlockMap & EpetraExt::BlockVector::GetBaseMap | ( | ) | const |
Return base map.
Definition at line 233 of file EpetraExt_BlockVector.cpp.
|
protected |
Definition at line 130 of file EpetraExt_BlockVector.h.
|
protected |
Definition at line 132 of file EpetraExt_BlockVector.h.