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 Member Functions | Private Attributes | List of all members
fei_trilinos::Aztec_LinSysCore Class Reference

#include <fei_Aztec_LinSysCore.hpp>

Inheritance diagram for fei_trilinos::Aztec_LinSysCore:
Inheritance graph
[legend]

Public Member Functions

 Aztec_LinSysCore (MPI_Comm comm)
 
virtual ~Aztec_LinSysCore ()
 
LinearSystemCoreclone ()
 
int parameters (int numParams, const char *const *params)
 
int setLookup (Lookup &lookup)
 
int setGlobalOffsets (int len, int *nodeOffsets, int *eqnOffsets, int *blkEqnOffsets)
 
int setConnectivities (GlobalID elemBlock, int numElements, int numNodesPerElem, const GlobalID *elemIDs, const int *const *connNodes)
 
int setStiffnessMatrices (GlobalID, int, const GlobalID *, const double *const *const *, int, const int *const *)
 
int setLoadVectors (GlobalID, int, const GlobalID *, const double *const *, int, const int *const *)
 
int setMatrixStructure (int **ptColIndices, int *ptRowLengths, int **blkColIndices, int *blkRowLengths, int *ptRowsPerBlkRow)
 
int setMultCREqns (int, int, int, int **, int **, int *, int *)
 
int setPenCREqns (int, int, int, int **, int **, int *)
 
int resetMatrixAndVector (double s)
 
int resetMatrix (double s)
 
int resetRHSVector (double s)
 
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)
 
int sumIntoSystemMatrix (int numPtRows, const int *ptRows, int numPtCols, const int *ptCols, const double *const *values)
 
int putIntoSystemMatrix (int numPtRows, const int *ptRows, int numPtCols, const int *ptCols, const double *const *values)
 
int getMatrixRowLength (int row, int &length)
 
int getMatrixRow (int row, double *coefs, int *indices, int len, int &rowLength)
 
int sumIntoRHSVector (int num, const double *values, const int *indices)
 
int putIntoRHSVector (int num, const double *values, const int *indices)
 
int getFromRHSVector (int num, double *values, const int *indices)
 
int matrixLoadComplete ()
 
int enforceEssentialBC (int *globalEqn, double *alpha, double *gamma, int len)
 
int enforceBlkEssentialBC (int *blkEqn, int *blkOffset, double *alpha, double *gamma, int len)
 
int enforceRemoteEssBCs (int numEqns, int *globalEqns, int **colIndices, int *colIndLen, double **coefs)
 
int enforceBlkRemoteEssBCs (int numEqns, int *blkEqns, int **blkColInds, int **blkColOffsets, int *blkColLens, double **remEssBCCoefs)
 
int getMatrixPtr (Data &data)
 
int copyInMatrix (double scalar, const Data &data)
 
int copyOutMatrix (double scalar, Data &data)
 
int sumInMatrix (double scalar, const Data &data)
 
int getRHSVectorPtr (Data &data)
 
int copyInRHSVector (double scalar, const Data &data)
 
int copyOutRHSVector (double scalar, Data &data)
 
int sumInRHSVector (double scalar, const Data &data)
 
int destroyMatrixData (Data &data)
 
int destroyVectorData (Data &data)
 
int setNumRHSVectors (int numRHSs, const int *rhsIDs)
 
int setRHSID (int rhsID)
 
int putInitialGuess (const int *eqnNumbers, const double *values, int len)
 
int getSolution (double *answers, int len)
 
int getSolnEntry (int eqnNumber, double &answer)
 
int formResidual (double *values, int len)
 
int launchSolver (int &solveStatus, int &iterations)
 
int putNodalFieldData (int, int, int *, int, const double *)
 
int writeSystem (const char *name)
 
double * getMatrixBeginPointer ()
 
int getMatrixOffset (int row, int col)
 
- Public Member Functions inherited from LinearSystemCore
 LinearSystemCore ()
 
virtual ~LinearSystemCore ()
 
virtual int getProperty (const char *, double &)
 

Private Member Functions

int createMiscStuff ()
 
