Epetra  Development
 All Classes Namespaces Files Functions Variables Enumerations Enumerator Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Epetra_FEVector Class Reference

#include <Epetra_FEVector.h>

Inheritance diagram for Epetra_FEVector:
Inheritance graph
[legend]
Collaboration diagram for Epetra_FEVector:
Collaboration graph
[legend]

Public Member Functions

 Epetra_FEVector (const Epetra_BlockMap &Map, int numVectors=1, bool ignoreNonLocalEntries=false)
 
 Epetra_FEVector (Epetra_DataAccess CV, const Epetra_BlockMap &Map, double *A, int MyLDA, int NumVectors, bool ignoreNonLocalEntries=false)
 Set multi-vector values from two-dimensional array. More...
 
 Epetra_FEVector (Epetra_DataAccess CV, const Epetra_BlockMap &Map, double **ArrayOfPointers, int NumVectors, bool ignoreNonLocalEntries=false)
 Set multi-vector values from array of pointers. More...
 
 Epetra_FEVector (const Epetra_FEVector &source)
 
virtual ~Epetra_FEVector ()
 
int SumIntoGlobalValues (int numIDs, const int *GIDs, const double *values, int vectorIndex=0)
 
int SumIntoGlobalValues (int numIDs, const long long *GIDs, const double *values, int vectorIndex=0)
 
int SumIntoGlobalValues (const Epetra_IntSerialDenseVector &GIDs, const Epetra_SerialDenseVector &values, int vectorIndex=0)
 
int SumIntoGlobalValues (const Epetra_LongLongSerialDenseVector &GIDs, const Epetra_SerialDenseVector &values, int vectorIndex=0)
 
int ReplaceGlobalValues (int numIDs, const int *GIDs, const double *values, int vectorIndex=0)
 
int ReplaceGlobalValues (int numIDs, const long long *GIDs, const double *values, int vectorIndex=0)
 
int ReplaceGlobalValues (const Epetra_IntSerialDenseVector &GIDs, const Epetra_SerialDenseVector &values, int vectorIndex=0)
 
int ReplaceGlobalValues (const Epetra_LongLongSerialDenseVector &GIDs, const Epetra_SerialDenseVector &values, int vectorIndex=0)
 
int SumIntoGlobalValues (int numIDs, const int *GIDs, const int *numValuesPerID, const double *values, int vectorIndex=0)
 
int SumIntoGlobalValues (int numIDs, const long long *GIDs, const int *numValuesPerID, const double *values, int vectorIndex=0)
 
int ReplaceGlobalValues (int numIDs, const int *GIDs, const int *numValuesPerID, const double *values, int vectorIndex=0)
 
int ReplaceGlobalValues (int numIDs, const long long *GIDs, const int *numValuesPerID, const double *values, int vectorIndex=0)
 
int GlobalAssemble (Epetra_CombineMode mode=Add, bool reuse_map_and_exporter=false)
 
void setIgnoreNonLocalEntries (bool flag)
 
Epetra_FEVectoroperator= (const Epetra_FEVector &source)
 
- Public Member Functions inherited from Epetra_MultiVector
int ReplaceMap (const Epetra_BlockMap &map)
 
int Reduce ()
 
 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.
 
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...
 
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...
 
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...
 
int SetSeed (unsigned int Seed_in)
 Set seed for Random function. More...
 
unsigned int Seed ()
 Get seed from Random function. More...
 
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...
 
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.
 
virtual void Print (std::ostream &os) const
 Print method.
 
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_BlockMapMap () const
 Returns the address of the Epetra_BlockMap for this multi-vector.
 
