FEI  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Classes | Public Member Functions | List of all members
fei::Matrix Class Referenceabstract

#include <fei_Matrix.hpp>

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

Classes

class  Factory
 

Public Member Functions

virtual ~Matrix ()
 
virtual const char * typeName ()=0
 
virtual int parameters (const fei::ParameterSet &paramset)=0
 
virtual fei::SharedPtr
< fei::MatrixGraph
getMatrixGraph () const =0
 
virtual void setMatrixGraph (fei::SharedPtr< fei::MatrixGraph > matrixGraph)=0
 
virtual int getGlobalNumRows () const =0
 
virtual int getLocalNumRows () const =0
 
virtual int getRowLength (int row, int &length) const =0
 
virtual int putScalar (double scalar)=0
 
virtual int copyOutRow (int row, int len, double *coefs, int *indices) const =0
 
virtual int sumIn (int numRows, const int *rows, int numCols, const int *cols, const double *const *values, int format=0)=0
 
virtual int copyIn (int numRows, const int *rows, int numCols, const int *cols, const double *const *values, int format=0)=0
 
virtual int sumInFieldData (int fieldID, int idType, int rowID, int colID, const double *const *data, int format=0)=0
 
virtual int sumInFieldData (int fieldID, int idType, int rowID, int colID, const double *data, int format=0)=0
 
virtual int sumIn (int blockID, int connectivityID, const double *const *values, int format=0)=0
 
virtual int globalAssemble ()=0
 
virtual int multiply (fei::Vector *x, fei::Vector *y)=0
 
virtual void setCommSizes ()=0
 
virtual int gatherFromOverlap (bool accumulate=true)=0
 
virtual int writeToFile (const char *filename, bool matrixMarketFormat=true)=0
 
virtual int writeToStream (FEI_OSTREAM &ostrm, bool matrixMarketFormat=true)=0
 
virtual bool usingBlockEntryStorage ()=0
 
virtual void markState ()=0
 
virtual bool changedSinceMark ()=0
 

Detailed Description

Abstract representation of an algebraic matrix. This representation does not require that data be accessed only on the 'owning' processor. In other words, this representation may be used with an overlapping data decomposition. In most cases the underlying library-specific matrix will have a non-overlapping data decomposition (each equation uniquely owned by a single processor). Overlapping data may be assembled into this abstract matrix locally, and will be funneled into the underlying non- overlapping matrix on the correct processor when the gatherFromOverlap() method is called. Conversely, if the user wants to retrieve overlapping data from the matrix locally, that data is not guaranteed to be available until the scatterToOverlap() method is called.

Definition at line 30 of file fei_Matrix.hpp.

Constructor & Destructor Documentation

virtual fei::Matrix::~Matrix ( )
inlinevirtual

Virtual destructor.

Definition at line 44 of file fei_Matrix.hpp.

Member Function Documentation

virtual const char* fei::Matrix::typeName ( )
pure virtual

Return an implementation-dependent name describing the run-time type of this object.

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::parameters ( const fei::ParameterSet paramset)
pure virtual

Method for supplying parameters

Implemented in fei::Matrix_Impl< T >.

virtual fei::SharedPtr<fei::MatrixGraph> fei::Matrix::getMatrixGraph ( ) const
pure virtual

Obtain the fei::MatrixGraph associated with this matrix

Implemented in fei::Matrix_Impl< T >.

virtual void fei::Matrix::setMatrixGraph ( fei::SharedPtr< fei::MatrixGraph matrixGraph)
pure virtual

Set the fei::MatrixGraph associated with this matrix

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::getGlobalNumRows ( ) const
pure virtual

Get the global number of rows in the matrix.

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::getLocalNumRows ( ) const
pure virtual

Get the local number of rows in the matrix.

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::getRowLength ( int  row,
int &  length 
) const
pure 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.

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::putScalar ( double  scalar)
pure virtual

Set a specified scalar throughout the matrix.

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::copyOutRow ( int  row,
int  len,
double *  coefs,
int *  indices 
) const
pure virtual

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

Parameters
rowGlobal 0-based equation number
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.)
lenLength of the caller-allocated coefs and indices arrays
Returns
error-code non-zero if any error occurs.

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::sumIn ( int  numRows,
const int *  rows,
int  numCols,
const int *  cols,
const double *const *  values,
int  format = 0 
)
pure 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

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::copyIn ( int  numRows,
const int *  rows,
int  numCols,
const int *  cols,
const double *const *  values,
int  format = 0 
)
pure 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

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::sumInFieldData ( int  fieldID,
int  idType,
int  rowID,
int  colID,
const double *const *  data,
int  format = 0 
)
pure 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

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::sumInFieldData ( int  fieldID,
int  idType,
int  rowID,
int  colID,
const double *  data,
int  format = 0 
)
pure 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

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::sumIn ( int  blockID,
int  connectivityID,
const double *const *  values,
int  format = 0 
)
pure 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

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::globalAssemble ( )
pure virtual

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

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::multiply ( fei::Vector x,
fei::Vector y 
)
pure virtual

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

Implemented in fei::Matrix_Impl< T >.

virtual void fei::Matrix::setCommSizes ( )
pure 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...)

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::gatherFromOverlap ( bool  accumulate = true)
pure 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.

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::writeToFile ( const char *  filename,
bool  matrixMarketFormat = true 
)
pure virtual

Write the matrix 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 matrix global dimensions on the first line, and all following lines will contain a space-separated triple with global row index first, global column index second and coefficient value third. Note also that if matrixMarketFormat is true, indices will be output in 1-based form, but if not true, indices will be 0-based.
Returns
error-code 0 if successful, -1 if some error occurs such as failure to open file.

Implemented in fei::Matrix_Impl< T >.

virtual int fei::Matrix::writeToStream ( FEI_OSTREAM &  ostrm,
bool  matrixMarketFormat = true 
)
pure virtual

Write the matrix contents into the specified ostream.

Parameters
ostrmostream to be written to.
matrixMarketFormatOptional argument, defaults to true. If true the data will be written in MatrixMarket real array format. If not true, the stream will receive the matrix global dimensions on the first line, and all following lines will contain a space-separated triple with global row index first, global column index second and coefficient value third. Note also that if matrixMarketFormat is true, indices will be output in 1-based form, but if not true, indices will be 0-based.
Returns
error-code 0 if successful, -1 if some error occurs.

Implemented in fei::Matrix_Impl< T >.

virtual bool fei::Matrix::usingBlockEntryStorage ( )
pure virtual

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

Implemented in fei::Matrix_Impl< T >.

virtual void fei::Matrix::markState ( )
pure virtual

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.

Implemented in fei::Matrix_Impl< T >.

virtual bool fei::Matrix::changedSinceMark ( )
pure virtual

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.

Implemented in fei::Matrix_Impl< T >.


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