FEI Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
fei::BlockLinearProblemManager Class Referenceabstract

#include <fei_BlockLinearProblemManager.hpp>

virtual ~BlockLinearProblemManager ()
 
virtual void setRowDistribution (const std::vector< int > &ownedIDs, const std::vector< int > &dofPerOwnedID, const std::vector< int > &ghostIDs, const std::vector< int > &dofPerGhostID)=0
 
virtual void setMatrixGraph (fei::SharedPtr< fei::SparseRowGraph > matrixGraph)=0
 
virtual void setMatrixValues (double scalar)=0
 
virtual int getNumOwnedIDs ()=0
 
virtual int getRowPointLength (int ownedID)=0
 
virtual int getRowBlockLength (int ownedID)=0
 
virtual int copyOutMatrixRow (int ownedID, int dofOffset, int numColIDs, int numCoefs, int *colIDs, int *dofPerColID, double *coefs)
 
virtual int insertMatrixValues (int rowID, int numRowDof, int colID, int numColDof, const double *const *values, bool sum_into)=0
 
virtual int insertMatrixValues (int rowID, int rowDofOffset, int colID, int colDofOffset, double value, bool sum_into)=0
 
virtual void setVectorValues (double scalar, bool soln_vector)=0
 
virtual int insertVectorValues (int ID, int numDof, const double *values, bool sum_into, bool soln_vector, int vectorIndex=0)=0
 
virtual int copyOutVectorValues (int ID, int numDof, double *values, bool soln_vector, int vectorIndex=0)=0
 
virtual double * getLocalVectorValuesPtr (bool soln_vector, int vectorIndex=0)=0
 
virtual int globalAssemble ()=0
 
virtual int solve (const fei::ParameterSet &parameters)=0
 

Detailed Description

Linear system assembly and solution manager.

This interface assumes that the parallel distribution of the linear system across processors is 1-D in the row-dimension, meaning that each processor owns a collection of complete rows of the matrix, and corresponding entries in the rhs and soln vectors.

Definition at line 27 of file fei_BlockLinearProblemManager.hpp.

Constructor & Destructor Documentation

virtual fei::BlockLinearProblemManager::~BlockLinearProblemManager ( )
inlinevirtual

Destructor.

Definition at line 32 of file fei_BlockLinearProblemManager.hpp.

Member Function Documentation

virtual void fei::BlockLinearProblemManager::setRowDistribution ( const std::vector< int > &  ownedIDs,
const std::vector< int > &  dofPerOwnedID,
const std::vector< int > &  ghostIDs,
const std::vector< int > &  dofPerGhostID 
)
pure virtual

Set the linear-system's global row distribution.

Parameters
ownedIDsList of IDs to be owned by local processor.
dofPerOwnedIDList giving num-degrees-of-freedom for each ID in the 'ownedGlobalIDs' list.
ghostIDsList of IDs that the local processor will need to know about, but which are not in the 'ownedIDs' list. ghostIDs correspond to matrix columns for which the corresponding row is owned by a remote processor.
dofPerGhostIDList giving num-degrees-of-freedom for each ID in the 'ownedGlobalIDs' list.
virtual void fei::BlockLinearProblemManager::setMatrixGraph ( fei::SharedPtr< fei::SparseRowGraph matrixGraph)
pure virtual

Set the matrix-graph structure. This is the nonzero structure for locally-owned matrix rows.

virtual void fei::BlockLinearProblemManager::setMatrixValues ( double  scalar)
pure virtual

Set a specified scalar value throughout the matrix.

virtual int fei::BlockLinearProblemManager::getNumOwnedIDs ( )
pure virtual

Query the number of locally-owned matrix block-rows.

virtual int fei::BlockLinearProblemManager::getRowPointLength ( int  ownedID)
pure virtual

Given a locally-owned matrix row, query the point-entry-length (number of scalar nonzeros) of that row. Note that although the block-row 'ownedID' corresponds to multiple point-rows, each of those point-rows has the same length.

virtual int fei::BlockLinearProblemManager::getRowBlockLength ( int  ownedID)
pure virtual

Given a locally-owned matrix row, query the block-entry-length (number of block columns) of that row.