int allocateMatrix (int **ptColIndices, int *ptRowLengths, int **blkColIndices, int *blkRowLengths, int *ptRowsPerBlkRow)
 
int VBRmatPlusScaledMat (AztecDVBR_Matrix *A, double scalar, AztecDVBR_Matrix *source)
 
int MSRmatPlusScaledMat (AztecDMSR_Matrix *A, double scalar, AztecDMSR_Matrix *source)
 
int createBlockMatrix (int **blkColIndices, int *blkRowLengths, int *ptRowsPerBlkRow)
 
int sumIntoBlockRow (int numBlkRows, const int *blkRows, int numBlkCols, const int *blkCols, const double *const *values, int numPtCols, bool overwriteInsteadOfAccumulate)
 
int copyBlockRow (int i, const int *blkRows, int numBlkCols, const int *blkCols, const double *const *values, double *coefs)
 
int modifyRHSforBCs ()
 
int explicitlySetDirichletBCs ()
 
int blockRowToPointRow (int blkRow)
 
int getBlockRow (int blkRow, double *&val, int &valLen, int *&blkColInds, int &blkColIndLen, int &numNzBlks, int &numNNZ)
 
int getBlkEqnsAndOffsets (int *ptEqns, int *blkEqns, int *blkOffsets, int numEqns)
 
int getBlockSize (int blkInd)
 
int sumIntoPointRow (int numPtRows, const int *ptRows, int numPtCols, const int *ptColIndices, const double *const *values, bool overwriteInsteadOfAccumulate)
 
int sumPointIntoBlockRow (int blkRow, int rowOffset, int blkCol, int colOffset, double value)
 
int setMatrixType (const char *name)
 
int selectSolver (const char *name)
 
int selectPreconditioner (const char *name)
 
void setSubdomainSolve (const char *name)
 
void setScalingOption (const char *param)
 
void setConvTest (const char *param)
 
void setPreCalc (const char *param)
 
void setTypeOverlap (const char *param)
 
void setOverlap (const char *param)
 
void setOrthog (const char *param)
 
void setAuxVec (const char *param)
 
void setAZ_output (const char *param)
 
void recordUserParams ()
 
void checkForParam (const char *paramName, int numParams_, char **paramStrings, double &param)
 
void recordUserOptions ()
 
void checkForOption (const char *paramName, int numParams_, char **paramStrings, int &param)
 
int blkRowEssBCMod (int blkEqn, int blkOffset, double *val, int *blkCols, int numCols, int numPtNNZ, double alpha, double gamma)
 
int blkColEssBCMod (int blkRow, int blkEqn, int blkOffset, double *val, int *blkCols, int numCols, int numPtNNZ, double alpha, double gamma)
 
void setDebugOutput (const char *path, const char *name)
 
void debugOutput (const char *msg) const
 
int writeA (const char *name)
 
int writeVec (Aztec_LSVector *v, const char *name)
 
int messageAbort (const char *msg) const
 

Private Attributes

MPI_Comm comm_
 
Lookuplookup_
 
bool haveLookup_
 
int numProcs_
 
int thisProc_
 
int masterProc_
 
int * update_
 
fei::SharedPtr< Aztec_Mapmap_
 
AztecDMSR_Matrix * A_
 
AztecDMSR_Matrix * A_ptr_
 
Aztec_LSVectorx_
 
Aztec_LSVector ** b_
 
Aztec_LSVectorbc_
 
int * essBCindices_
 
int numEssBCs_
 
bool bcsLoaded_
 
bool explicitDirichletBCs_
 
bool BCenforcement_no_column_mod_
 
Aztec_LSVectorb_ptr_
 
bool matrixAllocated_
 
bool vectorsAllocated_
 
bool blkMatrixAllocated_
 
bool matrixLoaded_
 
bool rhsLoaded_
 
bool needNewPreconditioner_
 
bool tooLateToChooseBlock_
 
bool blockMatrix_
 
fei::SharedPtr< Aztec_BlockMapblkMap_
 
AztecDVBR_MatrixblkA_
 
AztecDVBR_MatrixblkA_ptr_
 
int * blkUpdate_
 
AZ_MATRIX * azA_
 
