FEI Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
fei::Vector_Local Class Reference

#include <fei_Vector_Local.hpp>

Inheritance diagram for fei::Vector_Local:
Inheritance graph
[legend]

Public Member Functions

 Vector_Local (fei::SharedPtr< fei::VectorSpace > vecSpace)
 
virtual ~Vector_Local ()
 
const char * typeName () const
 
int update (double a, const fei::Vector *x, double b)
 
int scatterToOverlap ()
 
void setCommSizes ()
 
int gatherFromOverlap (bool accumulate=true)
 
int putScalar (double scalar)
 
int sumIn (int numValues, const int *indices, const double *values, int vectorIndex=0)
 
int copyIn (int numValues, const int *indices, const double *values, int vectorIndex=0)
 
fei::SharedPtr< fei::VectorSpacegetVectorSpace () const
 
void setVectorSpace (fei::SharedPtr< fei::VectorSpace > vecSpace)
 
int sumInFieldData (int fieldID, int idType, int numIDs, const int *IDs, const double *data, int vectorIndex=0)
 
int copyInFieldData (int fieldID, int idType, int numIDs, const int *IDs, const double *data, int vectorIndex=0)
 
int copyInFieldDataLocalIDs (int fieldID, int idType, int numIDs, const int *localIDs, const double *data, int vectorIndex=0)
 
int copyOutFieldData (int fieldID, int idType, int numIDs, const int *IDs, double *data, int vectorIndex=0)
 
int copyOut (int numValues, const int *indices, double *values, int vectorIndex=0) const
 
int writeToFile (const char *filename, bool matrixMarketFormat=true)
 
int writeToStream (FEI_OSTREAM &ostrm, bool matrixMarketFormat=true)
 
std::vector< double > & getCoefs ()
 
- Public Member Functions inherited from fei::Vector
virtual ~Vector ()
 

Private Member Functions

int giveToVector (int numValues, const int *indices, const double *values, bool sumInto, int vectorIndex)
 
int assembleFieldData (int fieldID, int idType, int numIDs, const int *IDs, const double *data, bool sumInto, int vectorIndex)
 
int assembleFieldDataLocalIDs (int fieldID, int idType, int numIDs, const int *localIDs, const double *data, bool sumInto, int vectorIndex)
 

Private Attributes

fei::SharedPtr< fei::VectorSpacevecSpace_
 
std::vector< double > coefs_
 
std::map< int, int > global_to_local_
 
std::vector< int > work_indices_
 

Detailed Description

Definition at line 19 of file fei_Vector_Local.hpp.

Constructor & Destructor Documentation

fei::Vector_Local::Vector_Local ( fei::SharedPtr< fei::VectorSpace vecSpace)
fei::Vector_Local::~Vector_Local ( )
virtual

Definition at line 40 of file fei_Vector_Local.cpp.

Member Function Documentation

const char* fei::Vector_Local::typeName ( ) const
inlinevirtual
  Return an implementation-dependent name describing the run-time type

of this object.

Implements fei::Vector.

Definition at line 25 of file fei_Vector_Local.hpp.

int fei::Vector_Local::update ( double  a,
const fei::Vector x,
double  b 
)
virtual

Update 'this' = b*'this' + a*x

Implements fei::Vector.

Definition at line 45 of file fei_Vector_Local.cpp.

References fei::console_out(), and FEI_ENDL.

int fei::Vector_Local::scatterToOverlap ( )
virtual

Use data in the underlying non-overlapping decomposition to update any shared data in the overlapping decomposition.

If any data is already held for the shared positions, that data will be replaced by the data from the 'owning' processor.

Implements fei::Vector.

Definition at line 54 of file fei_Vector_Local.cpp.

void fei::Vector_Local::setCommSizes ( )
virtual

perform initial communication to establish message sizes that will be needed for exchanging shared-node data. Called from within gatherFromOverlap usually, doesn't usually need to be explicitly called by client code. (Power users only...)

Implements fei::Vector.

Definition at line 58 of file fei_Vector_Local.cpp.

int fei::Vector_Local::gatherFromOverlap ( bool  accumulate = true)
virtual

Move any shared data from the overlapping decomposition to the underlying non-overlapping decomposition.

Implements fei::Vector.

Definition at line 63 of file fei_Vector_Local.cpp.

int fei::Vector_Local::putScalar ( double  scalar)
virtual

Set a specified scalar throughout the vector.

Implements fei::Vector.

Definition at line 67 of file fei_Vector_Local.cpp.

References coefs_.

int fei::Vector_Local::sumIn ( int  numValues,
const int *  indices,
const double *  values,
int  vectorIndex = 0 
)
virtual

Sum values into the vector, adding to any that may already exist at the specified indices.

Implements fei::Vector.

Definition at line 104 of file fei_Vector_Local.cpp.

References giveToVector().

int fei::Vector_Local::copyIn ( int  numValues,
const int *  indices,
const double *  values,
int  vectorIndex = 0 
)
virtual

Copy values into the vector, overwriting any that may already exist at the specified indices.

Implements fei::Vector.

Definition at line 111 of file fei_Vector_Local.cpp.

References giveToVector().

fei::SharedPtr< fei::VectorSpace > fei::Vector_Local::getVectorSpace ( ) const
virtual

Obtain the VectorSpace associated with this vector.

Implements fei::Vector.

Definition at line 118 of file fei_Vector_Local.cpp.

References vecSpace_.

void fei::Vector_Local::setVectorSpace ( fei::SharedPtr< fei::VectorSpace vecSpace)
virtual

