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

#include <fei_FEI_Impl.hpp>

Inheritance diagram for fei::FEI_Impl:
Inheritance graph
[legend]

Public Member Functions

 FEI_Impl (fei::SharedPtr< LibraryWrapper > wrapper, MPI_Comm comm, int masterRank=0)
 
 FEI_Impl (const fei::Factory *factory, MPI_Comm comm, int masterRank=0)
 
virtual ~FEI_Impl ()
 
fei::SharedPtr< fei::LinearSystemgetLinearSystem ()
 
int parameters (int numParams, const char *const *paramStrings)
 
int setIDLists (int numMatrices, const int *matrixIDs, int numRHSs, const int *rhsIDs)
 
int setSolveType (int solveType)
 
int initFields (int numFields, const int *fieldSizes, const int *fieldIDs, const int *fieldTypes=NULL)
 
int initElemBlock (GlobalID elemBlockID, int numElements, int numNodesPerElement, const int *numFieldsPerNode, const int *const *nodalFieldIDs, int numElemDofFieldsPerElement, const int *elemDOFFieldIDs, int interleaveStrategy)
 
int initElem (GlobalID elemBlockID, GlobalID elemID, const GlobalID *elemConn)
 
int initSlaveVariable (GlobalID slaveNodeID, int slaveFieldID, int offsetIntoSlaveField, int numMasterNodes, const GlobalID *masterNodeIDs, const int *masterFieldIDs, const double *weights, double rhsValue)
 
int deleteMultCRs ()
 
int initSharedNodes (int numSharedNodes, const GlobalID *sharedNodeIDs, const int *numProcsPerNode, const int *const *sharingProcIDs)
 
int initCRMult (int numCRNodes, const GlobalID *CRNodeIDs, const int *CRFieldIDs, int &CRID)
 
int initCRPen (int numCRNodes, const GlobalID *CRNodeIDs, const int *CRFieldIDs, int &CRID)
 
int initComplete ()
 
int setCurrentMatrix (int matrixID)
 
int setCurrentRHS (int rhsID)
 
int resetSystem (double s=0.0)
 
int resetMatrix (double s=0.0)
 
int resetRHSVector (double s=0.0)
 
int resetInitialGuess (double s=0.0)
 
int loadNodeBCs (int numNodes, const GlobalID *nodeIDs, int fieldID, const int *offsetsIntoField, const double *prescribedValues)
 
int loadElemBCs (int numElems, const GlobalID *elemIDs, int fieldID, const double *const *alpha, const double *const *beta, const double *const *gamma)
 
int sumInElem (GlobalID elemBlockID, GlobalID elemID, const GlobalID *elemConn, const double *const *elemStiffness, const double *elemLoad, int elemFormat)
 
int sumInElemMatrix (GlobalID elemBlockID, GlobalID elemID, const GlobalID *elemConn, const double *const *elemStiffness, int elemFormat)
 
int sumInElemRHS (GlobalID elemBlockID, GlobalID elemID, const GlobalID *elemConn, const double *elemLoad)
 
int loadCRMult (int CRMultID, int numCRNodes, const GlobalID *CRNodeIDs, const int *CRFieldIDs, const double *CRWeights, double CRValue)
 
int loadCRPen (int CRPenID, int numCRNodes, const GlobalID *CRNodeIDs, const int *CRFieldIDs, const double *CRWeights, double CRValue, double penValue)
 
int putIntoRHS (int IDType, int fieldID, int numIDs, const GlobalID *IDs, const double *coefficients)
 
int sumIntoRHS (int IDType, int fieldID, int numIDs, const GlobalID *IDs, const double *coefficients)
 
int setMatScalars (int numScalars, const int *IDs, const double *scalars)
 
int setRHSScalars (int numScalars, const int *IDs, const double *scalars)
 
int loadComplete (bool applyBCs=true, bool globalAssemble=true)
 
int residualNorm (int whichNorm, int numFields, int *fieldIDs, double *norms)
 
int solve (int &status)
 
int iterations (int &itersTaken) const
 
int version (const char *&versionString)
 
int cumulative_cpu_times (double &initTime, double &loadTime, double &solveTime, double &solnReturnTime)
 
int getBlockNodeSolution (GlobalID elemBlockID, int numNodes, const GlobalID *nodeIDs, int *offsets, double *results)
 