virtual int fei::BlockLinearProblemManager::copyOutMatrixRow ( int  ownedID,
int  dofOffset,
int  numColIDs,
int  numCoefs,
int *  colIDs,
int *  dofPerColID,
double *  coefs 
)
virtual

Given a locally-owned point-entry row, pass out a copy of the contents of that row.

Parameters
ownedIDGlobal row id.
dofOffsetOffset into block-row 'ownedID'.
numColIDsLength of caller-allocated colIDs and dofPerColID arrays.
numCoefsLength of caller-allocated coefs array.
colIDsOutput array of block-column-indices.
dofPerColIDOutput array which will hold the number of scalar dof per column-id. i.e., the number of scalar coefficients that correspond to colIDs[i] is given by dofPerID[i].
coefsOutput array of matrix coefficients.
Returns
error-code 0 if successful. Non-zero return-value may indicate that the specified row is not locally owned, or that the lengths of the user-allocated arrays aren't consistent.
virtual int fei::BlockLinearProblemManager::insertMatrixValues ( int  rowID,
int  numRowDof,
int  colID,
int  numColDof,
const double *const *  values,
bool  sum_into 
)
pure virtual

Put a C-style table (array of pointers) of coefficient data into the matrix. This is a rectangular array of coefficients for a single block-entry specified by rowID and colID. If the sum_into argument is true, values should be added to any that already exist at the specified locations. Otherwise (if sum_into is false) incoming values should overwrite already-existing values.

virtual int fei::BlockLinearProblemManager::insertMatrixValues ( int  rowID,
int  rowDofOffset,
int  colID,
int  colDofOffset,
double  value,
bool  sum_into 
)
pure virtual

Put a single coefficient into the matrix. If the sum_into argument is true, value should be added to any that already exists at the specified location. Otherwise (if sum_into is false) incoming value should overwrite the already-existing value.

virtual void fei::BlockLinearProblemManager::setVectorValues ( double  scalar,
bool  soln_vector 
)
pure virtual

Set a specified scalar value throughout the vector.

Parameters
scalarValue to be used.
soln_vectorIf true, scalar should be set in solution vector, otherwise set rhs vector.
virtual int fei::BlockLinearProblemManager::insertVectorValues ( int  ID,
int  numDof,
const double *  values,
bool  sum_into,
bool  soln_vector,
int  vectorIndex = 0 
)
pure virtual

Put coefficient data into a vector at the specified global indices. If any specified indices are out of range (negative or too large), a positive warning code will be returned and the corresponding positions in the values array will not be referenced.

Parameters
IDGlobal id of block-entry position at which incoming values are to be stored.
numDofNumber of scalar coefficients that correspond to ID.
valuesInput list of coefficients.
sum_intoIf true, incoming values should be added to values that may already be in the specified locations. If sum_into is false, then incoming values should overwrite existing values.
soln_vectorIf true, incoming values should be placed in the solution vector. Otherwise, they should be placed in the rhs vector.
vectorIndexIf the linear system has multiple rhs/soln vectors, then this parameter specifies which vector the incoming values should be put into.
virtual int fei::BlockLinearProblemManager::copyOutVectorValues ( int  ID,
int  numDof,
double *  values,
bool  soln_vector,
int  vectorIndex = 0 
)
pure virtual

Copy values for the specified vector id into the caller's 'values' array.

virtual double* fei::BlockLinearProblemManager::getLocalVectorValuesPtr ( bool  soln_vector,
int  vectorIndex = 0 
)
pure virtual

Dangerous, high-performance vector access. Return a pointer to local vector values. Implementations that can't support this may return NULL, in which case the caller will revert to using the copyOutVectorValues method.

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

Perform any necessary internal communications/synchronizations or other operations appropriate at the end of data input. For some implementations this may be a no-op. (Trilinos/Epetra implementations would call 'FillComplete' on the matrix at this point.)

virtual int fei::BlockLinearProblemManager::solve ( const fei::ParameterSet parameters)
pure virtual

Solve the linear-system. (Most implementations will first internally perform the globalAssemble() operation if the caller hasn't called it explicitly.)

Returns
Error/status code. At the implementation's discretion, the return-code should be a value indicating the success/error status of the solve. Generally a return-value of 0 will indicate a successful solve.

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