FEI  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | List of all members
fei::Matrix_Impl< T > Class Template Reference

#include <fei_Matrix_Impl.hpp>

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

Public Member Functions

 Matrix_Impl (fei::SharedPtr< T > matrix, fei::SharedPtr< fei::MatrixGraph > matrixGraph, int numLocalEqns, bool zeroSharedRows=true)
 
virtual ~Matrix_Impl ()
 
const char * typeName ()
 
int parameters (const fei::ParameterSet &paramset)
 
fei::SharedPtr< T > getMatrix ()
 
fei::SharedPtr< fei::MatrixGraphgetMatrixGraph () const
 
void setMatrixGraph (fei::SharedPtr< fei::MatrixGraph > matrixGraph)
 
int getGlobalNumRows () const
 
int getLocalNumRows () const
 
int putScalar (double scalar)
 
int getRowLength (int row, int &length) const
 
int copyOutRow (int row, int len, double *coefs, int *indices) const
 
int sumIn (int numRows, const int *rows, int numCols, const int *cols, const double *const *values, int format=0)
 
int copyIn (int numRows, const int *rows, int numCols, const int *cols, const double *const *values, int format=0)
 
int sumInFieldData (int fieldID, int idType, int rowID, int colID, const double *const *data, int format=0)
 
int sumInFieldData (int fieldID, int idType, int rowID, int colID, const double *data, int format=0)
 
int sumIn (int blockID, int connectivityID, const double *const *values, int format=0)
 
int globalAssemble ()
 
int multiply (fei::Vector *x, fei::Vector *y)
 
void setCommSizes ()
 
int gatherFromOverlap (bool accumulate=true)
 
int writeToFile (const char *filename, bool matrixMarketFormat=true)
 
int writeToStream (FEI_OSTREAM &ostrm, bool matrixMarketFormat=true)
 
bool usingBlockEntryStorage ()
 
int giveToUnderlyingMatrix (int numRows, const int *rows, int numCols, const int *cols, const double *const *values, bool sumInto, int format)
 
int giveToUnderlyingBlockMatrix (int row, int rowDim, int numCols, const int *cols, const int *LDAs, const int *colDims, const double *const *values, bool sumInto)
 
void markState ()
 
bool changedSinceMark ()
 
- Public Member Functions inherited from fei::Matrix
virtual ~Matrix ()
 
- Public Member Functions inherited from fei::Matrix_core
void setRHS (fei::SharedPtr< fei::Vector > rhsvector)
 
void setSlaveInfo (fei::SharedPtr< fei::MatrixGraph > matrixGraph)
 

Additional Inherited Members

- Static Public Member Functions inherited from fei::Matrix_core
static void copyTransposeToWorkArrays (int numRows, int numCols, const double *const *values, std::vector< double > &work_1D, std::vector< const double * > &work_2D)
 
- Protected Member Functions inherited from fei::Matrix_core
int copyPointRowsToBlockRow (int numPtRows, int numPtCols, const double *const *ptValues, int numBlkCols, const int *blkColDims, double **blkValues)
 
- Protected Member Functions inherited from fei::Logger
 Logger ()
 
virtual ~Logger ()
 
void setOutputLevel (OutputLevel olevel)
 
- Protected Attributes inherited from fei::Logger
OutputLevel output_level_
 
FEI_OSTREAM * output_stream_
 

Detailed Description

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

To be used for 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 above-mentioned gather 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) matrix.

Definition at line 53 of file fei_Matrix_Impl.hpp.

Constructor & Destructor Documentation

template<typename T >
fei::Matrix_Impl< T >::Matrix_Impl ( fei::SharedPtr< T >  matrix,
fei::SharedPtr< fei::MatrixGraph matrixGraph,
int  numLocalEqns,
bool  zeroSharedRows = true 
)

Constructor

Definition at line 436 of file fei_Matrix_Impl.hpp.

template<typename T >
fei::Matrix_Impl< T >::~Matrix_Impl ( )
virtual

Destructor

Definition at line 484 of file fei_Matrix_Impl.hpp.

Member Function Documentation

template<typename T>
const char* fei::Matrix_Impl< T >::typeName ( )
inlinevirtual

Return a name describing the run-time type of this object.

Implements fei::Matrix.