int getNodalSolution (int numNodes, const GlobalID *nodeIDs, int *offsets, double *results)
 
int getBlockFieldNodeSolution (GlobalID elemBlockID, int fieldID, int numNodes, const GlobalID *nodeIDs, double *results)
 
int getBlockElemSolution (GlobalID elemBlockID, int numElems, const GlobalID *elemIDs, int &numElemDOFPerElement, double *results)
 
int getNumCRMultipliers (int &numMultCRs)
 
int getCRMultIDList (int numMultCRs, int *multIDs)
 
int getCRMultipliers (int numCRs, const int *CRIDs, double *multipliers)
 
int putBlockNodeSolution (GlobalID elemBlockID, int numNodes, const GlobalID *nodeIDs, const int *offsets, const double *estimates)
 
int putBlockFieldNodeSolution (GlobalID elemBlockID, int fieldID, int numNodes, const GlobalID *nodeIDs, const double *estimates)
 
int putBlockElemSolution (GlobalID elemBlockID, int numElems, const GlobalID *elemIDs, int dofPerElem, const double *estimates)
 
int putCRMultipliers (int numMultCRs, const int *CRIDs, const double *multEstimates)
 
int getBlockNodeIDList (GlobalID elemBlockID, int numNodes, GlobalID *nodeIDs)
 
int getBlockElemIDList (GlobalID elemBlockID, int numElems, GlobalID *elemIDs)
 
int getNumSolnParams (GlobalID nodeID, int &numSolnParams) const
 
int getNumElemBlocks (int &numElemBlocks) const
 
int getNumBlockActNodes (GlobalID blockID, int &numNodes) const
 
int getNumBlockActEqns (GlobalID blockID, int &numEqns) const
 
int getNumNodesPerElement (GlobalID blockID, int &nodesPerElem) const
 
int getNumEqnsPerElement (GlobalID blockID, int &numEqns) const
 
int getNumBlockElements (GlobalID blockID, int &numElems) const
 
int getNumBlockElemDOF (GlobalID blockID, int &DOFPerElem) const
 
int getParameters (int &numParams, char **&paramStrings)
 
int getFieldSize (int fieldID, int &numScalars)
 
int getEqnNumbers (GlobalID ID, int idType, int fieldID, int &numEqns, int *eqnNumbers)
 
int getNodalFieldSolution (int fieldID, int numNodes, const GlobalID *nodeIDs, double *results)
 
int getNumLocalNodes (int &numNodes)
 
int getLocalNodeIDList (int &numNodes, GlobalID *nodeIDs, int lenNodeIDs)
 
int putNodalFieldData (int fieldID, int numNodes, const GlobalID *nodeIDs, const double *nodeData)
 
- Public Member Functions inherited from FEI
virtual ~FEI ()
 

Detailed Description

Implementation of the original abstract FEI interface. This FEI_Impl class plays the same role as the FEI_Implementation class, except that this class internally is based on the newer more modular fei:: classes. The result should be the same functionality as FEI_Implementation but with improved performance.

The ultimate goal is for this class to replace the FEI_Implementation class, thereby greatly reducing the total body of code that constitutes the fei code-distribution.

Definition at line 33 of file fei_FEI_Impl.hpp.

Constructor & Destructor Documentation

fei::FEI_Impl::FEI_Impl ( fei::SharedPtr< LibraryWrapper >  wrapper,
MPI_Comm  comm,
int  masterRank = 0 
)

Constructor

Definition at line 32 of file fei_FEI_Impl.cpp.

fei::FEI_Impl::FEI_Impl ( const fei::Factory factory,
MPI_Comm  comm,
int  masterRank = 0 
)

Constructor

Definition at line 90 of file fei_FEI_Impl.cpp.

fei::FEI_Impl::~FEI_Impl ( )
virtual

Destructor

Definition at line 154 of file fei_FEI_Impl.cpp.

Member Function Documentation

fei::SharedPtr< fei::LinearSystem > fei::FEI_Impl::getLinearSystem ( )

This function isn't in the original FEI 2.1 interface. But it is useful to users who may be mixing old and new fei objects and need access to underlying stuff.

Definition at line 247 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::parameters ( int  numParams,
const char *const *  paramStrings 
)
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.

Parameters
numParamsNumber of parameters being supplied.
paramStringsList of 'numParams' strings. Each string usually contains a key-value pair, separated by a space.