AZ_PRECOND * azP_
 
bool precondCreated_
 
AZ_SCALING * azS_
 
bool scalingCreated_
 
int * aztec_options_
 
double * aztec_params_
 
double * aztec_status_
 
double * tmp_x_
 
bool tmp_x_touched_
 
double ** tmp_b_
 
double * tmp_bc_
 
bool tmp_b_allocated_
 
bool ML_Vanek_
 
int numLevels_
 
int * rhsIDs_
 
int numRHSs_
 
int currentRHS_
 
int numGlobalEqns_
 
int localOffset_
 
int numLocalEqns_
 
int numGlobalEqnBlks_
 
int localBlkOffset_
 
int numLocalEqnBlks_
 
int * localBlockSizes_
 
int numNonzeroBlocks_
 
int outputLevel_
 
int numParams_
 
char ** paramStrings_
 
std::string name_
 
int debugOutput_
 
int debugFileCounter_
 
char * debugPath_
 
char * debugFileName_
 
FILE * debugFile_
 
std::map< std::string, unsigned > & named_solve_counter_
 

Detailed Description

Definition at line 70 of file fei_Aztec_LinSysCore.hpp.

Constructor & Destructor Documentation

fei_trilinos::Aztec_LinSysCore::Aztec_LinSysCore ( MPI_Comm  comm)
virtual fei_trilinos::Aztec_LinSysCore::~Aztec_LinSysCore ( )
virtual

Member Function Documentation

LinearSystemCore* fei_trilinos::Aztec_LinSysCore::clone ( )
virtual

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

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::parameters ( int  numParams,
const char *const *  params 
)
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"

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::setLookup ( Lookup lookup)
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

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::setGlobalOffsets ( int  len,
int *  nodeOffsets,
int *  eqnOffsets,
int *  blkEqnOffsets 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::setConnectivities ( GlobalID  elemBlock,
int  numElements,
int  numNodesPerElem,
const GlobalID elemIDs,
const int *const *  connNodes 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::setStiffnessMatrices ( GlobalID  elemBlock,
int  numElems,
const GlobalID elemIDs,
const double *const *const *  stiff,
int  numEqnsPerElem,
const int *const *  eqnIndices 
)
inlinevirtual

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).

Implements LinearSystemCore.

Definition at line 95 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::setLoadVectors ( GlobalID  elemBlock,
int  numElems,
const GlobalID elemIDs,
const double *const *  load,
int  numEqnsPerElem,
const int *const *  eqnIndices 
)
inlinevirtual

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).

Implements LinearSystemCore.

Definition at line 103 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::setMatrixStructure ( int **  ptColIndices,
int *  ptRrowLengths,
int **  blkColIndices,
int *  blkRowLengths,
int *  ptRowsPerBlkRow 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::setMultCREqns ( int  multCRSetID,
int  numCRs,
int  numNodesPerCR,
int **  nodeNumbers,
int **  eqnNumbers,
int *  fieldIDs,
int *  multiplierEqnNumbers 
)
inlinevirtual

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.

Implements LinearSystemCore.

Definition at line 117 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::setPenCREqns ( int  penCRSetID,
int  numCRs,
int  numNodesPerCR,
int **  nodeNumbers,
int **  eqnNumbers,
int *  fieldIDs 
)
inlinevirtual

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).

Implements LinearSystemCore.

