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

#include <fei_LinearSystemCore.hpp>

Inheritance diagram for LinearSystemCore:
Inheritance graph
[legend]

Public Member Functions

 LinearSystemCore ()
 
virtual ~LinearSystemCore ()
 
virtual LinearSystemCoreclone ()=0
 
virtual int parameters (int numParams, const char *const *params)=0
 
virtual int setLookup (Lookup &lookup)=0
 
virtual int getProperty (const char *, double &)
 
virtual int setGlobalOffsets (int len, int *nodeOffsets, int *eqnOffsets, int *blkEqnOffsets)=0
 
virtual int setConnectivities (GlobalID elemBlock, int numElements, int numNodesPerElem, const GlobalID *elemIDs, const int *const *connNodes)=0
 
virtual int setStiffnessMatrices (GlobalID elemBlock, int numElems, const GlobalID *elemIDs, const double *const *const *stiff, int numEqnsPerElem, const int *const *eqnIndices)=0
 
virtual int setLoadVectors (GlobalID elemBlock, int numElems, const GlobalID *elemIDs, const double *const *load, int numEqnsPerElem, const int *const *eqnIndices)=0
 
virtual int setMatrixStructure (int **ptColIndices, int *ptRrowLengths, int **blkColIndices, int *blkRowLengths, int *ptRowsPerBlkRow)=0
 
virtual int setMultCREqns (int multCRSetID, int numCRs, int numNodesPerCR, int **nodeNumbers, int **eqnNumbers, int *fieldIDs, int *multiplierEqnNumbers)=0
 
virtual int setPenCREqns (int penCRSetID, int numCRs, int numNodesPerCR, int **nodeNumbers, int **eqnNumbers, int *fieldIDs)=0
 
virtual int sumIntoSystemMatrix (int numPtRows, const int *ptRows, int numPtCols, const int *ptCols, int numBlkRows, const int *blkRows, int numBlkCols, const int *blkCols, const double *const *values)=0
 
virtual int sumIntoSystemMatrix (int numPtRows, const int *ptRows, int numPtCols, const int *ptCols, const double *const *values)=0
 
virtual int putIntoSystemMatrix (int numPtRows, const int *ptRows, int numPtCols, const int *ptCols, const double *const *values)=0
 
virtual int getMatrixRowLength (int row, int &length)=0
 
virtual int getMatrixRow (int row, double *coefs, int *indices, int len, int &rowLength)=0
 
virtual int sumIntoRHSVector (int num, const double *values, const int *indices)=0
 
virtual int putIntoRHSVector (int num, const double *values, const int *indices)=0
 
virtual int getFromRHSVector (int num, double *values, const int *indices)=0
 
virtual int matrixLoadComplete ()=0
 
virtual int putNodalFieldData (int fieldID, int fieldSize, int *nodeNumbers, int numNodes, const double *data)=0
 
virtual int resetMatrixAndVector (double s)=0
 
virtual int resetMatrix (double s)=0
 
virtual int resetRHSVector (double s)=0
 
virtual int enforceEssentialBC (int *globalEqn, double *alpha, double *gamma, int len)=0
 
virtual int enforceRemoteEssBCs (int numEqns, int *globalEqns, int **colIndices, int *colIndLen, double **coefs)=0
 
virtual int getMatrixPtr (Data &data)=0
 
virtual int copyInMatrix (double scalar, const Data &data)=0
 
virtual int copyOutMatrix (double scalar, Data &data)=0
 
virtual int sumInMatrix (double scalar, const Data &data)=0
 
virtual int getRHSVectorPtr (Data &data)=0
 
virtual int copyInRHSVector (double scalar, const Data &data)=0
 
virtual int copyOutRHSVector (double scalar, Data &data)=0
 
virtual int sumInRHSVector (double scalar, const Data &data)=0
 
virtual int destroyMatrixData (Data &data)=0
 
virtual int destroyVectorData (Data &data)=0
 