Implements FEI.

Definition at line 252 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::setIDLists ( int  numMatrices,
const int *  matrixIDs,
int  numRHSs,
const int *  rhsIDs 
)
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).

Parameters
numMatriceslength of matrixIDs parameter
matrixIDslist of user-defined identifiers for separate matrices to be assembled
numRHSslength of rhsIDs parameter
rhsIDslist of user-defined identifiers for separate rhs vectors to be assembled

Implements FEI.

Definition at line 206 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::setSolveType ( int  solveType)
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).

Parameters
solveTypecurrently supported values for this are: FEI_SINGLE_SOLVE, FEI_AGGREGATE_SUM

Implements FEI.

Definition at line 277 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::initFields ( int  numFields,
const int *  fieldSizes,
const int *  fieldIDs,
const int *  fieldTypes = NULL 
)
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.

Parameters
numFieldsGlobal number of fields in the entire problem, on all processors. (This is the length of the fieldSizes and fieldIDs lists.)
fieldSizesNumber of scalars contained in each field.
fieldIDsUser-supplied identifiers for each field.

Implements FEI.

Definition at line 285 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::initElemBlock ( GlobalID  elemBlockID,
int  numElements,
int  numNodesPerElement,
const int *  numFieldsPerNode,
const int *const *  nodalFieldIDs,
int  numElemDofFieldsPerElement,
const int *  elemDOFFieldIDs,
int  interleaveStrategy 
)
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.

Parameters
elemBlockIDThe user-defined identifier for this element-block.
numElementsThe number of elements in this block.
numNodesPerElementLength of the numFieldsPerNode list.
numFieldsPerNodeLengths of the rows of the nodalFieldIDs table.
nodalFieldIDsTable where row 'i' is the list of field ids for the ith node on every element in this element-block.
numElemDofFieldsPerElementLength of the elemDOFFieldIDs list.
elemDOFFieldIDslist of field identifiers for the element- centered degrees-of-freedom in the elements in this block.
interleaveStrategyIndicates 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, ...)

Implements FEI.

Definition at line 295 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::initElem ( GlobalID  elemBlockID,
GlobalID  elemID,
const GlobalID *  elemConn 
)
virtual

Initialize an element's connectivity. Provide a nodal connectivity list for inclusion in the sparse matrix structure being constructed.

Parameters
elemBlockIDWhich element-block this element belongs to.
elemIDA user-provided identifier for this element.
elemConnList of nodeIDs connected to this element. Length of this list must be 'numNodesPerElement' provided to the function 'initElemBlock' for this elemBlockID.

Implements FEI.

Definition at line 346 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::initSlaveVariable ( GlobalID  slaveNodeID,
int  slaveFieldID,
int  offsetIntoSlaveField,
int  numMasterNodes,
const GlobalID *  masterNodeIDs,
const int *  masterFieldIDs,
const double *  weights,
double  rhsValue 
)
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.

Parameters
slaveNodeIDNode identifier of the node containing the slave eqn.
slaveFieldIDField identifier corresponding to the slave eqn.
offsetIntoSlaveFieldDenotes location, within the field, of the slave eqn.
numMasterNodesNumber of nodes containing variables on which the slave depends.
masterNodeIDsNode identifiers of the master nodes.
masterFieldIDsList, length numMasterNodes, of the field at each master-node which contains the scalar variable(s) on which the slave depends.
weightsList, length sum-of-master-field-sizes, containing the weighting coefficients described above.
rhsValue

Implements FEI.

Definition at line 378 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::deleteMultCRs ( )
virtual

Request that any existing Lagrange-Multiplier constraints be deleted. (Intended to be called in preparation for loading new/different constraints.)

Implements FEI.

Definition at line 391 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::initSharedNodes ( int  numSharedNodes,
const GlobalID *  sharedNodeIDs,
const int *  numProcsPerNode,
const int *const *  sharingProcIDs 
)
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.

Parameters
numSharedNodesLength of the sharedNodeIDs and numProcsPerNode lists, and the number of rows in the sharingProcIDs table.
sharedNodeIDsList of nodes that are shared by this processor and one or more other processors.
numProcsPerNodeLengths of the rows of the sharingProcIDs table.
sharingProcIDsList, for each sharedNode, of processors which share that node.

