FEI
Version of the Day
|
#include <fei_FiniteElementData.hpp>
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 |
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.
|
pure virtual |
For setting argc/argv style parameters.
numParams | Number of strings in the params argument |
params | A list of strings which will usually contain space-separated key-value pairs. Example: "debugOutput /usr/users/me/work_dir" |
Implemented in FEData.
|
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.
lookup | Input. Reference to an implementation of the Lookup interface |
Implemented in FEData.
|
pure virtual |
For describing the general structure of the finite-element problem that is to be assembled.
numElemBlocks | Number 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. |
numElemsPerBlock | List of length numElemBlocks. |
numNodesPerElem | List of length numElemBlocks. |
elemMatrixSizePerBlock | List of length numElemBlocks, i-th entry gives the number of scalar degrees-of-freedom in element-matrices for the i-th element-block. |
totalNumNodes | Total number of local nodes. |
numSharedNodes | Number of nodes that are shared with other processors. |
Implemented in FEData.
|
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.
elemBlockID | Identifier for the element-block that this element belongs to. |
elemID | Locally zero-based identifier for this element. |
numNodes | Number of nodes for this element. |
nodeNumbers | List of length numNodes. |
numDofPerNode | List of length numNodes. |
dof_ids | List of length sum(numDofPerNode[i]) |
Implemented in FEData.
|
pure virtual |
For passing element-stiffness arrays.
elemBlockID | Identifier for the element-block that these elements belong to. |
elemID | Locally zero-based identifier for this element. |
numNodes | Number of nodes on this element. |
nodeNumbers | List of length numNodes |
numDofPerNode | List of length numNodes. |
dof_ids | List of length sum(numDofPerNode[i]) |
coefs | C-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.
|
pure virtual |
For passing element-load vectors.
elemBlockID | Identifier for the element-block that this element belongs to. |
elemID | Locally zero-based identifier for this element. |
numNodes | Number of nodes on this element. |
nodeNumbers | |
numDofPerNode | |
dof_ids | |
coefs | Packed list, length sum(dofPerNode[i]). |
Implemented in FEData.
|
pure virtual |
Specify dirichlet boundary-condition values.
numBCs | Number of boundary-condition values. |
nodeNumbers | List of length numBCs. |
dof_ids | List of length numBCs. |
values | List of length numBCs. |
Implemented in FEData.
|
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.
numRowNodes | |
rowNodeNumbers | List of length numRowNodes |
row_dof_ids | List of length numRowNodes |
numColNodesPerRow | List of length numRowNodes |
colNodeNumbers | List of length sum(numColNodesPerRow[i]) |
col_dof_ids | List of length sum(numColNodesPerRow[i]) |
coefs | List of length sum(numColNodesPerRow[i]) |
Implemented in FEData.
|
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.
Implemented in FEData.
|
pure virtual |
Function called to request the launching of the linear solver.
solveStatus | Output, should indicate the status of the solve. A successful solve is usually indicated by a value of 0. |
iterations | Output, how many iterations were performed. |
Implemented in FEData.
|
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.
|
pure virtual |
Function to signal that lagrange multiplier constraints should be deleted.
Implemented in FEData.
|
pure virtual |
Get solution value from the underlying solver.
nodeNumber | |
dof_id | |
value |
Implemented in FEData.
|
pure virtual |
Get solution value of a Lagrange Multiplier from the solver.
Implemented in FEData.
|
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.
fieldID | Identifier 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. |
nodeNumbers | List of nodes for which data is being supplied. |
numNodes | |
data | List of length numNodes * (size of field 'fieldID') |
Implemented in FEData.
|
pure virtual |
Specify lagrange-multipler constraint-relation.
Implemented in FEData.
|
pure virtual |
Specify penalty constraint-relation.
Implemented in FEData.