virtual int setNumRHSVectors (int numRHSs, const int *rhsIDs)=0
 
virtual int setRHSID (int rhsID)=0
 
virtual int putInitialGuess (const int *eqnNumbers, const double *values, int len)=0
 
virtual int getSolution (double *answers, int len)=0
 
virtual int getSolnEntry (int eqnNumber, double &answer)=0
 
virtual int formResidual (double *values, int len)=0
 
virtual int launchSolver (int &solveStatus, int &iterations)=0
 
virtual int writeSystem (const char *name)=0
 

Detailed Description

This is the original internal FEI interface to solver-libraries – the destination for the data being assembled into a linear system by the FEI implementation.

When creating a specific FEI implementation, i.e., a version that supports a specific underlying linear solver library, the main task that must be performed is the implementation of this interface, LinearSystemCore (Note: the LinearSystemCore interface is being replaced by the ESI_Broker interface.).

To date (as of August 2001), implementations of this interface exist for coupling the following solver libraries to the FEI implementation:

An implementation of LinearSystemCore holds and manipulates all solver-library-specific stuff, such as matrices/vectors, solvers/preconditioners, etc. An instance of this class is owned and used by the class that implements the public FEI spec. i.e., when element contributions, etc., are received from the finite-element application, the data is ultimately passed to this class for assembly into the sparse matrix and associated vectors. This class will also be asked to launch any underlying solver, and finally to return the solution.

Terminology notes:

NOTES:

  1. The LinearSystemCore object is given a 'Lookup' interface, from which it can obtain information about the structure of the finite-element problem. e.g., it can look up the number of fields, the fields' sizes, the number of element-blocks, etc., etc. For details on the Lookup interface, see the file Lookup.h.
  2. Except for making calls on the Lookup object, LinearSystemCore is reactive. All of the functions below are member functions of LinearSystemCore, which are called by FEI code. 'set' and 'put' functions are for information to be provided TO LinearSystemCore by the FEI implementation layer, while 'get' functions are for information to be requested FROM LinearSystemCore by the FEI implementation layer.
  3. Node-numbers that are given to LinearSystemCore (through setConnectivities, etc.,) are 0-based and contiguous. Each processor owns a contiguous block of node-numbers. However, non-local node-numbers may of course appear in connectivity lists. These node-numbers are not the same as the 'nodeIDs' that appear in the FEI public interface, and throughout the FEI implementation code. nodeIDs are provided by the application, and may not be assumed to be 0-based or contiguous. They are arbitrary node identifiers which need only be globally unique. nodeIDs never appear in LinearSystemCore functions.
  4. Not all LinearSystemCore functions are necessary for assembling a linear system. For instance, many implementations will ignore the data that is passed in by 'setConnectivities', 'setStiffnessMatrices', etc.
  5. Some data is passed into LinearSystemCore redundantly. i.e., matrix coefficient data appears in the 'setStiffnessMatrices' function as well as the 'sumIntoSystemMatrix/sumIntoSystemBlkMatrix' function. The difference is, only local data is supplied through sumIntoSystemMatrix. Additionally, the FEI implementation layer performs the necessary communication to ensure that data from shared finite-element nodes is moved onto the 'owning' processor before being passed to LinearSystemCore::sumIntoSystemMatrix. The finite-element application doesn't have the concept of node ownership, instead processors can equally share nodes. The FEI implementation designates a single processor as the owner of these nodes, and moves any shared stiffness data onto the appropriate owning processor. The data supplied through setStiffnessMatrices is the un-modified stiffness data supplied by the application, including any shared-non-local portions.

    Similarly for setLoadVectors and sumIntoRHSVector. sumIntoRHSVector only provides local data, while setLoadVectors supplies the un-modified element- load vectors from the application, including any shared portions.

Definition at line 124 of file fei_LinearSystemCore.hpp.

Constructor & Destructor Documentation

LinearSystemCore::LinearSystemCore ( )
inline

