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 | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
fei::Matrix_Local Class Reference

#include <fei_Matrix_Local.hpp>

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

Public Member Functions

 Matrix_Local (fei::SharedPtr< fei::MatrixGraph > matrixGraph, fei::SharedPtr< fei::SparseRowGraph > sparseRowGraph)
 
virtual ~Matrix_Local ()
 
const char * typeName ()
 
int parameters (const fei::ParameterSet &paramset)
 
int parameters (int numParams, const char *const *paramStrings)
 
fei::SharedPtr< fei::MatrixGraphgetMatrixGraph () const
 
void setMatrixGraph (fei::SharedPtr< fei::MatrixGraph > matrixGraph)
 
int getGlobalNumRows () const
 
int getLocalNumRows () const
 
int getRowLength (int row, int &length) const
 
int putScalar (double scalar)
 
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 ()
 
void markState ()
 
bool changedSinceMark ()
 
const std::vector< int > & getRowNumbers () const
 
const std::vector< int > & getRowOffsets () const
 
const std::vector< int > & getColumnIndices () const
 
const std::vector< double > & getCoefs () const
 
- Public Member Functions inherited from fei::Matrix
virtual ~Matrix ()
 
virtual double * getBeginPointer ()
 
virtual int getOffset (int row, int col)
 

Static Public Member Functions

static fei::SharedPtr
< fei::Matrix
create_Matrix_Local (fei::SharedPtr< fei::MatrixGraph > matrixGraph, bool blockEntry)
 

Private Member Functions

int getRowIndex (int rowNumber) const
 
int giveToMatrix (int numRows, const int *rows, int numCols, const int *cols, const double *const *values, bool sumInto, int format)
 

Private Attributes

fei::SharedPtr< fei::MatrixGraphmatrixGraph_
 
fei::SharedPtr
< fei::SparseRowGraph
sparseRowGraph_
 
std::vector< double > coefs_
 
bool stateChanged_
 
std::vector< double > work_data1D_
 
std::vector< const double * > work_data2D_
 

Detailed Description

Definition at line 20 of file fei_Matrix_Local.hpp.

Constructor & Destructor Documentation

fei::Matrix_Local::Matrix_Local ( fei::SharedPtr< fei::MatrixGraph matrixGraph,
fei::SharedPtr< fei::SparseRowGraph sparseRowGraph 
)

Definition at line 17 of file fei_Matrix_Local.cpp.

fei::Matrix_Local::~Matrix_Local ( )
virtual

Definition at line 28 of file fei_Matrix_Local.cpp.

Member Function Documentation

fei::SharedPtr< fei::Matrix > fei::Matrix_Local::create_Matrix_Local ( fei::SharedPtr< fei::MatrixGraph matrixGraph,
bool  blockEntry 
)
static

Definition at line 34 of file fei_Matrix_Local.cpp.

References fei::MatrixGraph::createGraph().

const char * fei::Matrix_Local::typeName ( )
virtual

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

Implements fei::Matrix.

Definition at line 45 of file fei_Matrix_Local.cpp.

int fei::Matrix_Local::parameters ( const fei::ParameterSet paramset)
virtual

Method for supplying parameters

Implements fei::Matrix.

Definition at line 50 of file fei_Matrix_Local.cpp.

int fei::Matrix_Local::parameters ( int  numParams,
const char *const *  paramStrings 
)

Method for supplying parameters

Definition at line 57 of file fei_Matrix_Local.cpp.

fei::SharedPtr< fei::MatrixGraph > fei::Matrix_Local::getMatrixGraph ( ) const
virtual

Obtain the fei::MatrixGraph associated with this matrix

Implements fei::Matrix.

Definition at line 63 of file fei_Matrix_Local.cpp.

References matrixGraph_.

void fei::Matrix_Local::setMatrixGraph ( fei::SharedPtr< fei::MatrixGraph matrixGraph)
virtual

Set the fei::MatrixGraph associated with this matrix

Implements fei::Matrix.

Definition at line 67 of file fei_Matrix_Local.cpp.

References matrixGraph_.

int fei::Matrix_Local::getGlobalNumRows ( ) const
virtual

Get the global number of rows in the matrix.

Implements fei::Matrix.

Definition at line 71 of file fei_Matrix_Local.cpp.

References fei::SparseRowGraph::rowNumbers, and sparseRowGraph_.

Referenced by getLocalNumRows().

int fei::Matrix_Local::getLocalNumRows ( ) const
virtual

Get the local number of rows in the matrix.

Implements fei::Matrix.

Definition at line 75 of file fei_Matrix_Local.cpp.

References getGlobalNumRows().

Referenced by getRowIndex(), and writeToStream().

int fei::Matrix_Local::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 87 of file fei_Matrix_Local.cpp.

References getRowIndex(), fei::SparseRowGraph::rowOffsets, and sparseRowGraph_.

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

Set a specified scalar throughout the matrix.

Implements fei::Matrix.

Definition at line 98 of file fei_Matrix_Local.cpp.

References coefs_, and stateChanged_.

int fei::Matrix_Local::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
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.

Implements fei::Matrix.

Definition at line 106 of file fei_Matrix_Local.cpp.

References coefs_, getRowIndex(), fei::SparseRowGraph::packedColumnIndices, fei::SparseRowGraph::rowOffsets, and sparseRowGraph_.