Implements FEI.

Definition at line 397 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::initCRMult ( int  numCRNodes,
const GlobalID *  CRNodeIDs,
const int *  CRFieldIDs,
int &  CRID 
)
virtual

Constraint relation initialization, Lagrange Multiplier formulation.

Parameters
numCRNodesLength of the CRNodeIDs and CRFieldIDs lists.
CRNodeIDsNodes involved in this constraint relation.
CRFieldIDsList of the the field being constrained at each node.
CRIDOutput. 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.

Implements FEI.

Definition at line 411 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::initCRPen ( int  numCRNodes,
const GlobalID *  CRNodeIDs,
const int *  CRFieldIDs,
int &  CRID 
)
virtual

Constraint relation initialization, Penalty function formulation .

Parameters
numCRNodesLength of the CRNodeIDs and CRFieldIDs lists.
CRNodeIDsNodes involved in this constraint relation.
CRFieldIDsList of the the field being constrained at each node.
CRIDOutput. An identifier by which this constraint relation may be referred to later, when loading its weight data and penalty value.

Implements FEI.

Definition at line 430 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::initComplete ( )
virtual

indicate that overall initialization sequence is complete

Implements FEI.

Definition at line 449 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::setCurrentMatrix ( int  matrixID)
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.

Parameters
matrixIDOne of the identifiers provided earlier via the function 'setIDLists'.

Implements FEI.

Definition at line 560 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::setCurrentRHS ( int  rhsID)
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.

Parameters
rhsIDOne of the identifiers provided earlier via the function 'setIDLists'.

Implements FEI.

Definition at line 575 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::resetSystem ( double  s = 0.0)
virtual

Set a value (usually zero) througout the linear system.

Parameters
sThe value to be written into the linear system.

Implements FEI.

Definition at line 590 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::resetMatrix ( double  s = 0.0)
virtual

Set a value (usually zero) througout the matrix.

Parameters
sThe value to be written into the matrix.

Implements FEI.

Definition at line 599 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::resetRHSVector ( double  s = 0.0)
virtual

Set a value (usually zero) througout the rhs vector.

Parameters
sThe value to be written into the rhs vector.

Implements FEI.

Definition at line 604 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::resetInitialGuess ( double  s = 0.0)
virtual

Set a value (usually, if not always, 0.0) throughout the initial guess (solution) vector.

Parameters
sInput. Scalar value to use in filling the solution vector.
Returns
error-code 0 if successful

Implements FEI.

Definition at line 609 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::loadNodeBCs ( int  numNodes,
const GlobalID *  nodeIDs,
int  fieldID,
const int *  offsetsIntoField,
const double *  prescribedValues 
)
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.

Parameters
numNodesLength of the nodeIDs list.
nodeIDsList of nodes upon which a boundary condition is to be imposed.
fieldIDThe solution field that will receive the boundary condition.
offsetsIntoFieldArray, length numNodes.
prescribedValuesArray, length numNodes.

Implements FEI.

Definition at line 614 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::loadElemBCs ( int  numElems,
const GlobalID *  elemIDs,
int  fieldID,
const double *const *  alpha,
const double *const *  beta,
const double *const *  gamma 
)
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.

Parameters
numElemsLength of the elemIDs list.
elemIDsList of elements for which a boundary condition is to be specified.
fieldIDThe solution field for which to apply the boundary condition.
alphaTable, as in 'loadNodeBCs', but with 'numElems' number-of- rows.
betaTable, same dimensions as alpha.
gammaTable, same dimensions as alpha.

Implements FEI.

Definition at line 629 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::sumInElem ( GlobalID  elemBlockID,
GlobalID  elemID,
const GlobalID *  elemConn,
const double *const *  elemStiffness,
const double *  elemLoad,
int  elemFormat 
)
virtual

Element-stiffness/load data loading. This function accumulates element stiffness and load data into the underlying matrix and rhs vector.

Parameters
elemBlockIDWhich element-block this element belongs to.
elemIDUser-supplied identifier for this element.
elemConnConnectivity list of nodes that are connected to this element.
elemStiffnessTable 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'.)
elemLoadElement-load vector.
elemFormatDesignates 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.

Implements FEI.

Definition at line 640 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::sumInElemMatrix ( GlobalID  elemBlockID,
GlobalID  elemID,
const GlobalID *  elemConn,
const double *const *  elemStiffness,
int  elemFormat 
)
virtual