Default constructor, typically overridden by the implementing class, which probably requires an MPI Communicator, and possibly other arguments.

Definition at line 129 of file fei_LinearSystemCore.hpp.

virtual LinearSystemCore::~LinearSystemCore ( )
inlinevirtual

Destructor, should call through to the implementation's destructor and destroy all allocated memory, internal objects, etc. Exceptions: objects created in reponse to calls to the functions 'copyOutMatrix' and 'copyOutRHSVector' are not destroyed here. The caller is assumed to have taken responsibility for those matrix/vector copies.

Definition at line 137 of file fei_LinearSystemCore.hpp.

Member Function Documentation

virtual LinearSystemCore* LinearSystemCore::clone ( )
pure virtual

For cloning a LinearSystemCore instance. Caller recieves a pointer to a new instantiation of the implementing object.

virtual int LinearSystemCore::parameters ( int  numParams,
const char *const *  params 
)
pure virtual

For setting argc/argv style parameters.

Parameters
numParamsNumber of strings in the params argument
paramsA list of strings which will usually contain space-separated key-value pairs. Example: "debugOutput /usr/users/me/work_dir"
virtual int LinearSystemCore::setLookup ( Lookup lookup)
pure virtual

Supply the LinearSystemCore implementation with an object (created and owned by the caller) that can be used to obtain various information about problem layout, shared finite-element nodes, etc. For details, see the documentation for the Lookup interface.

Parameters
lookupInput. Reference to an implementation of the Lookup interface
virtual int LinearSystemCore::getProperty ( const char *  ,
double &   
)
inlinevirtual

Query a named property (such as timing statistics, etc.) from the solver library.

Parameters
nameInput. Name of the property for which a value is being requested. value Output. Requested property's value.
Returns
error-code 0 if successful. -1 probably indicates that the named property is not recognized.

Definition at line 171 of file fei_LinearSystemCore.hpp.

virtual int LinearSystemCore::setGlobalOffsets ( int  len,
int *  nodeOffsets,
int *  eqnOffsets,
int *  blkEqnOffsets 
)
pure virtual

Supply LinearSystemCore with global offset information for the problem being assembled.

Parameters
lenLength of the following list arguments. This will be numProcs+1
nodeOffsetsThe FEI implementation assigns a global 0-based numbering to the finite-element nodes in the problem. Each processor is given ownership of a contiguous subset of these node-numbers. nodeOffsets[i] gives the first local node-number for processor i. nodeOffsets[len-1] gives the total number of nodes.
eqnOffsetseqnOffsets[i] gives the first local equation number for processor i, eqnOffsets[len-1] gives the global number of equations.
blkEqnOffsetsContains the same kind of information as eqnOffsets, but for 'block-equations'. A block-equation contains all of the point-equations present at a finite-element node. Special case: if this problem contains Lagrange Multiplier constraints, they will be equations that don't correspond to any node, and there will only be one of these equations mapped to a block-equation.
virtual int LinearSystemCore::setConnectivities ( GlobalID  elemBlock,
int  numElements,
int  numNodesPerElem,
const GlobalID *  elemIDs,
const int *const *  connNodes 
)
pure virtual

For passing element-connectivity arrays.

Parameters
elemBlockIdentifier for the element-block that these elements belong to.
numElementsLength of the elemIDs list.
numNodesPerElemLength of each row in the connNodes table.
elemIDsIdentifiers for each element for which connectivities are being supplied.
connNodesTable, with one row for each element. Each row is a list of the nodes that are connected to that element.
virtual int LinearSystemCore::setStiffnessMatrices ( GlobalID  elemBlock,
int  numElems,
const GlobalID *  elemIDs,
const double *const *const *  stiff,
int  numEqnsPerElem,
const int *const *  eqnIndices 
)
pure virtual

For passing element-stiffness arrays.

