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_Impl< T > Class Template Reference

#include <fei_Vector_Impl.hpp>

Inheritance diagram for fei::Vector_Impl< T >:
Inheritance graph
[legend]

Public Member Functions

 Vector_Impl (fei::SharedPtr< fei::VectorSpace > vecSpace, T *vector, int numLocalEqns, bool isSolutionVector=false, bool deleteVector=false)
 
virtual ~Vector_Impl ()
 
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 writeToFile (const char *filename, bool matrixMarketFormat=true)
 
int writeToStream (FEI_OSTREAM &ostrm, bool matrixMarketFormat=true)
 
void setUnderlyingVector (T *vec)
 
T * getUnderlyingVector ()
 
const T * getUnderlyingVector () const
 
int copyOut (int numValues, const int *indices, double *values, int vectorIndex=0) const
 
int copyOut_FE (int nodeNumber, int dofOffset, double &value)
 
- Public Member Functions inherited from fei::Vector
virtual ~Vector ()
 
- Public Member Functions inherited from fei::Vector_core
 Vector_core (fei::SharedPtr< fei::VectorSpace > vecSpace, int numLocalEqns)
 
virtual ~Vector_core ()
 
int copyOut (int numValues, const int *indices, double *values, int vectorIndex=0) const
 
int giveToVector (int numValues, const int *indices, const double *values, bool sumInto=true, int vectorIndex=0)
 
void setOverlap (int numRemoteEqns=0, const int *remoteEqns=NULL)
 

Private Member Functions

int giveToUnderlyingVector (int numValues, const int *indices, const double *values, bool sumInto=true, int vectorIndex=0)
 
int copyOutOfUnderlyingVector (int numValues, const int *indices, double *values, int vectorIndex=0) const
 
int sumIntoFEVector (int blockID, int connOffset, int numNodes, const int *nodeNumbers, const int *numIndicesPerNode, const int *dof_ids, const double *values)
 

Private Attributes

T * vector_
 
bool isSolution_
 
bool deleteVector_
 
int localProc_
 
int numProcs_
 
std::string dbgprefix_
 

Additional Inherited Members

- Protected Member Functions inherited from fei::Vector_core
int assembleFieldData (int fieldID, int idType, int numIDs, const int *IDs, const double *data, bool sumInto=true, int vectorIndex=0)
 
int assembleFieldDataLocalIDs (int fieldID, int idType, int numIDs, const int *localIDs, const double *data, bool sumInto=true, int vectorIndex=0)
 
void setCommSizes ()
 
fei::SharedPtr< fei::VectorSpaceget_vector_space () const
 
void set_vector_space (fei::SharedPtr< fei::VectorSpace > vspace)
 
int firstLocalOffset () const
 
int lastLocalOffset () const
 
std::vector< int > & work_indices ()
 
std::vector< int > & work_indices2 ()
 
bool haveFEVector ()
 
void setFEVector (bool flag)
 
std::vector< CSVec * > & remotelyOwned ()
 
const std::vector< CSVec * > & remotelyOwned () const
 
std::vector< int > & remotelyOwnedProcs ()
 
const std::vector< int > & remotelyOwnedProcs () const
 
fei::CSVecgetRemotelyOwned (int proc)
 
const fei::CSVecgetRemotelyOwned (int proc) const
 
- Protected Member Functions inherited from fei::Logger
 Logger ()
 
virtual ~Logger ()
 
void setOutputLevel (OutputLevel olevel)
 
void addLogID (int ID)
 
void addLogEqn (int eqn)
 
bool isLogID (int ID)
 
bool isLogEqn (int eqn)
 
std::set< int > & getLogIDs ()
 
std::set< int > & getLogEqns ()
 
- Protected Attributes inherited from fei::Vector_core
fei::SharedPtr< fei::EqnCommeqnComm_
 
- Protected Attributes inherited from fei::Logger
OutputLevel output_level_
 
FEI_OSTREAMoutput_stream_
 
std::set< int > logIDs_
 
std::set< int > logEqns_
 

Detailed Description

template<typename T>
class fei::Vector_Impl< T >

To be used for vector data assembly, including local assembly of shared data. Provides operations for gathering the overlapped data (locally- stored shared data) to a non-overlapped data distribution (e.g., send shared data to owning processor) and vice-versa for scattering non- overlapped data to the overlapped distribution.

