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

#include <fei_FiniteElementData.hpp>

Inheritance diagram for FiniteElementData:
Inheritance graph
[legend]

Public Member Functions

virtual int parameters (int numParams, char **params)=0
 
virtual int setLookup (Lookup &lookup)=0
 
virtual int describeStructure (int numElemBlocks, const int *numElemsPerBlock, const int *numNodesPerElem, const int *elemMatrixSizePerBlock, int totalNumNodes, int numSharedNodes, int numMultCRs)=0
 
virtual int setConnectivity (int elemBlockID, int elemID, int numNodes, const int *nodeNumbers, const int *numDofPerNode, const int *dof_ids)=0
 
virtual int setElemMatrix (int elemBlockID, int elemID, int numNodes, const int *nodeNumbers, const int *numDofPerNode, const int *dof_ids, const double *const *coefs)=0
 
virtual int setElemVector (int elemBlockID, int elemID, int numNodes, const int *nodeNumbers, const int *numDofPerNode, const int *dof_ids, const double *coefs)=0
 
virtual int setDirichletBCs (int numBCs, const int *nodeNumbers, const int *dof_ids, const double *values)=0
 
virtual int sumIntoMatrix (int numRowNodes, const int *rowNodeNumbers, const int *row_dof_ids, const int *numColNodesPerRow, const int *colNodeNumbers, const int *col_dof_ids, const double *coefs)=0
 
virtual int loadComplete ()=0
 
virtual int launchSolver (int &solveStatus, int &iterations)=0
 
virtual int reset ()=0
 
virtual int deleteConstraints ()=0
 
virtual int getSolnEntry (int nodeNumber, int dof_id, double &value)=0
 
virtual int getMultiplierSoln (int CRID, double &lagrangeMultiplier)=0
 
virtual int putNodalFieldData (int fieldID, int fieldSize, int numNodes, const int *nodeNumbers, const double *coefs)=0
 
virtual int setMultiplierCR (int CRID, int numNodes, const int *nodeNumbers, const int *dof_ids, const double *coefWeights, double rhsValue)=0
 
virtual int setPenaltyCR (int CRID, int numNodes, const int *nodeNumbers, const int *dof_ids, const double *coefWeights, double penaltyValue, double rhsValue)=0
 

Detailed Description

This interface is used to pass finite-element data from the FEI implementation to a solver library. This is primarily intended to be an internal interface in the FEI implementation, acting as an abstraction between the FEI and underlying solvers. This interface is an alternative to the LinearSystemCore and other purely algebraic interfaces. The purely algebraic interfaces accept almost entirely algebraic data, e.g., equation-numbered matrix/vector contributions, whereas this FiniteElementData interface accepts data in terms of nodeNumbers and degree-of-freedom ids. Equation numbers don't appear in this interface.

elemBlockID's and elemID's are not necessarily globally unique. In certain problems (e.g., multi-physics, in parallel) element- blocks may span processor boundaries, in which case the implementer of FiniteElementData may assume that elemBlockID's are globally unique. Elements are never shared between processors, so elemID's don't need to be globally unique. They may be assumed to be locally zero-based numbers.

nodeNumbers are globally unique, and are globally contiguous and zero-based.

dof_ids are contiguous and zero-based. The total set of solution-fields for the problem is mapped to a contiguous, zero-based set of dof_ids. There is a unique dof_id for each scalar component of each solution-field.

Definition at line 33 of file fei_FiniteElementData.hpp.

Member Function Documentation