Parameters
elemBlockIdentifier for the element-block that these elements belong to.
numElemsLength of the elemIDs list.
elemIDsIdentifiers for each element for which a stiffness array is being supplied.
stiffList of 'numElems' tables, each table is of size 'numEqnsPerElem' X 'numEqnsPerElem'.
numEqnsPerElem
eqnIndicesTable, with 'numElems' rows, each row being a list of 'numEqnsPerElem' scatter indices (0-based global matrix row/column indices).
virtual int LinearSystemCore::setLoadVectors ( GlobalID  elemBlock,
int  numElems,
const GlobalID *  elemIDs,
const double *const *  load,
int  numEqnsPerElem,
const int *const *  eqnIndices 
)
pure virtual

For passing element-load vectors.

Parameters
elemBlockIdentifier for the element-block that these elements belong to.
numElemsLength of the elemIDs list.
elemIDsIdentifiers for each element for which a load vector is being supplied.
loadTable with 'numElems' rows, each row is of length 'numEqnsPerElem'.
numEqnsPerElem
eqnIndicesTable, with 'numElems' rows, each row being a list of 'numEqnsPerElem' scatter indices (0-based global equation numbers).
virtual int LinearSystemCore::setMatrixStructure ( int **  ptColIndices,
int *  ptRrowLengths,
int **  blkColIndices,
int *  blkRowLengths,
int *  ptRowsPerBlkRow 
)
pure virtual

Supply LinearSystemCore with information defining the structure of the sparse matrix to be assembled. Implementers of LinearSystemCore may safely assume that this function will not be called until after the function 'setGlobalOffsets' has been called. Using the information provided via setGlobalOffsets, the number-of-local-equations can be trivially calculated. After setMatrixStructure has been called, there should be enough information to instantiate internal linear-algebra entities, such as vectors, matrix, etc.

Parameters
ptColIndicesTable, with num-local-eqns rows, and the i-th row is of length ptRowLengths[i].
ptRowLengths
blkColIndicesTable, with num-local-blkEqns rows, and the i-th row is of length blkRowLengths[i].
blkRowLengths
ptRowsPerBlkRowThe i-th local block-equation corresponds to ptRowsPerBlkRow[i] point-equations.
virtual int LinearSystemCore::setMultCREqns ( int  multCRSetID,
int  numCRs,
int  numNodesPerCR,
int **  nodeNumbers,
int **  eqnNumbers,
int *  fieldIDs,
int *  multiplierEqnNumbers 
)
pure virtual

Specify which global equation numbers correspond to Lagrange Multiplier equations. This function won't be called if there are no Lagrange Multiplier constraints in the problem. If this function is called, it is guaranteed to be called after 'setGlobalOffsets' and before 'setMatrixStructure'. The primary purpose of this function is to give LinearSystemCore implementers the opportunity to deal with constraints in a special way, rather than assembling everything into one matrix. If the problem being assembled does have Lagrange constraints, then the FEI implementation will request an ESI_MatrixRowWriteAccess interface for "C_Matrix" from LinearSystemCore. If that is not available, then the FEI implementation will request "A_Matrix" and assemble everything into that.

Parameters
numCRsnumber of constraint relations
numNodesPerCRnumber of constrained node in each constraint relation
nodeNumbersTable of constrained nodes. 'numCRs' rows, with the i-th row being of length numNodesPerCR[i].
eqnNumbersTable, same dimensions as 'nodeNumbers'. These are the global 0-based matrix column indices of the constraint coefficients.
multiplierEqnNumbersEquation numbers that the Lagrange Multipliers correspond to.
virtual int LinearSystemCore::setPenCREqns ( int  penCRSetID,
int  numCRs,
int  numNodesPerCR,
int **  nodeNumbers,
int **  eqnNumbers,
int *  fieldIDs 
)
pure virtual

Specify which nodes and equation numbers correspond to penalty constraints. This function is included for completeness, but hasn't yet been proven to be useful or necessary, and will probably not be included in the successor to LinearSystemCore (ESI_LSManager).