When shared data that is not locally-owned is assembled into this object, it will be held locally until the gatherFromOverlap() operation is performed. When data that is locally-owned is assembled into this object, it will be passed directly to the underlying algebraic (non-overlapping) vector.

When non-locally-owned shared data is requested from this vector object, the operation is only guaranteed to succeed if scatterToOverlap() has been called. If non-local data has recently been input to this vector, but the vector has not been 'synchronized' using gatherFromOverlap() and scatterToOverlap(), then that data may be used to answer the request but the values may be erroneous due to not including contributions from other processors.

Definition at line 57 of file fei_Vector_Impl.hpp.

Constructor & Destructor Documentation

template<typename T >
fei::Vector_Impl< T >::Vector_Impl ( fei::SharedPtr< fei::VectorSpace vecSpace,
T *  vector,
int  numLocalEqns,
bool  isSolutionVector = false,
bool  deleteVector = false 
)
template<typename T >
fei::Vector_Impl< T >::~Vector_Impl ( )
virtual

Destructor

Definition at line 294 of file fei_Vector_Impl.hpp.

Member Function Documentation

template<typename T>
const char* fei::Vector_Impl< T >::typeName ( ) const
inlinevirtual
  Return a name describing the run-time type

of this object.

Implements fei::Vector.

Definition at line 72 of file fei_Vector_Impl.hpp.

Referenced by fei::Vector_Impl< T >::Vector_Impl().

template<typename T >
int fei::Vector_Impl< T >::update ( double  a,
const fei::Vector x,
double  b 
)
virtual

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

Implements fei::Vector.

Definition at line 323 of file fei_Vector_Impl.hpp.

References fei::Vector_Impl< T >::getUnderlyingVector().

template<typename T >
int fei::Vector_Impl< T >::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 339 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, FEI_ENDL, FEI_OSTREAM, and fei::Vector_core::scatterToOverlap().

template<typename T >
void fei::Vector_Impl< T >::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 351 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, FEI_ENDL, FEI_OSTREAM, and fei::Vector_core::setCommSizes().

Referenced by fei::Vector_Impl< T >::Vector_Impl().

template<typename T >
int fei::Vector_Impl< T >::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 363 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, FEI_ENDL, FEI_OSTREAM, and fei::Vector_core::gatherFromOverlap().

template<typename T >
int fei::Vector_Impl< T >::putScalar ( double  scalar)
virtual

Set a specified scalar throughout the vector.

Implements fei::Vector.

Definition at line 301 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, CHK_ERR, FEI_ENDL, FEI_OSTREAM, and fei::set_values().

template<typename T >
int fei::Vector_Impl< T >::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 375 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, FEI_ENDL, and FEI_OSTREAM.

Referenced by fei::Vector_Impl< T >::Vector_Impl().

template<typename T >
int fei::Vector_Impl< T >::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 389 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, FEI_ENDL, and FEI_OSTREAM.

Referenced by fei::Vector_Impl< T >::Vector_Impl().

template<typename T>
fei::SharedPtr<fei::VectorSpace> fei::Vector_Impl< T >::getVectorSpace ( ) const
inlinevirtual

Obtain the VectorSpace associated with this vector.

Implements fei::Vector.

Definition at line 112 of file fei_Vector_Impl.hpp.

References fei::Vector_core::get_vector_space().

template<typename T>
void fei::Vector_Impl< T >::setVectorSpace ( fei::SharedPtr< fei::VectorSpace vecSpace)
inlinevirtual

Set the VectorSpace associated with this vector.

Implements fei::Vector.

Definition at line 117 of file fei_Vector_Impl.hpp.

References fei::Vector_core::set_vector_space().

template<typename T >
int fei::Vector_Impl< T >::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 446 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, FEI_ENDL, and FEI_OSTREAM.

template<typename T >
int fei::Vector_Impl< T >::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 463 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, FEI_ENDL, and FEI_OSTREAM.

template<typename T >
int fei::Vector_Impl< T >::copyInFieldDataLocalIDs ( int  fieldID,
int  idType,
int  numIDs,
const int *  localIDs,
const double *  data,
int  vectorIndex = 0 
)
virtual

Implements fei::Vector.