Definition at line 67 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::parameters ( const fei::ParameterSet paramset)
virtual

Parameters method

Implements fei::Matrix.

Definition at line 490 of file fei_Matrix_Impl.hpp.

template<typename T>
fei::SharedPtr<T> fei::Matrix_Impl< T >::getMatrix ( )
inline
Obtain the underlying matrix object. Note that this will generally be

only the locally-owned portion of the matrix, not including any data corresponding to shared-but-not-owned nodes, etc.

Definition at line 87 of file fei_Matrix_Impl.hpp.

template<typename T>
fei::SharedPtr<fei::MatrixGraph> fei::Matrix_Impl< T >::getMatrixGraph ( ) const
inlinevirtual

Obtain the fei::MatrixGraph associated with this matrix

Implements fei::Matrix.

Definition at line 90 of file fei_Matrix_Impl.hpp.

template<typename T >
void fei::Matrix_Impl< T >::setMatrixGraph ( fei::SharedPtr< fei::MatrixGraph matrixGraph)
inlinevirtual

Set the fei::MatrixGraph associated with this matrix

Implements fei::Matrix.

Definition at line 334 of file fei_Matrix_Impl.hpp.

template<typename T>
int fei::Matrix_Impl< T >::getGlobalNumRows ( ) const
inlinevirtual

Get the global number of rows in the matrix.

Implements fei::Matrix.

Definition at line 98 of file fei_Matrix_Impl.hpp.

template<typename T>
int fei::Matrix_Impl< T >::getLocalNumRows ( ) const
inlinevirtual

Get the local number of rows in the matrix.

Implements fei::Matrix.

Definition at line 106 of file fei_Matrix_Impl.hpp.

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

Set a specified scalar throughout the matrix.

Implements fei::Matrix.

Definition at line 523 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::getRowLength ( int  row,
int &  length 
) const
virtual

Get the length of a row of the matrix.

Parameters
rowGlobal 0-based equation number
lengthOutput. Length of the row.
Returns
error-code non-zero if any error occurs.

Implements fei::Matrix.

Definition at line 498 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::copyOutRow ( int  row,
int  len,
double *  coefs,
int *  indices 
) const
virtual

Obtain a copy of the coefficients and indices for a row of the matrix.

Parameters
rowGlobal 0-based equation number
lenLength of the caller-allocated coefs and indices arrays
coefsCaller-allocated array, length 'len', to be filled with coefficients
indicesCaller-allocated array, length 'len', to be filled with indices. (These indices will be global 0-based equation numbers.)
Returns
error-code non-zero if any error occurs.

Implements fei::Matrix.

Definition at line 552 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::sumIn ( int  numRows,
const int *  rows,
int  numCols,
const int *  cols,
const double *const *  values,
int  format = 0 
)
virtual

Sum coefficients into the matrix, adding them to any coefficients that may already exist at the specified row/column locations.

Parameters
numRows
rows
numCols
cols
values
formatFor compatibility with old FEI elemFormat... 0 means row-wise or row-major, 3 means column-major. Others not recognized

Implements fei::Matrix.

Definition at line 587 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::copyIn ( int  numRows,
const int *  rows,
int  numCols,
const int *  cols,
const double *const *  values,
int  format = 0 
)
virtual

Copy coefficients into the matrix, overwriting any coefficients that may already exist at the specified row/column locations.

Parameters
numRows
rows
numCols
cols
values
formatFor compatibility with old FEI elemFormat... 0 means row-wise or row-major, 3 means column-major. Others not recognized

Implements fei::Matrix.

Definition at line 610 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::sumInFieldData ( int  fieldID,
int  idType,
int  rowID,
int  colID,
const double *const *  data,
int  format = 0 
)
virtual

Sum coefficients into the matrix, specifying row/column locations by identifier/fieldID pairs.

Parameters
fieldIDInput. field-identifier for which data is being input.
idTypeInput. The identifier-type of the identifiers.
rowIDInput. Identifier in row-space, for which data is being input.
colIDInput. Identifier in column-space, for which data is being input.
dataInput. C-style table of data. num-rows is the field-size (i.e., number of scalar components that make up the field) of 'fieldID', as is num-columns.
formatFor compatibility with old FEI elemFormat... 0 means row-wise or row-major, 3 means column-major. Others not recognized
Returns
error-code 0 if successful