virtual int LinearSystemCore::sumIntoSystemMatrix ( int  numPtRows,
const int *  ptRows,
int  numPtCols,
const int *  ptCols,
int  numBlkRows,
const int *  blkRows,
int  numBlkCols,
const int *  blkCols,
const double *const *  values 
)
pure virtual
Provides point-entry data, as well as block-entry data. This is the

primary assembly function, through which the FEI implementation provides the local equation contributions of all element contributions.

virtual int LinearSystemCore::sumIntoSystemMatrix ( int  numPtRows,
const int *  ptRows,
int  numPtCols,
const int *  ptCols,
const double *const *  values 
)
pure virtual
 Purely point-entry version for accumulating coefficient data into the 
matrix. This will be called when a matrix contribution fills only part of
a block-equation. e.g., when a penalty constraint is being applied to a

single solution field on a node that has several solution fields. (A block-equation contains all solution field equations at a node.)

virtual int LinearSystemCore::putIntoSystemMatrix ( int  numPtRows,
const int *  ptRows,
int  numPtCols,
const int *  ptCols,
const double *const *  values 
)
pure virtual

Point-entry matrix data as for 'sumIntoSystemMatrix', but in this case the data should be "put" into the matrix (i.e., overwrite any coefficients already present) rather than being "summed" into the matrix.

virtual int LinearSystemCore::getMatrixRowLength ( int  row,
int &  length 
)
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, e.g., row is not local.
virtual int LinearSystemCore::getMatrixRow ( int  row,
double *  coefs,
int *  indices,
int  len,
int &  rowLength 
)
pure virtual

Obtain 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
rowLengthOutput. Actual length of this row. Not referenced if row is not in the local portion of the matrix.
Returns
error-code non-zero if any error occurs, e.g., row is not local.
virtual int LinearSystemCore::sumIntoRHSVector ( int  num,
const double *  values,
const int *  indices 
)
pure virtual

For accumulating coefficients into the rhs vector

virtual int LinearSystemCore::putIntoRHSVector ( int  num,
const double *  values,
const int *  indices 
)
pure virtual

For putting coefficients into the rhs vector

virtual int LinearSystemCore::getFromRHSVector ( int  num,
double *  values,
const int *  indices 
)
pure virtual

For getting coefficients out of the rhs vector

virtual int LinearSystemCore::matrixLoadComplete ( )
pure virtual

The FEI implementation calls this function to signal the linsyscore object that data-loading is finished.

virtual int LinearSystemCore::putNodalFieldData ( int  fieldID,
int  fieldSize,
int *  nodeNumbers,
int  numNodes,
const double *  data 
)
pure virtual

Pass nodal data that probably doesn't mean anything to the FEI implementation, but may mean something to the linear solver. Examples: geometric coordinates, nullspace data, etc.

Parameters
fieldIDIdentifier for the field that describes this data. Lists of field identifiers and field sizes defined for the finite-element problem may be obtained from the Lookup interface that is supplied to the LinearSystemCore by the FEI implementation.
nodeNumbersList of nodes for which data is being supplied.
numNodes
dataList of length numNodes * (size of field 'fieldID')
virtual int LinearSystemCore::resetMatrixAndVector ( double  s)
pure virtual

For setting the scalar 's' (usually 0.0) throughout the matrix rhs vector.

virtual int LinearSystemCore::resetMatrix ( double  s)
pure virtual

For setting the scalar 's' (usually 0.0) throughout the matrix.

virtual int LinearSystemCore::resetRHSVector ( double  s)
pure virtual

For setting the scalar 's' (usually 0.0) throughout the rhs vector.

virtual int LinearSystemCore::enforceEssentialBC ( int *  globalEqn,
double *  alpha,
double *  gamma,
int  len 
)
pure virtual