Definition at line 480 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, FEI_ENDL, and FEI_OSTREAM.

template<typename T >
int fei::Vector_Impl< T >::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 504 of file fei_Vector_Impl.hpp.

References fei::Vector_core::copyOutFieldData().

template<typename T >
int fei::Vector_Impl< T >::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 517 of file fei_Vector_Impl.hpp.

References fei::Vector_core::writeToFile().

template<typename T >
int fei::Vector_Impl< T >::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 525 of file fei_Vector_Impl.hpp.

References fei::Vector_core::writeToStream().

template<typename T>
void fei::Vector_Impl< T >::setUnderlyingVector ( T *  vec)
inline
  Set the library-specific underlying vector object that this

snl_fei::Vector is filtering data in and out of.

Definition at line 177 of file fei_Vector_Impl.hpp.

References fei::Vector_Impl< T >::vector_.

template<typename T>
T* fei::Vector_Impl< T >::getUnderlyingVector ( )
inline
  Get the library-specific underlying vector object that this

snl_fei::Vector is filtering data in and out of.

Definition at line 185 of file fei_Vector_Impl.hpp.

References fei::Vector_Impl< T >::vector_.

Referenced by fei::MatrixTraits< FillableMat >::matvec(), and fei::Vector_Impl< T >::update().

template<typename T>
const T* fei::Vector_Impl< T >::getUnderlyingVector ( ) const
inline

Definition at line 186 of file fei_Vector_Impl.hpp.

References fei::Vector_Impl< T >::vector_.

template<typename T >
int fei::Vector_Impl< T >::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 533 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, fei::Vector_core::copyOut(), FEI_ENDL, and FEI_OSTREAM.

template<typename T >
int fei::Vector_Impl< T >::copyOut_FE ( int  nodeNumber,
int  dofOffset,
double &  value 
)
virtual

please ignore

Implements fei::Vector_core.

Definition at line 497 of file fei_Vector_Impl.hpp.

template<typename T >
int fei::Vector_Impl< T >::giveToUnderlyingVector ( int  numValues,
const int *  indices,
const double *  values,
bool  sumInto = true,
int  vectorIndex = 0 
)
privatevirtual

Review this function. Is it redundant with other functions?

Implements fei::Vector_core.

Definition at line 403 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, FEI_ENDL, FEI_OSTREAM, and fei::VectorTraits< T >::putValuesIn().

template<typename T >
int fei::Vector_Impl< T >::copyOutOfUnderlyingVector ( int  numValues,
const int *  indices,
double *  values,
int  vectorIndex = 0 
) const
privatevirtual

Review this function. Is it redundant with other functions?

Implements fei::Vector_core.

Definition at line 429 of file fei_Vector_Impl.hpp.

References fei::BRIEF_LOGS, FEI_ENDL, and FEI_OSTREAM.

template<typename T >
int fei::Vector_Impl< T >::sumIntoFEVector ( int  blockID,
int  connOffset,
int  numNodes,
const int *  nodeNumbers,
const int *  numIndicesPerNode,
const int *  dof_ids,
const double *  values 
)
privatevirtual

Sum in data for FiniteElementData-specific structure. Power users only.

Implements fei::Vector_core.

Definition at line 548 of file fei_Vector_Impl.hpp.

Member Data Documentation

template<typename T>
T* fei::Vector_Impl< T >::vector_
private
template<typename T>
bool fei::Vector_Impl< T >::isSolution_
private

Definition at line 218 of file fei_Vector_Impl.hpp.

template<typename T>
bool fei::Vector_Impl< T >::deleteVector_
private

Definition at line 219 of file fei_Vector_Impl.hpp.

template<typename T>
int fei::Vector_Impl< T >::localProc_
private

Definition at line 221 of file fei_Vector_Impl.hpp.

Referenced by fei::Vector_Impl< T >::Vector_Impl().

template<typename T>
int fei::Vector_Impl< T >::numProcs_
private

Definition at line 222 of file fei_Vector_Impl.hpp.

Referenced by fei::Vector_Impl< T >::Vector_Impl().

template<typename T>
std::string fei::Vector_Impl< T >::dbgprefix_
private

Definition at line 223 of file fei_Vector_Impl.hpp.

Referenced by fei::Vector_Impl< T >::Vector_Impl().


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