Element-stiffness data loading. This function is the same as 'sumInElem' but only accepts stiffness data, not the load data for the rhs.

Parameters
elemBlockIDWhich element-block this element belongs to.
elemIDUser-supplied identifier for this element.
elemConnConnectivity list of nodes that are connected to this element.
elemStiffnessTable 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'.)
elemFormatDesignates 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.

Implements FEI.

Definition at line 660 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::sumInElemRHS ( GlobalID  elemBlockID,
GlobalID  elemID,
const GlobalID *  elemConn,
const double *  elemLoad 
)
virtual

Element-load data loading. This function is the same as 'sumInElem', but only accepts the load for the rhs, not the stiffness matrix.

Parameters
elemBlockIDWhich element-block this element belongs to.
elemIDUser-supplied identifier for this element.
elemConnConnectivity list of nodes that are connected to this element.
elemLoadElement-load vector.

Implements FEI.

Definition at line 673 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::loadCRMult ( int  CRMultID,
int  numCRNodes,
const GlobalID *  CRNodeIDs,
const int *  CRFieldIDs,
const double *  CRWeights,
double  CRValue 
)
virtual

Load weight/value data for a Lagrange Multiplier constraint relation.

Parameters
CRMultIDIdentifier returned from an earlier call to 'initCRMult'.
numCRNodesLength of CRNodeIDs and CRFieldIDs lists.
CRNodeIDsList of nodes in this constraint relation.
CRFieldIDsList of fields, one per node, to be constrained.
CRWeightsWeighting coefficients. This length of this list is the sum of the sizes associated with the fields identified in CRFieldIDs.
CRValueThe constraint's rhs value. Often (always?) zero.

Implements FEI.

Definition at line 689 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::loadCRPen ( int  CRPenID,
int  numCRNodes,
const GlobalID *  CRNodeIDs,
const int *  CRFieldIDs,
const double *  CRWeights,
double  CRValue,
double  penValue 
)
virtual

Load weight/value data for a Penalty constraint relation.

Parameters
CRPenIDIdentifier returned from an earlier call to 'initCRPen'.
numCRNodesLength of CRNodeIDs and CRFieldIDs lists.
CRNodeIDsList of nodes in this constraint relation.
CRFieldIDsList of fields, one per node, to be constrained.
CRWeightsWeighting coefficients. This length of this list is the sum of the sizes associated with the fields identified in CRFieldIDs.
CRValueThe constraint's rhs value. Often (always?) zero.
penValueThe penalty value.

Implements FEI.

Definition at line 703 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::putIntoRHS ( int  IDType,
int  fieldID,
int  numIDs,
const GlobalID *  IDs,
const double *  coefficients 
)
virtual

Put a copy of coefficient data into the rhs vector.

Implements FEI.

Definition at line 718 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::sumIntoRHS ( int  IDType,
int  fieldID,
int  numIDs,
const GlobalID *  IDs,
const double *  coefficients 
)
virtual

Sum a copy of coefficient data into the rhs vector.

Implements FEI.

Definition at line 731 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::setMatScalars ( int  numScalars,
const int *  IDs,
const double *  scalars 
)
virtual

Set scalars by which to multiply matrices, in cases where a linear- combination of several matrices is to be solved.

Parameters
numScalarsLength of the IDs and scalars lists.
IDsMatrix identifiers which must be a subset of those supplied via an earlier call to 'setIDLists'.
scalarsThe coefficients by which to multiply the matrices.

Implements FEI.

Definition at line 772 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::setRHSScalars ( int  numScalars,
const int *  IDs,
const double *  scalars 
)
virtual

Set scalars by which to multiply rhs vectors, in cases where a linear- combination of several rhs vectors is to be solved.

Parameters
numScalarsLength of the IDs and scalars lists.
IDsRHS-vector identifiers which must be a subset of those supplied via an earlier call to 'setIDLists'.
scalarsThe coefficients by which to multiply the rhs vectors.

Implements FEI.

Definition at line 794 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::loadComplete ( bool  applyBCs = true,
bool  globalAssemble = true 
)
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.

Implements FEI.

Definition at line 816 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::residualNorm ( int  whichNorm,
int  numFields,
int *  fieldIDs,
double *  norms 
)
virtual