The FEI implementation calls this function to inform LinearSystemCore of equations that need to have essential (Dirichlet) boundary conditions enforced on them. The intent is that the LinearSystemCore implementation will perform the column-modification b[i] = gamma[i]/alpha[i] if i == globalEqn[i], otherwise b[i] -= gamma[i]/alpha[i] * A(i,globalEqn[i]) if i != globalEqn[i]. After this operation is performed, all of row globalEqn[i] and column globalEqn[i] should be set to 0.0, except for the diagonal position. (Naturally the implementer is free to enforce the boundary condition another way if they wish.)

Parameters
globalEqnList, of length 'len', of global 0-based equation numbers.
alphaList, of length 'len', of coefficients. When the solution to the linear system is later requested, the solution value for globalEqn[i] should be gamma[i]/alpha[i].
gamma
len
virtual int LinearSystemCore::enforceRemoteEssBCs ( int  numEqns,
int *  globalEqns,
int **  colIndices,
int *  colIndLen,
double **  coefs 
)
pure virtual

The FEI implementation calls this function to inform LinearSystemCore that certain local equations need column-modifications made due to essential boundary-conditions being enforced on other processors. The column modification is roughly this: b(globalEqns[i]) -= A(globalEqns[i], colIndices[i][j]) * coefs[i][j], for i in [0..numEqns-1] and j in [0.. colIndLen[i]-1]. (Note that A(globalEqns[i], colIndices[i][j]) should be set = 0.0 after the appropriate value has been accumulated into b also.)

Parameters
numEqnsLength of 'globalEqns'
globalEqnsEquations that are local to this processor.
colIndicesTable, with 'numEqns' rows, and the i-th row is of length colIndLen[i]. The i-th row contains column indices in equation globalEqns[i]. These column indices correspond to equations (rows) that are owned by other processors, and those other processors are imposing essential boundary conditions on those equations.
colIndLenList of length 'numEqns'.
coefsThis table holds the gamma/alpha coeficients that are the value of the boundary conditions being enforced on each of the remote equations in the 'colIndices' table.
virtual int LinearSystemCore::getMatrixPtr ( Data data)
pure virtual

The FEI implementation calls this function to request a pointer to the internal 'A-matrix' data.

Parameters
dataSee Data class documentation.
virtual int LinearSystemCore::copyInMatrix ( double  scalar,
const Data data 
)
pure virtual

LinearSystemCore's internal 'A-matrix' should be replaced with a scaled copy of the incoming data.

Parameters
scalarcoefficient by which to scale the incoming data.
dataSee documentation for Data class.
virtual int LinearSystemCore::copyOutMatrix ( double  scalar,
Data data 
)
pure virtual

The FEI implementation calls this function to request a scaled copy of the internal 'A-matrix' data. The FEI implementation will then be responsible for deciding when this matrix data should be destroyed. The LinearSystemCore implementation should not keep a reference to the pointer that was handed out.

Parameters
scalar
dataSee documentation for Data class.
virtual int LinearSystemCore::sumInMatrix ( double  scalar,
const Data data 
)
pure virtual

A scaled copy of the incoming data should be added to the internal 'A-matrix'.

Parameters
scalar
dataSee documentation for Data class.
virtual int LinearSystemCore::getRHSVectorPtr ( Data data)
pure virtual

Same semantics as getMatrixPtr, but applied to rhs vector.

virtual int LinearSystemCore::copyInRHSVector ( double  scalar,
const Data data 
)
pure virtual

Same semantics as copyInMatrix, but applied to rhs vector.

virtual int LinearSystemCore::copyOutRHSVector ( double  scalar,
Data data 
)
pure virtual

Same semantics as copyOutMatrix, but applied to rhs vector.

virtual int LinearSystemCore::sumInRHSVector ( double  scalar,
const Data data 
)
pure virtual

Same semantics as sumInMatrix, but applied to rhs vector.

virtual int LinearSystemCore::destroyMatrixData ( Data data)
pure virtual