virtual int FiniteElementData::parameters ( int  numParams,
char **  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"

Implemented in FEData.

virtual int FiniteElementData::setLookup ( Lookup lookup)
pure virtual

Supply the FiniteElementData implementation with an object (created and owned by the FEI layer) 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

Implemented in FEData.

virtual int FiniteElementData::describeStructure ( int  numElemBlocks,
const int *  numElemsPerBlock,
const int *  numNodesPerElem,
const int *  elemMatrixSizePerBlock,
int  totalNumNodes,
int  numSharedNodes,
int  numMultCRs 
)
pure virtual

For describing the general structure of the finite-element problem that is to be assembled.

Parameters
numElemBlocksNumber of element-blocks. An element-block is a collection of elements having the same number of nodes and the same pattern of num-degrees-of-freedom-per-node.
numElemsPerBlockList of length numElemBlocks.
numNodesPerElemList of length numElemBlocks.
elemMatrixSizePerBlockList of length numElemBlocks, i-th entry gives the number of scalar degrees-of-freedom in element-matrices for the i-th element-block.
totalNumNodesTotal number of local nodes.
numSharedNodesNumber of nodes that are shared with other processors.
Returns
error-code 0 if successful

Implemented in FEData.

virtual int FiniteElementData::setConnectivity ( int  elemBlockID,
int  elemID,
int  numNodes,
const int *  nodeNumbers,
const int *  numDofPerNode,
const int *  dof_ids 
)
pure virtual

For passing element-connectivity arrays. dof_ids are contiguous and zero-based. The total set of solution-fields for the problem is mapped to a contiguous, zero-based set of dof_ids. There is a unique dof_id for each scalar component of each solution-field.

Parameters
elemBlockIDIdentifier for the element-block that this element belongs to.
elemIDLocally zero-based identifier for this element.
numNodesNumber of nodes for this element.
nodeNumbersList of length numNodes.
numDofPerNodeList of length numNodes.
dof_idsList of length sum(numDofPerNode[i])

Implemented in FEData.

virtual int FiniteElementData::setElemMatrix ( int  elemBlockID,
int  elemID,
int  numNodes,
const int *  nodeNumbers,
const int *  numDofPerNode,
const int *  dof_ids,
const double *const *  coefs 
)
pure virtual

For passing element-stiffness arrays.

Parameters
elemBlockIDIdentifier for the element-block that these elements belong to.
elemIDLocally zero-based identifier for this element.
numNodesNumber of nodes on this element.
nodeNumbersList of length numNodes
numDofPerNodeList of length numNodes.
dof_idsList of length sum(numDofPerNode[i])
coefsC-style table (list of pointers). Each row of the table is of length sum(dofPerNode[i]), and that is also the number of rows.

Implemented in FEData.

virtual int FiniteElementData::setElemVector ( int  elemBlockID,
int  elemID,
int  numNodes,
const int *  nodeNumbers,
const int *  numDofPerNode,
const int *  dof_ids,
const double *  coefs 
)
pure virtual

For passing element-load vectors.

Parameters
elemBlockIDIdentifier for the element-block that this element belongs to.
elemIDLocally zero-based identifier for this element.
numNodesNumber of nodes on this element.
nodeNumbers
numDofPerNode
dof_ids
coefsPacked list, length sum(dofPerNode[i]).

Implemented in FEData.

virtual int FiniteElementData::setDirichletBCs ( int  numBCs,
const int *  nodeNumbers,
const int *  dof_ids,
const double *  values 
)
pure virtual

Specify dirichlet boundary-condition values.

Parameters
numBCsNumber of boundary-condition values.
nodeNumbersList of length numBCs.
dof_idsList of length numBCs.
valuesList of length numBCs.

Implemented in FEData.

virtual int FiniteElementData::sumIntoMatrix ( int  numRowNodes,
const int *  rowNodeNumbers,
const int *  row_dof_ids,
const int *  numColNodesPerRow,
const int *  colNodeNumbers,
const int *  col_dof_ids,
const double *  coefs 
)
pure virtual

Sum coefficients into the matrix. This may be used to pass a dense rectangular block of coefs, or a diagonal band of coefs, or an arbitrary sparse block of coefs.

Parameters
numRowNodes
rowNodeNumbersList of length numRowNodes
row_dof_idsList of length numRowNodes
numColNodesPerRowList of length numRowNodes
colNodeNumbersList of length sum(numColNodesPerRow[i])
col_dof_idsList of length sum(numColNodesPerRow[i])
coefsList of length sum(numColNodesPerRow[i])

Implemented in FEData.

virtual int FiniteElementData::loadComplete ( )
pure virtual

Function called to signal to the FiniteElementData implementation that data-loading is complete and any synchronization or other final operations may be performed now. This is a collective function, must be called by all processors.

Returns
error-code 0 if successful

Implemented in FEData.

virtual int FiniteElementData::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.

Implemented in FEData.

virtual int FiniteElementData::reset ( )
pure virtual

Function to signal that all coefficient values should be zero'd in preparation for a new assemble/solve phase. This is to handle cases where no structural information is changing from one time-step to the next, but new coefficient values need to be loaded for the next solve.

Implemented in FEData.

virtual int FiniteElementData::deleteConstraints ( )
pure virtual

Function to signal that lagrange multiplier constraints should be deleted.

Implemented in FEData.

virtual int FiniteElementData::getSolnEntry ( int  nodeNumber,
int  dof_id,
double &  value 
)
pure virtual

Get solution value from the underlying solver.

Parameters
nodeNumber
dof_id
value

Implemented in FEData.

virtual int FiniteElementData::getMultiplierSoln ( int  CRID,
double &  lagrangeMultiplier 
)
pure virtual

Get solution value of a Lagrange Multiplier from the solver.

Implemented in FEData.

virtual int FiniteElementData::putNodalFieldData ( int  fieldID,
int  fieldSize,
int  numNodes,
const int *  nodeNumbers,
const double *  coefs 
)
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 by the FEI implementation.
nodeNumbersList of nodes for which data is being supplied.
numNodes
dataList of length numNodes * (size of field 'fieldID')

Implemented in FEData.

virtual int FiniteElementData::setMultiplierCR ( int  CRID,
int  numNodes,
const int *  nodeNumbers,
const int *  dof_ids,
const double *  coefWeights,
double  rhsValue 
)
pure virtual

Specify lagrange-multipler constraint-relation.

Implemented in FEData.

virtual int FiniteElementData::setPenaltyCR ( int  CRID,
int  numNodes,
const int *  nodeNumbers,
const int *  dof_ids,
const double *  coefWeights,
double  penaltyValue,
double  rhsValue 
)
pure virtual

Specify penalty constraint-relation.

Implemented in FEData.


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