get residual norms

Implements FEI.

Definition at line 920 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::solve ( int &  status)
virtual

launch underlying solver

Implements FEI.

Definition at line 1042 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::iterations ( int &  itersTaken) const
virtual

Query number of iterations taken for last solve.

Parameters
itersTakenIterations performed during any previous solve.

Implements FEI.

Definition at line 1062 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::version ( const char *&  versionString)
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.

Parameters
versionStringOutput 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.

Implements FEI.

Definition at line 1069 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::cumulative_cpu_times ( double &  initTime,
double &  loadTime,
double &  solveTime,
double &  solnReturnTime 
)
virtual

query for some accumulated timing information. Collective function.

Implements FEI.

Definition at line 1076 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getBlockNodeSolution ( GlobalID  elemBlockID,
int  numNodes,
const GlobalID *  nodeIDs,
int *  offsets,
double *  results 
)
virtual

return all nodal solution params on a block-by-block basis

Implements FEI.

Definition at line 1089 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNodalSolution ( int  numNodes,
const GlobalID *  nodeIDs,
int *  offsets,
double *  results 
)
virtual

return all nodal solution params for an arbitrary list of nodes

Implements FEI.

Definition at line 1099 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getBlockFieldNodeSolution ( GlobalID  elemBlockID,
int  fieldID,
int  numNodes,
const GlobalID *  nodeIDs,
double *  results 
)
virtual

return nodal solution for one field on a block-by-block basis

Implements FEI.

Definition at line 1145 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getBlockElemSolution ( GlobalID  elemBlockID,
int  numElems,
const GlobalID *  elemIDs,
int &  numElemDOFPerElement,
double *  results 
)
virtual

return element solution params on a block-by-block basis

Implements FEI.

Definition at line 1155 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNumCRMultipliers ( int &  numMultCRs)
virtual

Query number of lagrange-multiplier constraints on local processor

Implements FEI.

Definition at line 1215 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getCRMultIDList ( int  numMultCRs,
int *  multIDs 
)
virtual

Obtain list of lagrange-multiplier IDs

Implements FEI.

Definition at line 1221 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getCRMultipliers ( int  numCRs,
const int *  CRIDs,
double *  multipliers 
)
virtual

get Lagrange Multipliers

Implements FEI.

Definition at line 1229 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::putBlockNodeSolution ( GlobalID  elemBlockID,
int  numNodes,
const GlobalID *  nodeIDs,
const int *  offsets,
const double *  estimates 
)
virtual

put nodal-based solution guess on a block-by-block basis

Implements FEI.

Definition at line 1244 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::putBlockFieldNodeSolution ( GlobalID  elemBlockID,
int  fieldID,
int  numNodes,
const GlobalID *  nodeIDs,
const double *  estimates 
)
virtual

put nodal-based guess for one field on a block-by-block basis

Implements FEI.

Definition at line 1254 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::putBlockElemSolution ( GlobalID  elemBlockID,
int  numElems,
const GlobalID *  elemIDs,
int  dofPerElem,
const double *  estimates 
)
virtual

put element-based solution guess on a block-by-block basis

Implements FEI.

Definition at line 1264 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::putCRMultipliers ( int  numMultCRs,
const int *  CRIDs,
const double *  multEstimates 
)
virtual

put Lagrange solution to FE analysis on a constraint-set basis

Implements FEI.

Definition at line 1274 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getBlockNodeIDList ( GlobalID  elemBlockID,
int  numNodes,
GlobalID *  nodeIDs 
)
virtual

return info associated with blocked nodal solution

Implements FEI.

Definition at line 1282 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getBlockElemIDList ( GlobalID  elemBlockID,
int  numElems,
GlobalID *  elemIDs 
)
virtual

return info associated with blocked element solution

Implements FEI.

Definition at line 1325 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNumSolnParams ( GlobalID  nodeID,
int &  numSolnParams 
) const
virtual

Query number of degrees-of-freedom for specified node

Implements FEI.

Definition at line 1344 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNumElemBlocks ( int &  numElemBlocks) const
virtual

Query number of element-blocks on local processor

Implements FEI.

Definition at line 1350 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNumBlockActNodes ( GlobalID  blockID,
int &  numNodes 
) const
virtual

return the number of active nodes in a given element block

Implements FEI.

