FEI Package Browser (Single Doxygen Collection)
Version of the Day
|
#include <FEI.hpp>
Public Member Functions | |
virtual | ~FEI () |
virtual int | parameters (int numParams, const char *const *paramStrings)=0 |
virtual int | setIDLists (int numMatrices, const int *matrixIDs, int numRHSs, const int *rhsIDs)=0 |
virtual int | setSolveType (int solveType)=0 |
virtual int | initFields (int numFields, const int *fieldSizes, const int *fieldIDs, const int *fieldTypes=NULL)=0 |
virtual int | initElemBlock (GlobalID elemBlockID, int numElements, int numNodesPerElement, const int *numFieldsPerNode, const int *const *nodalFieldIDs, int numElemDofFieldsPerElement, const int *elemDOFFieldIDs, int interleaveStrategy)=0 |
virtual int | initElem (GlobalID elemBlockID, GlobalID elemID, const GlobalID *elemConn)=0 |
virtual int | initSharedNodes (int numSharedNodes, const GlobalID *sharedNodeIDs, const int *numProcsPerNode, const int *const *sharingProcIDs)=0 |
virtual int | initCRMult (int numCRNodes, const GlobalID *CRNodeIDs, const int *CRFieldIDs, int &CRID)=0 |
virtual int | initCRPen (int numCRNodes, const GlobalID *CRNodeIDs, const int *CRFieldIDs, int &CRID)=0 |
virtual int | initSlaveVariable (GlobalID slaveNodeID, int slaveFieldID, int offsetIntoSlaveField, int numMasterNodes, const GlobalID *masterNodeIDs, const int *masterFieldIDs, const double *weights, double rhsValue)=0 |
virtual int | initComplete ()=0 |
virtual int | setCurrentMatrix (int matrixID)=0 |
virtual int | setCurrentRHS (int rhsID)=0 |
virtual int | resetSystem (double s=0.0)=0 |
virtual int | resetMatrix (double s=0.0)=0 |
virtual int | resetRHSVector (double s=0.0)=0 |
virtual int | resetInitialGuess (double s)=0 |
virtual int | deleteMultCRs ()=0 |
virtual int | loadNodeBCs (int numNodes, const GlobalID *nodeIDs, int fieldID, const int *offsetsIntoField, const double *prescribedValues)=0 |
virtual int | loadElemBCs (int numElems, const GlobalID *elemIDs, int fieldID, const double *const *alpha, const double *const *beta, const double *const *gamma)=0 |
virtual int | sumInElem (GlobalID elemBlockID, GlobalID elemID, const GlobalID *elemConn, const double *const *elemStiffness, const double *elemLoad, int elemFormat)=0 |
virtual int | sumInElemMatrix (GlobalID elemBlockID, GlobalID elemID, const GlobalID *elemConn, const double *const *elemStiffness, int elemFormat)=0 |
virtual int | sumInElemRHS (GlobalID elemBlockID, GlobalID elemID, const GlobalID *elemConn, const double *elemLoad)=0 |
virtual int | loadCRMult (int CRMultID, int numCRNodes, const GlobalID *CRNodeIDs, const int *CRFieldIDs, const double *CRWeights, double CRValue)=0 |
virtual int | loadCRPen (int CRPenID, int numCRNodes, const GlobalID *CRNodeIDs, const int *CRFieldIDs, const double *CRWeights, double CRValue, double penValue)=0 |
virtual int | putIntoRHS (int IDType, int fieldID, int numIDs, const GlobalID *IDs, const double *coefficients)=0 |
virtual int | sumIntoRHS (int IDType, int fieldID, int numIDs, const GlobalID *IDs, const double *coefficients)=0 |
virtual int | sumIntoMatrixDiagonal (int, int, int, const GlobalID *, const double *) |
virtual int | setMatScalars (int numScalars, const int *IDs, const double *scalars)=0 |
virtual int | setRHSScalars (int numScalars, const int *IDs, const double *scalars)=0 |
virtual int | loadComplete (bool applyBCs=true, bool globalAssemble=true)=0 |
virtual int | residualNorm (int whichNorm, int numFields, int *fieldIDs, double *norms)=0 |
virtual int | solve (int &status)=0 |
virtual int | iterations (int &itersTaken) const =0 |
virtual int | getFieldSize (int fieldID, int &numScalars)=0 |
virtual int | getEqnNumbers (GlobalID ID, int idType, int fieldID, int &numEqns, int *eqnNumbers)=0 |
virtual int | getNodalFieldSolution (int fieldID, int numNodes, const GlobalID *nodeIDs, double *results)=0 |
virtual int | getNumLocalNodes (int &numNodes)=0 |
virtual int | getLocalNodeIDList (int &numNodes, GlobalID *nodeIDs, int lenNodeIDs)=0 |
virtual int | putNodalFieldData (int fieldID, int numNodes, const GlobalID *nodeIDs, const double *data)=0 |
virtual int | getBlockNodeSolution (GlobalID elemBlockID, int numNodes, const GlobalID *nodeIDs, int *offsets, double *results)=0 |
virtual int | getNodalSolution (int numNodes, const GlobalID *nodeIDs, int *offsets, double *results)=0 |
virtual int | getBlockFieldNodeSolution (GlobalID elemBlockID, int fieldID, int numNodes, const GlobalID *nodeIDs, double *results)=0 |
virtual int | getBlockElemSolution (GlobalID elemBlockID, int numElems, const GlobalID *elemIDs, int &numElemDOFPerElement, double *results)=0 |
virtual int | getNumCRMultipliers (int &numMultCRs)=0 |
virtual int | getCRMultIDList (int numMultCRs, int *multIDs)=0 |
virtual int | getCRMultipliers (int numCRs, const int *CRIDs, double *results)=0 |
virtual int | putBlockNodeSolution (GlobalID elemBlockID, int numNodes, const GlobalID *nodeIDs, const int *offsets, const double *estimates)=0 |
virtual int | putBlockFieldNodeSolution (GlobalID elemBlockID, int fieldID, int numNodes, const GlobalID *nodeIDs, const double *estimates)=0 |
virtual int | putBlockElemSolution (GlobalID elemBlockID, int numElems, const GlobalID *elemIDs, int dofPerElem, const double *estimates)=0 |
virtual int | putCRMultipliers (int numMultCRs, const int *CRMultIDs, const double *multEstimates)=0 |
virtual int | getBlockNodeIDList (GlobalID elemBlockID, int numNodes, GlobalID *nodeIDs)=0 |
virtual int | getBlockElemIDList (GlobalID elemBlockID, int numElems, GlobalID *elemIDs)=0 |
virtual int | version (const char *&versionString)=0 |
virtual int | cumulative_cpu_times (double &initPhase, double &loadPhase, double &solve, double &solnReturn)=0 |
virtual int | getNumSolnParams (GlobalID globalNodeID, int &numSolnParams) const =0 |
virtual int | getNumElemBlocks (int &numElemBlocks) const =0 |
virtual int | getNumBlockActNodes (GlobalID elemBlockID, int &numNodes) const =0 |
virtual int | getNumBlockActEqns (GlobalID elemBlockID, int &numEqns) const =0 |
virtual int | getNumNodesPerElement (GlobalID elemBlockID, int &nodesPerElem) const =0 |
virtual int | getNumEqnsPerElement (GlobalID elemBlockID, int &eqnsPerElem) const =0 |
virtual int | getNumBlockElements (GlobalID blockID, int &numElems) const =0 |
virtual int | getNumBlockElemDOF (GlobalID blockID, int &DOFPerElem) const =0 |
public Finite Element Interface specification, version 2.1, in C++.
Abstract base class.
Note: all FEI functions return an int error code. A value of 0 indicates that there was no error. Errors are usually indicated by -1, except where noted in the documentation for a particular function.
|
pure virtual |
Set parameters associated with solver choice, etc. This function may be called at any time after the FEI object is instantiated. This function may be called repeatedly with different parameters, which will accumulate those parameters into an internal 'master'-list of parameters.
numParams | Number of parameters being supplied. |
paramStrings | List of 'numParams' strings. Each string usually contains a key-value pair, separated by a space. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by beam_oldfei_main(), driverData::call_fei_method(), feiDriver_main(), main(), poisson_main(), and FEI_tester::testInitialization().
|
pure virtual |
Specify matrixIDs and rhsIDs to be used in cases where multiple matrices and/or rhs vectors are being assembled. This function does not need to be called if only one matrix and rhs are being assembled. Note: the values of the matrix and rhs identifiers must be non-negative. Important Note: If this function is called, it must be called BEFORE setSolveType is called. setSolveType must then be called with the parameter FEI_AGGREGATE_SUM (eigen-solves and product-solves aren't supported yet).
numMatrices | length of matrixIDs parameter |
matrixIDs | list of user-defined identifiers for separate matrices to be assembled |
numRHSs | length of rhsIDs parameter |
rhsIDs | list of user-defined identifiers for separate rhs vectors to be assembled |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method(), and FEI_tester::setIDlists().
|
pure virtual |
Set the type of solve to be performed. This distinguishes between a 'standard' single solve of Ax=b, an eigen-solve (not yet supported), an 'aggregate-sum' solve (a linear-combination of several separate A's and b's), and an 'aggregate-product' solve (not supported).
solveType | currently supported values for this are: FEI_SINGLE_SOLVE, FEI_AGGREGATE_SUM |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by beam_oldfei_main(), driverData::call_fei_method(), FEI_tester::initializationPhase(), main(), and poisson_main().
|
pure virtual |
Identify all the fields present in the analysis. A field may be a scalar such as temperature or pressure, or a 3-vector for velocity, etc. Non-solution fields may be denoted by a negative fieldID. This allows for situations where the application wants to pass data that the FEI doesn't need, (such as geometric coordinates) through to the underlying linear algebra library. This may be done via the various put*Solution functions.
numFields | Global number of fields in the entire problem, on all processors. (This is the length of the fieldSizes and fieldIDs lists.) |
fieldSizes | Number of scalars contained in each field. |
fieldIDs | User-supplied identifiers for each field. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by beam_oldfei_main(), driverData::call_fei_method(), test_Factory::factory_test1(), FEI_tester::initializationPhase(), main(), and poisson_main().
|
pure virtual |
Initialize the description of an element-block. This function informs the fei implementation of the defining characteristics for a block of elements. An element-block must be homogeneous – all elements in the block must have the same number of nodes, same number of solution fields per node, etc.
elemBlockID | The user-defined identifier for this element-block. |
numElements | The number of elements in this block. |
numNodesPerElement | Length of the numFieldsPerNode list. |
numFieldsPerNode | Lengths of the rows of the nodalFieldIDs table. |
nodalFieldIDs | Table where row 'i' is the list of field ids for the ith node on every element in this element-block. |
numElemDofFieldsPerElement | Length of the elemDOFFieldIDs list. |
elemDOFFieldIDs | list of field identifiers for the element- centered degrees-of-freedom in the elements in this block. |
interleaveStrategy | Indicates the ordering of solution-components in the element-wise (e.g., stiffness) contribution arrays. Valid values are FEI_NODE_MAJOR (all field-components for first node are followed by all field-components for second node, ...) or FEI_FIELD_MAJOR (first-field for all nodes, then second-field for all nodes, ...) |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method(), init_elem_connectivities(), HexBeam_Functions::init_elem_connectivities(), and FEI_tester::initializationPhase().
|
pure virtual |
Initialize an element's connectivity. Provide a nodal connectivity list for inclusion in the sparse matrix structure being constructed.
elemBlockID | Which element-block this element belongs to. |
elemID | A user-provided identifier for this element. |
elemConn | List of nodeIDs connected to this element. Length of this list must be 'numNodesPerElement' provided to the function 'initElemBlock' for this elemBlockID. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method(), init_elem_connectivities(), HexBeam_Functions::init_elem_connectivities(), and FEI_tester::initializationPhase().
|
pure virtual |
Identify a list of nodes that are shared by processors other than the local processor. This function must be called symmetrically on sharing processors – if a node is identified to processor 0 as being shared with processor 1, then that node must also be identified to processor 1 as being shared with processor 0. This function may be called repeatedly, and shared nodes may appear in multiple calls to this function.
numSharedNodes | Length of the sharedNodeIDs and numProcsPerNode lists, and the number of rows in the sharingProcIDs table. |
sharedNodeIDs | List of nodes that are shared by this processor and one or more other processors. |
numProcsPerNode | Lengths of the rows of the sharingProcIDs table. |
sharingProcIDs | List, for each sharedNode, of processors which share that node. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method(), HexBeam_Functions::init_shared_nodes(), FEI_tester::initializationPhase(), and set_shared_nodes().
|
pure virtual |
Constraint relation initialization, Lagrange Multiplier formulation.
numCRNodes | Length of the CRNodeIDs and CRFieldIDs lists. |
CRNodeIDs | Nodes involved in this constraint relation. |
CRFieldIDs | List of the the field being constrained at each node. |
CRID | Output. An identifier by which this constraint relation may be referred to later, when loading its weight data and recovering its Lagrange Multiplier after the solve. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method(), HexBeam_Functions::init_constraints(), and FEI_tester::initializationPhase().
|
pure virtual |
Constraint relation initialization, Penalty function formulation .
numCRNodes | Length of the CRNodeIDs and CRFieldIDs lists. |
CRNodeIDs | Nodes involved in this constraint relation. |
CRFieldIDs | List of the the field being constrained at each node. |
CRID | Output. An identifier by which this constraint relation may be referred to later, when loading its weight data and penalty value. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::initializationPhase().
|
pure virtual |
Advise the FEI that a nodal variable is slaved to a linear combination of other variables, plus a right-hand-side value (note that the rhsValue will often be zero). Since a field may contain more than one scalar component, the particular scalar equation that's being slaved must be specified by not only a nodeID and fieldID, but also an offset into the slave field.
The general form of the dependency being specified is: seqn = sum ( weight_i * meqn_i ) + rhsValue where 'seqn' means slave-equation and 'meqn' means master equation.
Example: to specify that a slave-equation is the average of two master- equations: seqn = 0.5*meqn_1 + 0.5*meqn_2 + 0.0 (Where 0.0 is the rhsValue in this case.)
The list of weights supplied will be assumed to be of length sum(masterFieldSizes). i.e., the slave equation may be dependent on more than one component of the specified master field. In cases where a master field contains more than one scalar component, but only one of those components is relevant to the dependency being specified, then positions in the weights list corresponding to the non-relevant field-components should contain zeros.
This mechanism can also be used as an alternative way to specify essential boundary conditions, where the rhsValue will be the boundary condition value, with no master nodes or weights.
Note: This is a new and experimental capability, and is not compatible with all other FEI capabilities. In particular, the following precaution should be taken: Don't identify both slave variables and constraint-relations for the same degree of freedom. They are mutually exclusive.
slaveNodeID | Node identifier of the node containing the slave eqn. |
slaveFieldID | Field identifier corresponding to the slave eqn. |
offsetIntoSlaveField | Denotes location, within the field, of the slave eqn. |
numMasterNodes | Number of nodes containing variables on which the slave depends. |
masterNodeIDs | Node identifiers of the master nodes. |
masterFieldIDs | List, length numMasterNodes, of the field at each master-node which contains the scalar variable(s) on which the slave depends. |
weights | List, length sum-of-master-field-sizes, containing the weighting coefficients described above. |
rhsValue |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::initializationPhase().
|
pure virtual |
Indicate that initialization phase is complete. This function will internally finish calculating the structure of the sparse matrix, and provide that structure to the underlying linear algebra library. At that point the linear algebra library may or may not go ahead and allocate the matrix and vectors for the linear system(s).
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by beam_oldfei_main(), driverData::call_fei_method(), FEI_tester::initializationPhase(), main(), and poisson_main().
|
pure virtual |
Set current matrix data 'context'. When assembling multiple matrices, this allows the application to specify which matrix should receive data provided by any subsequent data-loading function calls.
matrixID | One of the identifiers provided earlier via the function 'setIDLists'. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::aggregateLoadPhase(), and driverData::call_fei_method().
|
pure virtual |
Set current rhs data 'context'. When assembling multiple rhs vectors, this allows the application to specify which rhs should receive data provided by any subsequent data-loading function calls.
rhsID | One of the identifiers provided earlier via the function 'setIDLists'. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::aggregateLoadPhase().
|
pure virtual |
Set a value (usually zero) througout the linear system.
s | The value to be written into the linear system. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method(), and FEI_tester::testLoading().
|
pure virtual |
Set a value (usually zero) througout the matrix.
s | The value to be written into the matrix. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method(), and FEI_tester::testLoading().
|
pure virtual |
Set a value (usually zero) througout the rhs vector.
s | The value to be written into the rhs vector. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method(), and FEI_tester::testLoading().
|
pure virtual |
Set a value (usually, if not always, 0.0) throughout the initial guess (solution) vector.
s | Input. Scalar value to use in filling the solution vector. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method().
|
pure virtual |
Request that any existing Lagrange-Multiplier constraints be deleted. (Intended to be called in preparation for loading new/different constraints.)
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method().
|
pure virtual |
Load nodal boundary condition data. This allows the application to specify a boundary condition (dirichlet) on a list of nodes.
The boundary condition specified via this function applies to the same solution field on all nodes in the nodeIDs list.
The i-th entry in the offsetsIntoField array specifies which component of the specified field will be prescribed by the i-th entry in the prescribedValues array.
numNodes | Length of the nodeIDs list. |
nodeIDs | List of nodes upon which a boundary condition is to be imposed. |
fieldID | The solution field that will receive the boundary condition. |
offsetsIntoField | Array, length numNodes. |
prescribedValues | Array, length numNodes. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::aggregateLoadPhase(), load_BC_data(), HexBeam_Functions::load_BC_data(), and FEI_tester::normalLoadPhase().
|
pure virtual |
Load boundary condition data for element-dof. The form of these boundary conditions is as follows for a given field upon which the condition is to be imposed (description courtesy of Kim Mish). If the primary field solution unknown is denoted by u, and the dual of the solution (e.g., force as opposed to displacement) is denoted by q, then a generic boundary condition can be specified by alpha*u + beta*q = gamma. A dirichlet boundary condition is given when alpha is non-zero but beta is zero. A neumann boundary condition is given when alpha is zero but beta is non-zero. A mixed boundary condition is given when alpha and beta are both non-zero.
numElems | Length of the elemIDs list. |
elemIDs | List of elements for which a boundary condition is to be specified. |
fieldID | The solution field for which to apply the boundary condition. |
alpha | Table, as in 'loadNodeBCs', but with 'numElems' number-of- rows. |
beta | Table, same dimensions as alpha. |
gamma | Table, same dimensions as alpha. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Element-stiffness/load data loading. This function accumulates element stiffness and load data into the underlying matrix and rhs vector.
elemBlockID | Which element-block this element belongs to. |
elemID | User-supplied identifier for this element. |
elemConn | Connectivity list of nodes that are connected to this element. |
elemStiffness | Table of element-stiffness data. Dimensions of this table defined by the sum of the sizes of the fields associated with this element. (This information supplied earlier via 'initElemBlock'.) |
elemLoad | Element-load vector. |
elemFormat | Designates the way in which the 'elemStiffness' stiffness-matrix data is laid out. Valid values for this parameter can be found in the file fei_defs.h. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method().
|
pure virtual |
Element-stiffness data loading. This function is the same as 'sumInElem' but only accepts stiffness data, not the load data for the rhs.
elemBlockID | Which element-block this element belongs to. |
elemID | User-supplied identifier for this element. |
elemConn | Connectivity list of nodes that are connected to this element. |
elemStiffness | Table of element-stiffness data. Dimensions of this table defined by the sum of the sizes of the fields associated with this element. (This information supplied earlier via 'initElemBlock'.) |
elemFormat | Designates the way in which the 'elemStiffness' stiffness-matrix data is laid out. Valid values for this parameter can be found in the file fei_defs.h. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::aggregateLoadPhase(), driverData::call_fei_method(), load_elem_data(), HexBeam_Functions::load_elem_data(), load_elem_data_putrhs(), and FEI_tester::normalLoadPhase().
|
pure virtual |
Element-load data loading. This function is the same as 'sumInElem', but only accepts the load for the rhs, not the stiffness matrix.
elemBlockID | Which element-block this element belongs to. |
elemID | User-supplied identifier for this element. |
elemConn | Connectivity list of nodes that are connected to this element. |
elemLoad | Element-load vector. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::aggregateLoadPhase(), driverData::call_fei_method(), load_elem_data(), HexBeam_Functions::load_elem_data(), and FEI_tester::normalLoadPhase().
|
pure virtual |
Load weight/value data for a Lagrange Multiplier constraint relation.
CRMultID | Identifier returned from an earlier call to 'initCRMult'. |
numCRNodes | Length of CRNodeIDs and CRFieldIDs lists. |
CRNodeIDs | List of nodes in this constraint relation. |
CRFieldIDs | List of fields, one per node, to be constrained. |
CRWeights | Weighting coefficients. This length of this list is the sum of the sizes associated with the fields identified in CRFieldIDs. |
CRValue | The constraint's rhs value. Often (always?) zero. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method(), HexBeam_Functions::load_constraints(), and FEI_tester::normalLoadPhase().
|
pure virtual |
Load weight/value data for a Penalty constraint relation.
CRPenID | Identifier returned from an earlier call to 'initCRPen'. |
numCRNodes | Length of CRNodeIDs and CRFieldIDs lists. |
CRNodeIDs | List of nodes in this constraint relation. |
CRFieldIDs | List of fields, one per node, to be constrained. |
CRWeights | Weighting coefficients. This length of this list is the sum of the sizes associated with the fields identified in CRFieldIDs. |
CRValue | The constraint's rhs value. Often (always?) zero. |
penValue | The penalty value. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::normalLoadPhase().
|
pure virtual |
Put a copy of coefficient data into the rhs vector.
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by load_elem_data_putrhs().
|
pure virtual |
Sum a copy of coefficient data into the rhs vector.
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
inlinevirtual |
Reimplemented in FEI_Implementation.
|
pure virtual |
Set scalars by which to multiply matrices, in cases where a linear- combination of several matrices is to be solved.
numScalars | Length of the IDs and scalars lists. |
IDs | Matrix identifiers which must be a subset of those supplied via an earlier call to 'setIDLists'. |
scalars | The coefficients by which to multiply the matrices. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::aggregateLoadPhase().
|
pure virtual |
Set scalars by which to multiply rhs vectors, in cases where a linear- combination of several rhs vectors is to be solved.
numScalars | Length of the IDs and scalars lists. |
IDs | RHS-vector identifiers which must be a subset of those supplied via an earlier call to 'setIDLists'. |
scalars | The coefficients by which to multiply the rhs vectors. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::aggregateLoadPhase().
|
pure virtual |
Indicate that all data loading is complete, and that the underlying linear system can now be "finalized". e.g., boundary conditions enforced, shared contributions exchanged among processors, etc.
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by beam_oldfei_main(), driverData::call_fei_method(), load_elem_data_putrhs(), main(), and FEI_tester::testLoading().
|
pure virtual |
Request residual norms of the underlying linear-system, broken down by solution field. This function should obviously not be called until after the matrix and vector data has been fully assembled. Calculates the residual vector r = b - Ax, then takes a norm of r, separating out the components corresponding to the various fields in 'fieldIDs'. If the solution vector x contains zeros, then r == b.
whichNorm | Determines which norm is calculated. 1 -> 1-norm, 2 -> 2-norm, 0 -> infinity-norm. |
numFields | Length of the fieldIDs and norms lists. |
fieldIDs | The fields for which residual norms are to be returned. |
norms | A norm corresponding to each field in fieldIDs. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::exerciseResidualNorm(), and poisson_main().
|
pure virtual |
Launch the solver to solve the assembled linear system.
status | Indicates the status of the solve. Varies depending on the specific underlying solver library. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by beam_oldfei_main(), driverData::call_fei_method(), main(), poisson_main(), and FEI_tester::testSolve().
|
pure virtual |
Query number of iterations taken for last solve.
itersTaken | Iterations performed during any previous solve. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by poisson_main(), and FEI_tester::testSolve().
|
pure virtual |
Query the size of a field. This info is supplied to the FEI (initFields) by the application, but may not be readily available on the app side at all times. Thus, it would be nice if the FEI could answer this query.
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by test_Factory::factory_test1().
|
pure virtual |
Since the ultimate intent for matrix-access is to bypass the FEI and go straight to the underlying data objects, we need a translation function to map between the IDs that the FEI deals in, and equation numbers that linear algebra objects deal in.
ID | Identifier of either a node or an element. |
idType | Can take either of the values FEI_NODE or FEI_ELEMENT. |
fieldID | Identifies a particular field at this [node||element]. |
numEqns | Output. Number of equations associated with this node/field (or element/field) pair. |
eqnNumbers | Caller-allocated array. On exit, this is filled with the above-described equation-numbers. They are global 0-based numbers. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Get the solution data for a particular field, on an arbitrary set of nodes.
fieldID | Input. field identifier for which solution data is being requested. |
numNodes | Input. Length of the nodeIDs list. |
nodeIDs | Input. List specifying the nodes on which solution data is being requested. |
results | Allocated by caller, but contents are output. Solution data for the i-th node/element starts in position i*fieldSize, where fieldSize is the number of scalar components that make up 'fieldID'. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by poisson_main().
|
pure virtual |
Get the number of nodes that are local to this processor (includes nodes that are shared by other processors).
numNodes | Output. Number of local nodes. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::exercisePutFunctions(), poisson_main(), and FEI_tester::save_block_node_soln().
|
pure virtual |
Get a list of the nodeIDs that are local to this processor (includes nodes that are shared by other processors).
numNodes | Output. Same as the value output by 'getNumLocalNodes'. |
nodeIDs | Caller-allocated array, contents to be filled by this function. |
lenNodeIDs | Input. Length of the caller-allocated nodeIDs array. If lenNodeIDs is less than numNodes, then only 'lenNodeIDs' nodeIDs are provided, of course. If lenNodeIDs is greater than numNodes, then only 'numNodes' positions of the nodeIDs array are referenced. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::exercisePutFunctions(), poisson_main(), and FEI_tester::save_block_node_soln().
|
pure virtual |
Input data associated with the specified field on a specified list of nodes.
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::exercisePutFunctions().
|
pure virtual |
Return nodal solutions for an element-block. The user supplies the list of nodes for which solutions are required. This list may be a subset of the nodes in the element-block.
elemBlockID | Element-block for which the nodal solutions are desired. |
numNodes | Length of the nodeIDs list. |
nodeIDs | User-supplied list of nodes in this element-block. The list of all nodes in an element-block may be obtained from the FEI via the function 'getBlockNodeIDList'. |
offsets | List of offsets into the results list. The first solution value for node i is results[offsets[i]]. The number of solution values for node i is offsets[i+1]-offsets[i]. The offsets list is of length numNodes+1. |
results | List of the nodal solution values. The results list should be allocated with length 'getNumBlockActEqns', if 'nodeIDs' contains all nodes in the element-block. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Return nodal solutions for an arbitrary list of nodes. The user supplies the list of node-identifiers for which solutions are required. This list may be any subset of the nodes that reside on the local processor.
numNodes | Number of nodes |
nodeIDs | Node-identifiers, list of length numNodes |
offsets | List of length numNodes+1, allocated by caller. On exit, this list will contain offsets into the results array at which nodal solutions are located. The solution values for the nodeIDs[i] will begin at position offsets[i]. The last solution value for nodeIDs[i] will be located at position offsets[i+1]-1. |
results | List, allocated by caller, of length sum-of-num-dof-per-node. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::save_block_node_soln().
|
pure virtual |
Return nodal solutions for one field.
elemBlockID | Element-block identifier. |
fieldID | The field for which solution values are desired. |
numNodes | Length of the nodeIDs list. |
nodeIDs | User-supplied list of nodes for which solution values are desired. The list of all nodes in an element-block may be obtained from the FEI via the function 'getBlockNodeIDList'. |
results | List of solution values. The solution values for node i start at results[i*fieldSize] where fieldSize is the size of fieldID. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Return elem-dof solution params.
elemBlockID | Element-block identifier. |
numElems | Length of the elemIDs list. |
elemIDs | User-supplied list of elements for which solution values are desired. |
numElemDOFPerElement | Output. Number of elem-dof per element. |
results | List of solution values. The solution values for element i start at results[i*numElemDOFPerElement] |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::save_block_elem_soln().
|
pure virtual |
Number of constraint relation multipliers. The number of Lagrange Multiplier constraint relations in the problem.
numMultCRs | Output |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::save_multiplier_soln().
|
pure virtual |
List of identifiers for Lagrange Multipliers.
numMultCRs | Input. Value obtained from 'getNumCRMultipliers'. |
multIDs | Output. User-allocated list, to be filled, on exit, with the identifiers for the first 'numMultCRs' multiplier constraints in the problem. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::save_multiplier_soln().
|
pure virtual |
Return Lagrange Multiplier solutions.
numCRs | number of constraint relations for which multipliers are being requested. |
CRIDs | Identifiers of constraint-relations for which multipliers are being requested. |
results | List of Lagrange Multipliers, one per constraint relation. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::save_multiplier_soln().
|
pure virtual |
Put nodal-based solution estimates for an element-block.
elemBlockID | Element-block identifier. |
numNodes | Length of nodeIDs list. |
nodeIDs | Those nodes for which solutions are being supplied. |
offsets | List, length numNodes+1, of offsets into the estimates list. |
estimates | List of solution estimates. The solution estimates for node i will be assumed to start in estimates[offsets[i]]. The length of the estimates list should be offsets[numNodes]. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Put nodal-based solution estimates for one field for an element-block.
elemBlockID | Element-block identifier. |
fieldID | The field for which estimates are being supplied. |
numNodes | Length of the nodeIDs list. |
nodeIDs | The nodes for which solution estimates are being supplied. |
estimates | List of initial guesses. Should be of length numNodes * fieldSize, where fieldSize is the size of field 'fieldID'. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by driverData::call_fei_method().
|
pure virtual |
Put element-dof solution guesses for an element-block.
elemBlockID | Element-block identifier. |
numElems | Length of the elemIDs list. |
elemIDs | Those elements for which elem-dof initial guesses are being supplied. |
dofPerElem | Number of degrees-of-freedom per element. |
estimates | List of length numElems*dofPerElem. The estimates for element i should start in estimates[i*dofPerElem]. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Put Lagrange Multiplier solution guesses.
numMultCRs | Length of the CRMultIDs and multEstimates lists. |
CRMultIDs | Identifiers obtained from earlier calls to 'initCRMult'. |
multEstimates | Initial guesses for the Lagrange Multipliers. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Return list of nodes associated with an element-block.
elemBlockID | Element-block identifier. |
numNodes | Allocated length of the nodeIDs list. If this is less than the value obtained via 'getNumBlockActNodes', then the first numNodes node id's will be placed in the nodeIDs list, ordered in ascending value. |
nodeIDs | Output. All nodes associated with this element-block. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Return llist of all elements associated with an element-block.
elemBlockID | Element-block identifier. |
numElems | Length of the elemIDs list. |
elemIDs | Output. All elements associated with this element-block. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::save_block_elem_soln().
|
pure virtual |
Return a version string. This string is owned by the FEI implementation, the calling application should not delete/free it. This string will contain the FEI implementation's version number, and if possible, a build time/date.
versionString | Output reference to a char*. The C interface will have a char** here. This function is simply setting versionString to point to the internal version string. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by feiDriver_main(), poisson_main(), and FEI_tester::testInitialization().
|
pure virtual |
Return cumulative cpu-time spent in each of 4 major FEI phases.
initPhase | Time in seconds, spent in the initialization phase. |
loadPhase | Time in seconds, spent loading data, up until the solver was launched. |
solve | Time in seconds, spent in the call to the underlying solver. |
solnReturn | Time in seconds, spent in the solution-return functions. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by poisson_main().
|
pure virtual |
Return the number of scalar degrees of freedom associated with a node.
globalNodeID | Globally unique node identifier |
numSolnParams | Sum of the sizes of the solution fields associated with the node. This will be the union of the set of fields defined on this node over all element-blocks it is associated with. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Return the number of element-blocks.
numElemBlocks |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Return the number of active nodes in an element-block.
elemBlockID | Input. |
numNodes | Output. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::testInitialization().
|
pure virtual |
Return the number of active equations in an element block.
elemBlockID | Input. |
numEqns | Output. Includes both nodal equations and elem-dof equations. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Return the number of nodes associated with elements of an element-block.
elemBlockID | Input. |
nodesPerElem | Output. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Return the number of equations at each element in an element-block. Includes elem-dof equations.
elemBlockID | Input. |
eqnsPerElem | Output. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
|
pure virtual |
Return the number of elements in an element-block.
blockID | Input. |
numElems | Output. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::save_block_elem_soln().
|
pure virtual |
Return the number of element-dof at elements in this element-block.
blockID | Input. |
DOFPerElem | Output. |
Implemented in FEI_Implementation, and fei::FEI_Impl.
Referenced by FEI_tester::save_block_elem_soln().