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 | Private Attributes | List of all members
LinProbMgr_EpetraBasic Class Reference

#include <fei_LinProbMgr_EpetraBasic.hpp>

Inheritance diagram for LinProbMgr_EpetraBasic:
Inheritance graph
[legend]

Public Member Functions

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

Private Attributes

MPI_Comm comm_
 
std::vector< int > ownedRows_
 
fei::SharedPtr< Epetra_Commepetra_comm_
 
fei::SharedPtr< Epetra_Mapepetra_rowmap_
 
fei::SharedPtr
< fei::SparseRowGraph
fei_srgraph_
 
fei::SharedPtr< Epetra_CrsGraphcrsgraph_
 
fei::SharedPtr< Epetra_CrsMatrixA_
 
int numVectors_
 
fei::SharedPtr
< Epetra_MultiVector
x_
 
fei::SharedPtr
< Epetra_MultiVector
b_
 

Detailed Description

Definition at line 54 of file fei_LinProbMgr_EpetraBasic.hpp.

Constructor & Destructor Documentation

LinProbMgr_EpetraBasic::LinProbMgr_EpetraBasic ( MPI_Comm  comm)
virtual LinProbMgr_EpetraBasic::~LinProbMgr_EpetraBasic ( )
virtual

Member Function Documentation

void LinProbMgr_EpetraBasic::setRowDistribution ( const std::vector< int > &  ownedGlobalRows)
virtual

Set the linear-system's global row distribution.

Parameters
ownedGlobalRowsList of row-numbers to be owned by local processor.

Implements fei::LinearProblemManager.

void LinProbMgr_EpetraBasic::setMatrixGraph ( fei::SharedPtr< fei::SparseRowGraph matrixGraph)
virtual

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

Implements fei::LinearProblemManager.

void LinProbMgr_EpetraBasic::setMatrixValues ( double  scalar)
virtual

Set a specified scalar value throughout the matrix.

Implements fei::LinearProblemManager.

void LinProbMgr_EpetraBasic::setVectorValues ( double  scalar,
bool  soln_vector 
)
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.

Implements fei::LinearProblemManager.

int LinProbMgr_EpetraBasic::getLocalNumRows ( )
virtual

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

Implements fei::LinearProblemManager.

int LinProbMgr_EpetraBasic::getRowLength ( int  row)
virtual

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

Implements fei::LinearProblemManager.

int LinProbMgr_EpetraBasic::copyOutMatrixRow ( int  row,
int  len,
double *  coefs,
int *  indices 
)
virtual

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

Parameters
rowGlobal row number
lenLength of the user-allocated arrays coefs and indices.
coefsUser-allocated array which will hold matrix coefficients on output.
indicesUser-allocated array which will hold column-indices on output.
Returns
error-code 0 if successful. Non-zero return-value may indicate that the specified row is not locally owned.

Implements fei::LinearProblemManager.

int LinProbMgr_EpetraBasic::insertMatrixValues ( int  numRows,
const int *  rows,
int  numCols,
const int *  cols,
const double *const *  values,
bool  sum_into 
)
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.

Implements fei::LinearProblemManager.

int LinProbMgr_EpetraBasic::insertVectorValues ( int  numValues,
const int *  globalIndices,
const double *  values,
bool  sum_into,
bool  soln_vector,
int  vectorIndex = 0 
)
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
numValuesNumber of coefficient values being input.
globalIndicesList of global-indices specifying the locations in the vector for incoming values to be placed.
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.

Implements fei::LinearProblemManager.

int LinProbMgr_EpetraBasic::copyOutVectorValues ( int  numValues,
const int *  globalIndices,
double *  values,
bool  soln_vector,
int  vectorIndex = 0 
)
virtual

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

Implements fei::LinearProblemManager.

double* LinProbMgr_EpetraBasic::getLocalVectorValuesPtr ( bool  soln_vector,
int  vectorIndex = 0 
)
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.

Implements fei::LinearProblemManager.

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

Implements fei::LinearProblemManager.

fei::SharedPtr<Epetra_CrsMatrix> LinProbMgr_EpetraBasic::get_A_matrix ( )

Return the Epetra matrix.

fei::SharedPtr<Epetra_MultiVector> LinProbMgr_EpetraBasic::get_rhs_vector ( )

Return the rhs Epetra vector.

fei::SharedPtr<Epetra_MultiVector> LinProbMgr_EpetraBasic::get_solution_vector ( )

Return the solution Epetra vector.

Member Data Documentation

MPI_Comm LinProbMgr_EpetraBasic::comm_
private

Definition at line 186 of file fei_LinProbMgr_EpetraBasic.hpp.

std::vector<int> LinProbMgr_EpetraBasic::ownedRows_
private

Definition at line 187 of file fei_LinProbMgr_EpetraBasic.hpp.

fei::SharedPtr<Epetra_Comm> LinProbMgr_EpetraBasic::epetra_comm_
private

Definition at line 188 of file fei_LinProbMgr_EpetraBasic.hpp.

fei::SharedPtr<Epetra_Map> LinProbMgr_EpetraBasic::epetra_rowmap_
private

Definition at line 189 of file fei_LinProbMgr_EpetraBasic.hpp.

fei::SharedPtr<fei::SparseRowGraph> LinProbMgr_EpetraBasic::fei_srgraph_
private

Definition at line 190 of file fei_LinProbMgr_EpetraBasic.hpp.

fei::SharedPtr<Epetra_CrsGraph> LinProbMgr_EpetraBasic::crsgraph_
private

Definition at line 191 of file fei_LinProbMgr_EpetraBasic.hpp.

fei::SharedPtr<Epetra_CrsMatrix> LinProbMgr_EpetraBasic::A_
private

Definition at line 192 of file fei_LinProbMgr_EpetraBasic.hpp.

int LinProbMgr_EpetraBasic::numVectors_
private

Definition at line 193 of file fei_LinProbMgr_EpetraBasic.hpp.

fei::SharedPtr<Epetra_MultiVector> LinProbMgr_EpetraBasic::x_
private

Definition at line 194 of file fei_LinProbMgr_EpetraBasic.hpp.

fei::SharedPtr<Epetra_MultiVector> LinProbMgr_EpetraBasic::b_
private

Definition at line 195 of file fei_LinProbMgr_EpetraBasic.hpp.


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