int fei::Matrix_Local::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 177 of file fei_Matrix_Local.cpp.

References giveToMatrix().

int fei::Matrix_Local::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 187 of file fei_Matrix_Local.cpp.

References giveToMatrix().

int fei::Matrix_Local::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 197 of file fei_Matrix_Local.cpp.

References fei::MatrixGraph::getColSpace(), fei::VectorSpace::getFieldSize(), fei::VectorSpace::getGlobalIndex(), fei::MatrixGraph::getRowSpace(), giveToMatrix(), and matrixGraph_.

Referenced by sumInFieldData().

int fei::Matrix_Local::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 225 of file fei_Matrix_Local.cpp.

References fei::VectorSpace::getFieldSize(), fei::MatrixGraph::getRowSpace(), matrixGraph_, and sumInFieldData().

int fei::Matrix_Local::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 248 of file fei_Matrix_Local.cpp.

References fei::MatrixGraph::getConnectivityIndices(), fei::MatrixGraph::getConnectivityNumIndices(), giveToMatrix(), and matrixGraph_.

int fei::Matrix_Local::globalAssemble ( )
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.

Implements fei::Matrix.

Definition at line 263 of file fei_Matrix_Local.cpp.

int fei::Matrix_Local::multiply ( fei::Vector x,
fei::Vector y 
)
virtual

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

Implements fei::Matrix.

Definition at line 267 of file fei_Matrix_Local.cpp.

References FEI_COUT, and FEI_ENDL.

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

Definition at line 275 of file fei_Matrix_Local.cpp.

int fei::Matrix_Local::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 280 of file fei_Matrix_Local.cpp.

int fei::Matrix_Local::writeToFile ( const char *  filename,
bool  matrixMarketFormat = true 
)
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.

Implements fei::Matrix.

Definition at line 287 of file fei_Matrix_Local.cpp.

References FEI_OFSTREAM, FEI_OSTRINGSTREAM, fei::VectorSpace::getCommunicator(), fei::MatrixGraph::getRowSpace(), IOS_OUT, fei::localProc(), matrixGraph_, MPI_Comm, and writeToStream().

int fei::Matrix_Local::writeToStream ( FEI_OSTREAM ostrm,
bool  matrixMarketFormat = true 
)
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.

Implements fei::Matrix.

Definition at line 306 of file fei_Matrix_Local.cpp.

References coefs_, FEI_ENDL, fei::MatrixGraph::getColSpace(), fei::VectorSpace::getEqnNumbers(), getLocalNumRows(), fei::MatrixGraph::getRowSpace(), IOS_FLOATFIELD, IOS_SCIENTIFIC, matrixGraph_, fei::SparseRowGraph::packedColumnIndices, fei::SparseRowGraph::rowNumbers, fei::SparseRowGraph::rowOffsets, and sparseRowGraph_.

Referenced by writeToFile().

bool fei::Matrix_Local::usingBlockEntryStorage ( )
virtual

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

Implements fei::Matrix.

Definition at line 354 of file fei_Matrix_Local.cpp.

void fei::Matrix_Local::markState ( )
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.

Implements fei::Matrix.

Definition at line 358 of file fei_Matrix_Local.cpp.

References stateChanged_.

bool fei::Matrix_Local::changedSinceMark ( )
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.

Implements fei::Matrix.

Definition at line 364 of file fei_Matrix_Local.cpp.

References stateChanged_.

const std::vector< int > & fei::Matrix_Local::getRowNumbers ( ) const

Definition at line 368 of file fei_Matrix_Local.cpp.

References fei::SparseRowGraph::rowNumbers, and sparseRowGraph_.

const std::vector< int > & fei::Matrix_Local::getRowOffsets ( ) const

Definition at line 372 of file fei_Matrix_Local.cpp.

References fei::SparseRowGraph::rowOffsets, and sparseRowGraph_.

const std::vector< int > & fei::Matrix_Local::getColumnIndices ( ) const
const std::vector< double > & fei::Matrix_Local::getCoefs ( ) const

Definition at line 380 of file fei_Matrix_Local.cpp.

References coefs_.

int fei::Matrix_Local::getRowIndex ( int  rowNumber) const
private
int fei::Matrix_Local::giveToMatrix ( int  numRows,
const int *  rows,
int  numCols,
const int *  cols,
const double *const *  values,
bool  sumInto,
int  format 
)
private

Member Data Documentation

fei::SharedPtr<fei::MatrixGraph> fei::Matrix_Local::matrixGraph_
private
fei::SharedPtr<fei::SparseRowGraph> fei::Matrix_Local::sparseRowGraph_
private
std::vector<double> fei::Matrix_Local::coefs_
private

Definition at line 256 of file fei_Matrix_Local.hpp.

Referenced by copyOutRow(), getCoefs(), giveToMatrix(), putScalar(), and writeToStream().

bool fei::Matrix_Local::stateChanged_
private

Definition at line 257 of file fei_Matrix_Local.hpp.

Referenced by changedSinceMark(), giveToMatrix(), markState(), and putScalar().

std::vector<double> fei::Matrix_Local::work_data1D_
private

Definition at line 258 of file fei_Matrix_Local.hpp.

Referenced by giveToMatrix().

std::vector<const double*> fei::Matrix_Local::work_data2D_
private

Definition at line 259 of file fei_Matrix_Local.hpp.

Referenced by giveToMatrix().


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