Definition at line 1356 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNumBlockActEqns ( GlobalID  blockID,
int &  numEqns 
) const
virtual

return the number of active equations in a given element block

Implements FEI.

Definition at line 1366 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNumNodesPerElement ( GlobalID  blockID,
int &  nodesPerElem 
) const
virtual

return the number of nodes associated with elements of a given block ID

Implements FEI.

Definition at line 1372 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNumEqnsPerElement ( GlobalID  blockID,
int &  numEqns 
) const
virtual

return the number of equations (including element eqns) associated with elements of a given block ID

Implements FEI.

Definition at line 1378 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNumBlockElements ( GlobalID  blockID,
int &  numElems 
) const
virtual

return the number of elements associated with this blockID

Implements FEI.

Definition at line 1384 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNumBlockElemDOF ( GlobalID  blockID,
int &  DOFPerElem 
) const
virtual

return the number of elements eqns for elems w/ this blockID

Implements FEI.

Definition at line 1397 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getParameters ( int &  numParams,
char **&  paramStrings 
)

return the parameters that have been set so far. The caller should NOT delete the paramStrings pointer.

Definition at line 1406 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getFieldSize ( int  fieldID,
int &  numScalars 
)
virtual

Query the size of a field. This info is supplied to the FEI (initFields) by the application, but may not be easily obtainable on the app side at all times. Thus, it would be nice if the FEI could answer this query.

Implements FEI.

Definition at line 1413 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getEqnNumbers ( GlobalID  ID,
int  idType,
int  fieldID,
int &  numEqns,
int *  eqnNumbers 
)
virtual

Since the ultimate intent for matrix-access is to bypass the FEI and go straight to 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.

Parameters
IDIdentifier of either a node or an element.
idTypeCan take either of the values FEI_NODE or FEI_ELEMENT.
fieldIDIdentifies a particular field at this [node||element].
numEqnsOutput. Number of equations associated with this node/field (or element/field) pair.
eqnNumbersCaller-allocated array. On exit, this is filled with the above-described equation-numbers. They are global 0-based numbers.

Implements FEI.

Definition at line 1419 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNodalFieldSolution ( int  fieldID,
int  numNodes,
const GlobalID *  nodeIDs,
double *  results 
)
virtual

Get the solution data for a particular field, on an arbitrary set of nodes.

Parameters
fieldIDInput. field identifier for which solution data is being requested.
numNodesInput. Length of the nodeIDs list.
nodeIDsInput. List specifying the nodes on which solution data is being requested.
resultsAllocated 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'.
Returns
error-code 0 if successful

Implements FEI.

Definition at line 1430 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getNumLocalNodes ( int &  numNodes)
virtual

Get the number of nodes that are local to this processor (includes nodes that are shared by other processors).

Parameters
numNodesOutput. Number of local nodes.
Returns
error-code 0 if successful

Implements FEI.

Definition at line 1440 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::getLocalNodeIDList ( int &  numNodes,
GlobalID *  nodeIDs,
int  lenNodeIDs 
)
virtual

Get a list of the nodeIDs that are local to this processor (includes nodes that are shared by other processors).

Parameters
numNodesOutput. Same as the value output by 'getNumLocalNodes'.
nodeIDsCaller-allocated array, contents to be filled by this function.
lenNodeIDsInput. 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.
Returns
error-code 0 if successful

Implements FEI.

Definition at line 1446 of file fei_FEI_Impl.cpp.

int fei::FEI_Impl::putNodalFieldData ( int  fieldID,
int  numNodes,
const GlobalID *  nodeIDs,
const double *  nodeData 
)
virtual
Pass nodal data for a specified field through to the solver. Example

is geometric coordinates, etc.

Parameters
fieldIDfield identifier for the data to be passed. This is probably not one of the 'solution-field' identifiers. It should be an identifier that the solver is expecting and knows how to handle. This field id must previously have been provided to the fei implementation via the normal initFields method.
numNodesnumber of nodes for which data is being provided
nodeIDsList of length numNodes, giving node-identifiers for which data is being provided.
nodeDataList of length fieldSize*numNodes, where fieldSize is the size which was associated with fieldID when initFields was called. The data for nodeIDs[i] begins in position i*fieldSize of this array.

Implements FEI.

Definition at line 1454 of file fei_FEI_Impl.cpp.


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