Set the VectorSpace associated with this vector.

Implements fei::Vector.

Definition at line 122 of file fei_Vector_Local.cpp.

References vecSpace_.

int fei::Vector_Local::sumInFieldData ( int  fieldID,
int  idType,
int  numIDs,
const int *  IDs,
const double *  data,
int  vectorIndex = 0 
)
virtual

Sum field data into the vector, adding to any coefficients that may already exist at the specified locations. If the specified fieldID doesn't exist at one or more of the specified IDs, then the corresponding positions in the data array will simply not be used.

Implements fei::Vector.

Definition at line 170 of file fei_Vector_Local.cpp.

References assembleFieldData().

int fei::Vector_Local::copyInFieldData ( int  fieldID,
int  idType,
int  numIDs,
const int *  IDs,
const double *  data,
int  vectorIndex = 0 
)
virtual

Copy field data into the vector, overwriting any coefficients that may already exist at the specified locations. If the specified fieldID doesn't exist at one or more of the specified IDs, then the corresponding positions in the data array will simply not be used.

Implements fei::Vector.

Definition at line 182 of file fei_Vector_Local.cpp.

References assembleFieldData().

int fei::Vector_Local::copyInFieldDataLocalIDs ( int  fieldID,
int  idType,
int  numIDs,
const int *  localIDs,
const double *  data,
int  vectorIndex = 0 
)
virtual

Implements fei::Vector.

Definition at line 194 of file fei_Vector_Local.cpp.

References assembleFieldDataLocalIDs().

int fei::Vector_Local::copyOutFieldData ( int  fieldID,
int  idType,
int  numIDs,
const int *  IDs,
double *  data,
int  vectorIndex = 0 
)
virtual

Copy field data out of the vector, into the caller-allocated data array. If the specified fieldID doesn't exist at one or more of the specified IDs, then the corresponding positions in the data array will simply not be referenced.

Implements fei::Vector.

Definition at line 206 of file fei_Vector_Local.cpp.

References CHK_ERR, coefs_, fei::console_out(), FEI_ENDL, fei::VectorSpace::getFieldSize(), fei::VectorSpace::getGlobalIndices(), global_to_local_, vecSpace_, and work_indices_.

int fei::Vector_Local::copyOut ( int  numValues,
const int *  indices,
double *  values,
int  vectorIndex = 0 
) const
virtual
  Retrieve a copy of values from the vector for the specified indices.

Note that if the specified indices are not local in the underlying non-overlapping data decomposition, these values are not guaranteed to be correct until after the scatterToOverlap() method has been called.

Implements fei::Vector.

Definition at line 237 of file fei_Vector_Local.cpp.

References coefs_, fei::console_out(), FEI_ENDL, and global_to_local_.

int fei::Vector_Local::writeToFile ( const char *  filename,
bool  matrixMarketFormat = true 
)
virtual
  Write the vector's contents into the specified file.
Parameters
filenameText name of the file to be created or overwritten. If in a parallel environment, each processor will take turns writing into the file.
matrixMarketFormatOptional argument, defaults to true. If true the contents of the file will be MatrixMarket real array format. If not true, the contents of the file will contain the vector's global dimension on the first line, and all following lines will contain a space-separated pair with global index first and coefficient value second.
Returns
error-code 0 if successful, -1 if some error occurs such as failure to open file.

Implements fei::Vector.

Definition at line 267 of file fei_Vector_Local.cpp.

References FEI_OFSTREAM, FEI_OSTRINGSTREAM, fei::VectorSpace::getCommunicator(), IOS_OUT, fei::localProc(), vecSpace_, and writeToStream().

int fei::Vector_Local::writeToStream ( FEI_OSTREAM ostrm,
bool  matrixMarketFormat = true 
)
virtual
  Write the vector's contents to the specified ostream.
Parameters
ostrmostream to be written to.
matrixMarketFormatOptional argument, defaults to true. If true the contents of the vector will be written in MatrixMarket real array format. If not true, the stream will be given the vector's global dimension on the first line, and all following lines will contain a space-separated pair with global index first and coefficient value second.
Returns
error-code 0 if successful, -1 if some error occurs such as failure to open file.

Implements fei::Vector.

Definition at line 280 of file fei_Vector_Local.cpp.

References coefs_, FEI_ENDL, global_to_local_, IOS_FLOATFIELD, and IOS_SCIENTIFIC.

Referenced by writeToFile().

std::vector< double > & fei::Vector_Local::getCoefs ( )

Definition at line 261 of file fei_Vector_Local.cpp.

References coefs_.

int fei::Vector_Local::giveToVector ( int  numValues,
const int *  indices,
const double *  values,
bool  sumInto,
int  vectorIndex 
)
private
int fei::Vector_Local::assembleFieldData ( int  fieldID,
int  idType,
int  numIDs,
const int *  IDs,
const double *  data,
bool  sumInto,
int  vectorIndex 
)
private
int fei::Vector_Local::assembleFieldDataLocalIDs ( int  fieldID,
int  idType,
int  numIDs,
const int *  localIDs,
const double *  data,
bool  sumInto,
int  vectorIndex 
)
private

Member Data Documentation

fei::SharedPtr<fei::VectorSpace> fei::Vector_Local::vecSpace_
private
std::vector<double> fei::Vector_Local::coefs_
private
std::map<int,int> fei::Vector_Local::global_to_local_
private
std::vector<int> fei::Vector_Local::work_indices_
private

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