Implements fei::Matrix.

Definition at line 632 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::sumInFieldData ( int  fieldID,
int  idType,
int  rowID,
int  colID,
const double *  data,
int  format = 0 
)
virtual

Sum coefficients into the matrix, specifying row/column locations by identifier/fieldID pairs.

Parameters
fieldIDInput. field-identifier for which data is being input.
idTypeInput. The identifier-type of the identifiers.
rowIDInput. Identifier in row-space, for which data is being input.
colIDInput. Identifier in column-space, for which data is being input.
dataInput. 1-D list representing a packed table of data. Data may be backed in row-major or column-major order and this may be specified with the 'format' argument. The "table" of data is of size num-rows X num-columns and num-rows is the field-size (i.e., number of scalar components that make up the field) of 'fieldID', as is num-columns.
formatFor compatibility with old FEI elemFormat... 0 means row-wise or row-major, 3 means column-major. Others not recognized
Returns
error-code 0 if successful

Implements fei::Matrix.

Definition at line 668 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::sumIn ( int  blockID,
int  connectivityID,
const double *const *  values,
int  format = 0 
)
virtual

Sum coefficients, associated with a connectivity-block that was initialized on the MatrixGraph object, into this matrix.

Parameters
blockID
connectivityID
values
formatFor compatibility with old FEI elemFormat... 0 means row-wise or row-major, 3 means column-major. Others not recognized

Implements fei::Matrix.

Definition at line 695 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::globalAssemble ( )
inlinevirtual

Perform any necessary internal communications/synchronizations or other operations appropriate at end of data input. For some implementations this will be a no-op.

Implements fei::Matrix.

Definition at line 317 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::multiply ( fei::Vector x,
fei::Vector y 
)
virtual

Form a matrix-vector product y = 'this' * x

Implements fei::Matrix.

Definition at line 926 of file fei_Matrix_Impl.hpp.

template<typename T >
void fei::Matrix_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::Matrix.

Definition at line 934 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::gatherFromOverlap ( bool  accumulate = true)
virtual

After local overlapping data has been input, (e.g., element-data for a finite-element application) call this method to have data that corresponds to shared identifiers be communicated from sharing-but-not- owning processors, to owning processors.

Implements fei::Matrix.

Definition at line 945 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::writeToFile ( const char *  filename,
bool  matrixMarketFormat = true 
)
virtual

Implementation of fei::Matrix::writeToFile

Implements fei::Matrix.

Definition at line 1230 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::writeToStream ( FEI_OSTREAM &  ostrm,
bool  matrixMarketFormat = true 
)
virtual

Implementation of fei::Matrix::writeToStream

Implements fei::Matrix.

Definition at line 1328 of file fei_Matrix_Impl.hpp.

template<typename T>
bool fei::Matrix_Impl< T >::usingBlockEntryStorage ( )
inlinevirtual

Query whether the underlying matrix object is a block-entry matrix.

Implements fei::Matrix.

Definition at line 252 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::giveToUnderlyingMatrix ( int  numRows,
const int *  rows,
int  numCols,
const int *  cols,
const double *const *  values,
bool  sumInto,
int  format 
)

for experts only

Definition at line 356 of file fei_Matrix_Impl.hpp.

template<typename T >
int fei::Matrix_Impl< T >::giveToUnderlyingBlockMatrix ( int  row,
int  rowDim,
int  numCols,
const int *  cols,
const int *  LDAs,
const int *  colDims,
const double *const *  values,
bool  sumInto 
)

for experts only

Definition at line 389 of file fei_Matrix_Impl.hpp.

template<typename T >
void fei::Matrix_Impl< T >::markState ( )
inlinevirtual

Set a "mark" point on the current state of the matrix, so that later a query can be made to see if the matrix has changed since this mark was set.

Implements fei::Matrix.

Definition at line 341 of file fei_Matrix_Impl.hpp.

template<typename T >
bool fei::Matrix_Impl< T >::changedSinceMark ( )
inlinevirtual

Query whether the matrix has changed since markState() was called. If markState() hasn't been called since the matrix was constructed, then this query will return true.

Implements fei::Matrix.

Definition at line 349 of file fei_Matrix_Impl.hpp.


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