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

#include <fei_LinearProblemManager.hpp>

Inherited by LinProbMgr_EpetraBasic.

virtual ~LinearProblemManager ()
 
virtual void setRowDistribution (const std::vector< int > &ownedGlobalRows)=0
 
virtual void setMatrixGraph (fei::SharedPtr< fei::SparseRowGraph > matrixGraph)=0
 
virtual void setMatrixValues (double scalar)=0
 
virtual int getLocalNumRows ()=0
 
virtual int getRowLength (int row)=0
 
virtual int copyOutMatrixRow (int row, int len, double *coefs, int *indices)=0
 
virtual int insertMatrixValues (int numRows, const int *rows, int numCols, const int *cols, const double *const *values, bool sum_into)=0
 
virtual void setVectorValues (double scalar, bool soln_vector)=0
 
virtual int insertVectorValues (int numValues, const int *globalIndices, const double *values, bool sum_into, bool soln_vector, int vectorIndex=0)=0
 
virtual int copyOutVectorValues (int numValues, const int *globalIndices, double *values, bool soln_vector, int vectorIndex=0)=0
 
virtual double * getLocalVectorValuesPtr (bool soln_vector, int vectorIndex=0)=0
 
virtual int globalAssemble ()=0
 

Detailed Description

Linear system assembly and solution manager.

Definition at line 24 of file fei_LinearProblemManager.hpp.

Constructor & Destructor Documentation

virtual fei::LinearProblemManager::~LinearProblemManager ( )
inlinevirtual

Destructor.

Definition at line 29 of file fei_LinearProblemManager.hpp.

Member Function Documentation

virtual void fei::LinearProblemManager::setRowDistribution ( const std::vector< int > &  ownedGlobalRows)
pure virtual

Set the linear-system's global row distribution.

Parameters
ownedGlobalRowsList of row-numbers to be owned by local processor.
virtual void fei::LinearProblemManager::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::LinearProblemManager::setMatrixValues ( double  scalar)
pure virtual

Set a specified scalar value throughout the matrix.

virtual int fei::LinearProblemManager::getLocalNumRows ( )
pure virtual

Query the number of local rows. This is expected to be the number of point-entry rows on the local processor.

virtual int fei::LinearProblemManager::getRowLength ( int  row)
pure virtual

Given a locally-owned global row number, query the length (number of nonzeros) of that row.

virtual int fei::LinearProblemManager::copyOutMatrixRow ( int  row,
int  len,
double *  coefs,
int *  indices 
)
pure virtual

Given a locally-owned global row number, pass out a copy of the contents of that row.

Parameters
rowGlobal row number
lenLength of user-allocated 'coefs' and 'indices' arrays. if 'len' != 'getRowLength(row)', then the number of coefs/indices returned will be max(len,getRowLength(row)).
coefsOutput array of matrix coefficients.
indicesOutput array of column-indices.
Returns
error-code 0 if successful. Non-zero return-value may indicate that the specified row is not locally owned.
virtual int fei::LinearProblemManager::insertMatrixValues ( int  numRows,
const int *  rows,
int  numCols,
const int *  cols,
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 rows/columns defined by the 'rows' and 'cols' lists. 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 void fei::LinearProblemManager::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::LinearProblemManager::insertVectorValues ( int  numValues,
const int *  globalIndices,
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) the corresponding positions in the values array will not be referenced, and a positive warning code will be returned.

Parameters
numValuesLength of caller-allocated 'globalIndices' and 'values' arrays.
globalIndicesList of global-indices specifying the locations in the vector for incoming values to be placed.
valuesList of incoming values.
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::LinearProblemManager::copyOutVectorValues ( int  numValues,
const int *  globalIndices,
double *  values,
bool  soln_vector,
int  vectorIndex = 0 
)
pure virtual

Copy values for the specified vector indices into the caller-allocated 'values' array.

virtual double* fei::LinearProblemManager::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::LinearProblemManager::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.)


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