Definition at line 124 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::resetMatrixAndVector ( double  s)
virtual

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

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::resetMatrix ( double  s)
virtual

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

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::resetRHSVector ( double  s)
virtual

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

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::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 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::sumIntoSystemMatrix ( int  numPtRows,
const int *  ptRows,
int  numPtCols,
const int *  ptCols,
const double *const *  values 
)
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.)

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::putIntoSystemMatrix ( int  numPtRows,
const int *  ptRows,
int  numPtCols,
const int *  ptCols,
const double *const *  values 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::getMatrixRowLength ( int  row,
int &  length 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::getMatrixRow ( int  row,
double *  coefs,
int *  indices,
int  len,
int &  rowLength 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::sumIntoRHSVector ( int  num,
const double *  values,
const int *  indices 
)
virtual

For accumulating coefficients into the rhs vector

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::putIntoRHSVector ( int  num,
const double *  values,
const int *  indices 
)
virtual

For putting coefficients into the rhs vector

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::getFromRHSVector ( int  num,
double *  values,
const int *  indices 
)
virtual

For getting coefficients out of the rhs vector

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::matrixLoadComplete ( )
virtual

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

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::enforceEssentialBC ( int *  globalEqn,
double *  alpha,
double *  gamma,
int  len 
)
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

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::enforceBlkEssentialBC ( int *  blkEqn,
int *  blkOffset,
double *  alpha,
double *  gamma,
int  len 
)
int fei_trilinos::Aztec_LinSysCore::enforceRemoteEssBCs ( int  numEqns,
int *  globalEqns,
int **  colIndices,
int *  colIndLen,
double **  coefs 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::enforceBlkRemoteEssBCs ( int  numEqns,
int *  blkEqns,
int **  blkColInds,
int **  blkColOffsets,
int *  blkColLens,
double **  remEssBCCoefs 
)
int fei_trilinos::Aztec_LinSysCore::getMatrixPtr ( Data data)
virtual

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

Parameters
dataSee Data class documentation.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::copyInMatrix ( double  scalar,
const Data data 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::copyOutMatrix ( double  scalar,
Data data 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::sumInMatrix ( double  scalar,
const Data data 
)
virtual

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

Parameters
scalar
dataSee documentation for Data class.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::getRHSVectorPtr ( Data data)
virtual

Same semantics as getMatrixPtr, but applied to rhs vector.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::copyInRHSVector ( double  scalar,
const Data data 
)
virtual

Same semantics as copyInMatrix, but applied to rhs vector.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::copyOutRHSVector ( double  scalar,
Data data 
)
virtual

Same semantics as copyOutMatrix, but applied to rhs vector.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::sumInRHSVector ( double  scalar,
const Data data 
)
virtual

Same semantics as sumInMatrix, but applied to rhs vector.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::destroyMatrixData ( Data data)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::destroyVectorData ( Data data)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::setNumRHSVectors ( int  numRHSs,
const int *  rhsIDs 
)
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).

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::setRHSID ( int  rhsID)
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

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::putInitialGuess ( const int *  eqnNumbers,
const double *  values,
int  len 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::getSolution ( double *  answers,
int  len 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::getSolnEntry ( int  eqnNumber,
double &  answer 
)
virtual

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

Parameters
eqnNumberGlobal 0-based equation number.
answer

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::formResidual ( double *  values,
int  len 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::launchSolver ( int &  solveStatus,
int &  iterations 
)
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.

Implements LinearSystemCore.

int fei_trilinos::Aztec_LinSysCore::putNodalFieldData ( int  fieldID,
int  fieldSize,
int *  nodeNumbers,
int  numNodes,
const double *  data 
)
inlinevirtual

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')

Implements LinearSystemCore.

Definition at line 268 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::writeSystem ( const char *  name)
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.

Implements LinearSystemCore.

double* fei_trilinos::Aztec_LinSysCore::getMatrixBeginPointer ( )
inlinevirtual

Reimplemented from LinearSystemCore.

Definition at line 273 of file fei_Aztec_LinSysCore.hpp.

References A_ptr_.

int fei_trilinos::Aztec_LinSysCore::getMatrixOffset ( int  row,
int  col 
)
inlinevirtual

Reimplemented from LinearSystemCore.

Definition at line 275 of file fei_Aztec_LinSysCore.hpp.

References A_ptr_.

int fei_trilinos::Aztec_LinSysCore::createMiscStuff ( )
private
int fei_trilinos::Aztec_LinSysCore::allocateMatrix ( int **  ptColIndices,
int *  ptRowLengths,
int **  blkColIndices,
int *  blkRowLengths,
int *  ptRowsPerBlkRow 
)
private
int fei_trilinos::Aztec_LinSysCore::VBRmatPlusScaledMat ( AztecDVBR_Matrix A,
double  scalar,
AztecDVBR_Matrix source 
)
private
int fei_trilinos::Aztec_LinSysCore::MSRmatPlusScaledMat ( AztecDMSR_Matrix *  A,
double  scalar,
AztecDMSR_Matrix *  source 
)
private
int fei_trilinos::Aztec_LinSysCore::createBlockMatrix ( int **  blkColIndices,
int *  blkRowLengths,
int *  ptRowsPerBlkRow 
)
private
int fei_trilinos::Aztec_LinSysCore::sumIntoBlockRow ( int  numBlkRows,
const int *  blkRows,
int  numBlkCols,
const int *  blkCols,
const double *const *  values,
int  numPtCols,
bool  overwriteInsteadOfAccumulate 
)
private
int fei_trilinos::Aztec_LinSysCore::copyBlockRow ( int  i,
const int *  blkRows,
int  numBlkCols,
const int *  blkCols,
const double *const *  values,
double *  coefs 
)
private
int fei_trilinos::Aztec_LinSysCore::modifyRHSforBCs ( )
private
int fei_trilinos::Aztec_LinSysCore::explicitlySetDirichletBCs ( )
private
int fei_trilinos::Aztec_LinSysCore::blockRowToPointRow ( int  blkRow)
private
int fei_trilinos::Aztec_LinSysCore::getBlockRow ( int  blkRow,
double *&  val,
int &  valLen,
int *&  blkColInds,
int &  blkColIndLen,
int &  numNzBlks,
int &  numNNZ 
)
private
int fei_trilinos::Aztec_LinSysCore::getBlkEqnsAndOffsets ( int *  ptEqns,
int *  blkEqns,
int *  blkOffsets,
int  numEqns 
)
private
int fei_trilinos::Aztec_LinSysCore::getBlockSize ( int  blkInd)
private
int fei_trilinos::Aztec_LinSysCore::sumIntoPointRow ( int  numPtRows,
const int *  ptRows,
int  numPtCols,
const int *  ptColIndices,
const double *const *  values,
bool  overwriteInsteadOfAccumulate 
)
private
int fei_trilinos::Aztec_LinSysCore::sumPointIntoBlockRow ( int  blkRow,
int  rowOffset,
int  blkCol,
int  colOffset,
double  value 
)
private
int fei_trilinos::Aztec_LinSysCore::setMatrixType ( const char *  name)
private
int fei_trilinos::Aztec_LinSysCore::selectSolver ( const char *  name)
private
int fei_trilinos::Aztec_LinSysCore::selectPreconditioner ( const char *  name)
private
void fei_trilinos::Aztec_LinSysCore::setSubdomainSolve ( const char *  name)
private
void fei_trilinos::Aztec_LinSysCore::setScalingOption ( const char *  param)
private
void fei_trilinos::Aztec_LinSysCore::setConvTest ( const char *  param)
private
void fei_trilinos::Aztec_LinSysCore::setPreCalc ( const char *  param)
private
void fei_trilinos::Aztec_LinSysCore::setTypeOverlap ( const char *  param)
private
void fei_trilinos::Aztec_LinSysCore::setOverlap ( const char *  param)
private
void fei_trilinos::Aztec_LinSysCore::setOrthog ( const char *  param)
private
void fei_trilinos::Aztec_LinSysCore::setAuxVec ( const char *  param)
private
void fei_trilinos::Aztec_LinSysCore::setAZ_output ( const char *  param)
private
void fei_trilinos::Aztec_LinSysCore::recordUserParams ( )
private
void fei_trilinos::Aztec_LinSysCore::checkForParam ( const char *  paramName,
int  numParams_,
char **  paramStrings,
double &  param 
)
private
void fei_trilinos::Aztec_LinSysCore::recordUserOptions ( )
private
void fei_trilinos::Aztec_LinSysCore::checkForOption ( const char *  paramName,
int  numParams_,
char **  paramStrings,
int &  param 
)
private
int fei_trilinos::Aztec_LinSysCore::blkRowEssBCMod ( int  blkEqn,
int  blkOffset,
double *  val,
int *  blkCols,
int  numCols,
int  numPtNNZ,
double  alpha,
double  gamma 
)
private
int fei_trilinos::Aztec_LinSysCore::blkColEssBCMod ( int  blkRow,
int  blkEqn,
int  blkOffset,
double *  val,
int *  blkCols,
int  numCols,
int  numPtNNZ,
double  alpha,
double  gamma 
)
private
void fei_trilinos::Aztec_LinSysCore::setDebugOutput ( const char *  path,
const char *  name 
)
private
void fei_trilinos::Aztec_LinSysCore::debugOutput ( const char *  msg) const
private
int fei_trilinos::Aztec_LinSysCore::writeA ( const char *  name)
private
int fei_trilinos::Aztec_LinSysCore::writeVec ( Aztec_LSVector v,
const char *  name 
)
private
int fei_trilinos::Aztec_LinSysCore::messageAbort ( const char *  msg) const
private

Member Data Documentation

MPI_Comm fei_trilinos::Aztec_LinSysCore::comm_
private

Definition at line 372 of file fei_Aztec_LinSysCore.hpp.

Lookup* fei_trilinos::Aztec_LinSysCore::lookup_
private

Definition at line 374 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::haveLookup_
private

Definition at line 375 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::numProcs_
private

Definition at line 377 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::thisProc_
private

Definition at line 378 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::masterProc_
private

Definition at line 379 of file fei_Aztec_LinSysCore.hpp.

int* fei_trilinos::Aztec_LinSysCore::update_
private

Definition at line 381 of file fei_Aztec_LinSysCore.hpp.

fei::SharedPtr<Aztec_Map> fei_trilinos::Aztec_LinSysCore::map_
private

Definition at line 382 of file fei_Aztec_LinSysCore.hpp.

AztecDMSR_Matrix* fei_trilinos::Aztec_LinSysCore::A_
private

Definition at line 383 of file fei_Aztec_LinSysCore.hpp.

AztecDMSR_Matrix* fei_trilinos::Aztec_LinSysCore::A_ptr_
private

Definition at line 384 of file fei_Aztec_LinSysCore.hpp.

Referenced by getMatrixBeginPointer(), and getMatrixOffset().

Aztec_LSVector* fei_trilinos::Aztec_LinSysCore::x_
private

Definition at line 385 of file fei_Aztec_LinSysCore.hpp.

Aztec_LSVector ** fei_trilinos::Aztec_LinSysCore::b_
private

Definition at line 385 of file fei_Aztec_LinSysCore.hpp.

Aztec_LSVector * fei_trilinos::Aztec_LinSysCore::bc_
private

Definition at line 385 of file fei_Aztec_LinSysCore.hpp.

int* fei_trilinos::Aztec_LinSysCore::essBCindices_
private

Definition at line 386 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::numEssBCs_
private

Definition at line 387 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::bcsLoaded_
private

Definition at line 388 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::explicitDirichletBCs_
private

Definition at line 389 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::BCenforcement_no_column_mod_
private

Definition at line 390 of file fei_Aztec_LinSysCore.hpp.

Aztec_LSVector* fei_trilinos::Aztec_LinSysCore::b_ptr_
private

Definition at line 391 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::matrixAllocated_
private

Definition at line 392 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::vectorsAllocated_
private

Definition at line 393 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::blkMatrixAllocated_
private

Definition at line 394 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::matrixLoaded_
private

Definition at line 395 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::rhsLoaded_
private

Definition at line 396 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::needNewPreconditioner_
private

Definition at line 397 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::tooLateToChooseBlock_
private

Definition at line 399 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::blockMatrix_
private

Definition at line 400 of file fei_Aztec_LinSysCore.hpp.

fei::SharedPtr<Aztec_BlockMap> fei_trilinos::Aztec_LinSysCore::blkMap_
private

Definition at line 401 of file fei_Aztec_LinSysCore.hpp.

AztecDVBR_Matrix* fei_trilinos::Aztec_LinSysCore::blkA_
private

Definition at line 402 of file fei_Aztec_LinSysCore.hpp.

AztecDVBR_Matrix* fei_trilinos::Aztec_LinSysCore::blkA_ptr_
private

Definition at line 403 of file fei_Aztec_LinSysCore.hpp.

int* fei_trilinos::Aztec_LinSysCore::blkUpdate_
private

Definition at line 404 of file fei_Aztec_LinSysCore.hpp.

AZ_MATRIX* fei_trilinos::Aztec_LinSysCore::azA_
private

Definition at line 406 of file fei_Aztec_LinSysCore.hpp.

AZ_PRECOND* fei_trilinos::Aztec_LinSysCore::azP_
private

Definition at line 407 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::precondCreated_
private

Definition at line 408 of file fei_Aztec_LinSysCore.hpp.

AZ_SCALING* fei_trilinos::Aztec_LinSysCore::azS_
private

Definition at line 409 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::scalingCreated_
private

Definition at line 410 of file fei_Aztec_LinSysCore.hpp.

int* fei_trilinos::Aztec_LinSysCore::aztec_options_
private

Definition at line 412 of file fei_Aztec_LinSysCore.hpp.

double* fei_trilinos::Aztec_LinSysCore::aztec_params_
private

Definition at line 413 of file fei_Aztec_LinSysCore.hpp.

double* fei_trilinos::Aztec_LinSysCore::aztec_status_
private

Definition at line 414 of file fei_Aztec_LinSysCore.hpp.

double* fei_trilinos::Aztec_LinSysCore::tmp_x_
private

Definition at line 416 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::tmp_x_touched_
private

Definition at line 417 of file fei_Aztec_LinSysCore.hpp.

double** fei_trilinos::Aztec_LinSysCore::tmp_b_
private

Definition at line 418 of file fei_Aztec_LinSysCore.hpp.

double* fei_trilinos::Aztec_LinSysCore::tmp_bc_
private

Definition at line 419 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::tmp_b_allocated_
private

Definition at line 420 of file fei_Aztec_LinSysCore.hpp.

bool fei_trilinos::Aztec_LinSysCore::ML_Vanek_
private

Definition at line 422 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::numLevels_
private

Definition at line 423 of file fei_Aztec_LinSysCore.hpp.

int* fei_trilinos::Aztec_LinSysCore::rhsIDs_
private

Definition at line 425 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::numRHSs_
private

Definition at line 426 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::currentRHS_
private

Definition at line 428 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::numGlobalEqns_
private

Definition at line 430 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::localOffset_
private

Definition at line 431 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::numLocalEqns_
private

Definition at line 432 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::numGlobalEqnBlks_
private

Definition at line 434 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::localBlkOffset_
private

Definition at line 435 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::numLocalEqnBlks_
private

Definition at line 436 of file fei_Aztec_LinSysCore.hpp.

int* fei_trilinos::Aztec_LinSysCore::localBlockSizes_
private

Definition at line 437 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::numNonzeroBlocks_
private

Definition at line 439 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::outputLevel_
private

Definition at line 441 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::numParams_
private

Definition at line 442 of file fei_Aztec_LinSysCore.hpp.

char** fei_trilinos::Aztec_LinSysCore::paramStrings_
private

Definition at line 443 of file fei_Aztec_LinSysCore.hpp.

std::string fei_trilinos::Aztec_LinSysCore::name_
private

Definition at line 445 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::debugOutput_
private

Definition at line 446 of file fei_Aztec_LinSysCore.hpp.

int fei_trilinos::Aztec_LinSysCore::debugFileCounter_
private

Definition at line 447 of file fei_Aztec_LinSysCore.hpp.

char* fei_trilinos::Aztec_LinSysCore::debugPath_
private

Definition at line 448 of file fei_Aztec_LinSysCore.hpp.

char* fei_trilinos::Aztec_LinSysCore::debugFileName_
private

Definition at line 449 of file fei_Aztec_LinSysCore.hpp.

FILE* fei_trilinos::Aztec_LinSysCore::debugFile_
private

Definition at line 450 of file fei_Aztec_LinSysCore.hpp.

std::map<std::string,unsigned>& fei_trilinos::Aztec_LinSysCore::named_solve_counter_
private

Definition at line 452 of file fei_Aztec_LinSysCore.hpp.


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