Utility function for destroying the matrix in a Data container. The caller (owner of 'data') can't destroy the matrix because they don't know what concrete type it is and can't get to its destructor. The contents of 'data' is a matrix previously passed out via 'copyOutMatrix'.

Parameters
dataSee documentation for Data class.
virtual int LinearSystemCore::destroyVectorData ( Data data)
pure virtual
  Utility function for destroying the vector in a Data container. The 

caller (owner of 'data') can't destroy the vector because they don't know what concrete type it is and can't get to its destructor. The contents of 'data' is a vector previously passed out via 'copyOutRHSVector'.

Parameters
dataSee documentation for Data class.
virtual int LinearSystemCore::setNumRHSVectors ( int  numRHSs,
const int *  rhsIDs 
)
pure virtual

Indicate the number of rhs-vectors being assembled/solved for. This function will be called by the FEI implementation at or near the beginning of the problem assembly. If numRHSs is greater than 1, then calls to 'getMemberInterface' requesting an interface to an rhs vector will use the 'objName' argument "b_Vector_n", where n is in [0 .. numRHSs-1]. If there is only one rhs-vector, then 'objName' will be simply "b_Vector".

Parameters
numRHSsLength of the rhsIDs list.
rhsIDsCaller-supplied integer identifiers for the rhs vectors. This argument will probably be removed, as it is obsolete (a carry-over from LinearSystemCore).
virtual int LinearSystemCore::setRHSID ( int  rhsID)
pure virtual

Set the 'current context' to the rhs-vector corresponding to rhsID. Subsequent data received via 'sumIntoRHSVector' should be directed into this rhs-vector. Any other function-calls having to do with the rhs-vector should also effect this rhs-vector.

Parameters
rhsID
virtual int LinearSystemCore::putInitialGuess ( const int *  eqnNumbers,
const double *  values,
int  len 
)
pure virtual

The FEI implementation will call this function to supply initial-guess data that should be used as the starting 'x-vector' if an iterative solution is to be performed.

Parameters
eqnNumbersGlobal 0-based equation numbers for which the initial guess should be set.
valuesThe initial guess data.
lenNumber of equations for which an initial guess is being supplied.
virtual int LinearSystemCore::getSolution ( double *  answers,
int  len 
)
pure virtual

The FEI implementation will call this function to request all local solution values.

Parameters
answersSolution coefficients.
lenThis should equal the number of local equations. If it is less, the LinearSystemCore implementation should simply pass out the first 'len' local solution values. If it is more, then just pass out numLocalEqns solution values.
virtual int LinearSystemCore::getSolnEntry ( int  eqnNumber,
double &  answer 
)
pure virtual

The FEI implementation will call this function to request a single solution value.

Parameters
eqnNumberGlobal 0-based equation number.
answer
virtual int LinearSystemCore::formResidual ( double *  values,
int  len 
)
pure virtual

This will be called to request that LinearSystemCore form the residual vector r = b - A*x, and pass the coefficients for r back out in the 'values' list.

Parameters
values
lenThis should equal num-local-eqns.
virtual int LinearSystemCore::launchSolver ( int &  solveStatus,
int &  iterations 
)
pure virtual

Function called to request the launching of the linear solver.

Parameters
solveStatusOutput, should indicate the status of the solve. A successful solve is usually indicated by a value of 0.
iterationsOutput, how many iterations were performed.
Returns
error-code, 0 if convergence tolerance was achieved within the specified maximum number of iterations. If error return is non-zero, the calling application will be expected to check solveStatus, and consult the solver-library's documentation to figure out exactly what happened.
virtual int LinearSystemCore::writeSystem ( const char *  name)
pure virtual

This function's intent is to provide a file-name to be used by LinearSystemCore in writing the linear system into disk files. Format is not specified. Implementers may choose to augment this name in the style of writing 3 files: A_name, x_name, b_name, or some other convention. This function is ill-defined, obsolete, and will probably never be called by the FEI implementation.


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