const Epetra_CommComm () 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_CompObjectoperator= (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_FlopsGetFlopCounter () 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

template<typename int_type >
int inputValues (int numIDs, const int_type *GIDs, const double *values, bool suminto, int vectorIndex)
 
template<typename int_type >
int inputValues (int numIDs, const int_type *GIDs, const int *numValuesPerID, const double *values, bool suminto, int vectorIndex)
 
template<typename int_type >
int inputNonlocalValue (int_type GID, double value, bool suminto, int vectorIndex)
 
template<typename int_type >
int inputNonlocalValues (int_type GID, int numValues, const double *values, bool suminto, int vectorIndex)
 
template<typename int_type >
void createNonlocalMapAndExporter ()
 Allocate the Map, Export object, and MultiVector for nonlocal data.
 
void destroyNonlocalMapAndExporter ()
 Deallocate the Map, Export object, and MultiVector for nonlocal data.
 
template<typename int_type >
void zeroNonlocalData ()
 Make all the nonlocal multivector entries zero.
 
void destroyNonlocalData ()
 Deallocate storage for nonlocal data. More...
 
template<typename int_type >
std::vector< int_type > & nonlocalIDs ()
 
template<>
std::vector< int > & nonlocalIDs ()
 
template<>
std::vector< long long > & nonlocalIDs ()
 
- Protected Member Functions inherited from Epetra_MultiVector
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

long long myFirstID_
 
int myNumIDs_
 
std::vector< int > nonlocalIDs_int_
 
std::vector< long long > nonlocalIDs_LL_
 
std::vector< int > nonlocalElementSize_
 
std::vector< std::vector
< double > > 
nonlocalCoefs_
 Array of arrays (one per column) of nonlocal coefficients. More...
 
Epetra_BlockMapnonlocalMap_
 Map describing distribution of nonlocal data.
 
Epetra_Exportexporter_
 Export object that sums nonlocal data into a nonoverlapping distribution.
 
Epetra_MultiVectornonlocalVector_
 Multivector that holds nonlocal data; source for the Export operation.
 
bool ignoreNonLocalEntries_
 
- Protected Attributes inherited from Epetra_MultiVector
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_
 

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
 

Detailed Description

Epetra Finite-Element Vector. This class inherits Epetra_MultiVector and thus provides all Epetra_MultiVector functionality.

The added functionality provided by Epetra_FEVector is the ability to perform finite-element style vector assembly. It accepts sub-vector contributions, such as those that would come from element-load vectors, etc. These sub-vectors need not be owned by the local processor. In other words, the user can assemble overlapping data (e.g., corresponding to shared finite-element nodes). When the user is finished assembling their vector data, they then call the method Epetra_FEVector::GlobalAssemble() which gathers the overlapping data (all non-local data that was input on each processor) into the data-distribution specified by the map with which the Epetra_FEVector was constructed.

Constructor & Destructor Documentation

Epetra_FEVector::Epetra_FEVector ( const Epetra_BlockMap Map,
int  numVectors = 1,
bool  ignoreNonLocalEntries = false 
)

Constructor that requires a map specifying a non-overlapping data layout.

Parameters
MapMap describing a non-overlapping distribution for the underlying Epetra_MultiVector into which this Epetra_FEVector will funnel data.
numVectorsOptional argument, default value is 1. (See the documentation for Epetra_MultiVector for the meaning of this argument.
ignoreNonLocalEntriesOptional argument, default value is false. Under certain special circumstances it is desirable to have non-local contributions ignored rather than saving them for the GlobalAssemble step.
Epetra_FEVector::Epetra_FEVector ( Epetra_DataAccess  CV,
const Epetra_BlockMap Map,
double *  A,
int  MyLDA,
int  NumVectors,
bool  ignoreNonLocalEntries = false 
)

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.

Epetra_FEVector::Epetra_FEVector ( Epetra_DataAccess  CV,
const Epetra_BlockMap Map,
double **  ArrayOfPointers,
int  NumVectors,
bool  ignoreNonLocalEntries = false 
)

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.

Epetra_FEVector::Epetra_FEVector ( const Epetra_FEVector source)

Copy constructor.

virtual Epetra_FEVector::~Epetra_FEVector ( )
virtual

Destructor

Member Function Documentation

void Epetra_FEVector::destroyNonlocalData ( )
protected

Deallocate storage for nonlocal data.

This does not include whatever destroyNonlocalMapAndExporter() deallocates.

int Epetra_FEVector::GlobalAssemble ( Epetra_CombineMode  mode = Add,
bool  reuse_map_and_exporter = false 
)

Gather any overlapping/shared data into the non-overlapping partitioning defined by the Map that was passed to this vector at construction time. Data imported from other processors is stored on the owning processor with a "sumInto" or accumulate operation. This is a collective method – every processor must enter it before any will complete it.

Optimization for power-users: The optional parameter 'reuse_map_and_exporter' defaults to false. By default, a map that describes the non-local data is re-created at each call to GlobalAssemble, along with an exporter used to do the communication. This is expensive. If you know that the layout of your nonlocal data has not changed since your previous call to GlobalAssemble, you can set this flag to true and it will reuse the previously created map and exporter rather than creating new ones.

int Epetra_FEVector::ReplaceGlobalValues ( int  numIDs,
const int *  GIDs,
const double *  values,
int  vectorIndex = 0 
)

Copy values into the vector overwriting any values that already exist for the specified indices.

int Epetra_FEVector::ReplaceGlobalValues ( const Epetra_IntSerialDenseVector GIDs,
const Epetra_SerialDenseVector values,
int  vectorIndex = 0 
)

Copy values into the vector, replacing any values that already exist for the specified GIDs.

Parameters
GIDsList of global ids. Must be the same length as the accompanying list of values.
valuesList of coefficient values. Must be the same length as the accompanying list of GIDs.
void Epetra_FEVector::setIgnoreNonLocalEntries ( bool  flag)
inline

Set whether or not non-local data values should be ignored.

int Epetra_FEVector::SumIntoGlobalValues ( int  numIDs,
const int *  GIDs,
const double *  values,
int  vectorIndex = 0 
)

Accumulate values into the vector, adding them to any values that already exist for the specified indices.

int Epetra_FEVector::SumIntoGlobalValues ( const Epetra_IntSerialDenseVector GIDs,
const Epetra_SerialDenseVector values,
int  vectorIndex = 0 
)

Accumulate values into the vector, adding them to any values that already exist for the specified GIDs.

Parameters
GIDsList of global ids. Must be the same length as the accompanying list of values.
valuesList of coefficient values. Must be the same length as the accompanying list of GIDs.

Member Data Documentation

std::vector<std::vector<double> > Epetra_FEVector::nonlocalCoefs_
protected

Array of arrays (one per column) of nonlocal coefficients.

Entry j in this array is an array of nonlocal coefficients for column j of the multivector. nonlocalCoefs_[j] in turn contains the data for all the nonlocal block entries. Entries within a block are stored in column-major order (Fortran style), in contiguous chunks of Map().MaxElementSize() entries.

In summary: for k = 0 to numNonlocalIDs_ - 1, nonlocalCoefs_[j][j*M + i] contains the entries of nonlocal block entry k, for i = 0 to nonlocalElementSize_[k] - 1.


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