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

#include <fei_Lookup.hpp>

Inheritance diagram for Lookup:
Inheritance graph
[legend]

Public Member Functions

virtual ~Lookup ()
 
virtual int getNumFields ()=0
 
virtual int getFieldSize (int fieldID)=0
 
virtual const int * getFieldIDsPtr ()=0
 
virtual const int * getFieldSizesPtr ()=0
 
virtual int getNumElemBlocks ()=0
 
virtual const GlobalID * getElemBlockIDs ()=0
 
virtual void getElemBlockInfo (GlobalID blockID, int &interleaveStrategy, int &lumpingStrategy, int &numElemDOF, int &numElements, int &numNodesPerElem, int &numEqnsPerElem)=0
 
virtual const int * getNumFieldsPerNode (GlobalID blockID)=0
 
virtual const int *const * getFieldIDsTable (GlobalID blockID)=0
 
virtual int getEqnNumber (int nodeNumber, int fieldID)=0
 
virtual int getAssociatedNodeNumber (int eqnNumber)=0
 
virtual int getAssociatedFieldID (int eqnNumber)=0
 
virtual bool isInLocalElement (int nodeNumber)=0
 
virtual int getNumSubdomains (int nodeNumber)=0
 
virtual int * getSubdomainList (int nodeNumber)=0
 
virtual int getNumSharedNodes ()=0
 
virtual const int * getSharedNodeNumbers ()=0
 
virtual const int * getSharedNodeProcs (int nodeNumber)=0
 
virtual int getNumSharingProcs (int nodeNumber)=0
 
virtual bool isExactlyBlkEqn (int ptEqn)=0
 
virtual int ptEqnToBlkEqn (int ptEqn)=0
 
virtual int getOffsetIntoBlkEqn (int blkEqn, int ptEqn)=0
 
virtual int getBlkEqnSize (int blkEqn)=0
 

Detailed Description

This interface is intended to be used by a LinearSystemCore implementation or a FiniteElementData implementation to look up various information about the structure of the finite-element problem being assembled into LinearSystemCore by a FEI implementation.

Background: The finite element problem consists of a set of element-blocks, each of which contains a set of elements. Each element has a list of connected nodes, and each node has a set of fields. Each field consists of 1 to several scalar quantities. Each of those scalar quantities corresponds to an equation in the linear system. Exception: some fields may not be solution-fields. This is indicated by a negative fieldID. There are no equations corresponding to fields with negative fieldIDs. Data that is passed in associated with negative fieldIDs is probably coordinate or nullspace data, or other data passed from the application to the solver (note that this non-solution-field data won't appear in element-matrices, it will be passed separately through a special function for supplying nodal data).

elem-block IDs and field IDs are application-provided numbers, and no assumption may be made regarding their order, contiguousness, etc.

Equation numbers are assigned to node/field pairs by the FEI implementation, and are also globally unique, and non-shared. Each equation resides on exactly one processor. Equation numbers are 0-based.

NOTES: 1. functions that return an equation number, or a size (e.g., num-equations, num-fields, etc.) may indicate an error or 'not-found' condition by returning a negative number. Functions that return a pointer may indicate an error by returning NULL.

Definition at line 40 of file fei_Lookup.hpp.

Constructor & Destructor Documentation

virtual Lookup::~Lookup ( )
inlinevirtual

Destructor.

Definition at line 43 of file fei_Lookup.hpp.

Member Function Documentation

virtual int Lookup::getNumFields ( )
pure virtual

Get the (global) number of fields defined for the problem. A field may consist of 1 to several scalar quantities. Examples include pressure, displacement (a 3-vector in 3D), etc.

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::getFieldSize ( int  fieldID)
pure virtual

Given a fieldID, obtain the associated fieldSize (number of scalar components)

Parameters
fieldIDinteger identifier for a field

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual const int* Lookup::getFieldIDsPtr ( )
pure virtual

Return a pointer to a list (of length numFields) of the fieldIDs

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual const int* Lookup::getFieldSizesPtr ( )
pure virtual

Return a pointer to a list (of length numFields) of the fieldSizes

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::getNumElemBlocks ( )
pure virtual

Return the number of element-blocks in the (local) finite-element problem.

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual const GlobalID* Lookup::getElemBlockIDs ( )
pure virtual

Return a pointer to the list (of length numElemBlocks) containing the element-block identifiers for the (local) finite-element problem.

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual void Lookup::getElemBlockInfo ( GlobalID  blockID,
int &  interleaveStrategy,
int &  lumpingStrategy,
int &  numElemDOF,
int &  numElements,
int &  numNodesPerElem,
int &  numEqnsPerElem 
)
pure virtual

Given a blockID, provide several pieces of element-block information.

Parameters
interleaveStrategyelement-equation ordering: 0 => node-major, 1 => field-major
lumpingStrategyelement-matrices may be lumped if they're mass matrices, 0 => not lumped, 1 => lumped
numElemDOFnumber of element-dof at each element in this block
numElementsnumber of elements in this block
numNodesPerElemnumber of nodes connected to each element in this block
numEqnsPerElemnumber of scalar equations at each element in this block

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual const int* Lookup::getNumFieldsPerNode ( GlobalID  blockID)
pure virtual

Given a blockID, return a pointer to a list (of length numNodesPerElem) of numFieldsPerNode.

Parameters
blockIDidentifier of the elem-block in question

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual const int* const* Lookup::getFieldIDsTable ( GlobalID  blockID)
pure virtual

Given a blockID, return a pointer to a table, (num-rows == numNodesPerElem, row-length[i] == fieldsPerNode[i]) containing the fieldIDs at each node of elements in that element-block.

Parameters
blockIDidentifier of the elem-block in question

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::getEqnNumber ( int  nodeNumber,
int  fieldID 
)
pure virtual

Given a nodeNumber/fieldID pair, this function returns the first global (0-based) equation number associated with that nodeNumber/fieldID pair.

Parameters
nodeNumber
fieldID

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::getAssociatedNodeNumber ( int  eqnNumber)
pure virtual

Given a global (0-based) equation number, return the node-number.

Parameters
eqnNumber

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::getAssociatedFieldID ( int  eqnNumber)
pure virtual

Given a global (0-based) equation number, return the fieldID.

Parameters
eqnNumber

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual bool Lookup::isInLocalElement ( int  nodeNumber)
pure virtual

Given a nodeNumber, determine whether that node is connected to a local element.

Parameters
nodeNumber
Returns
true if nodeNumber is connected to a local element, false otherwise.

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::getNumSubdomains ( int  nodeNumber)
pure virtual

Given a nodeNumber, return the number of subdomains that contain this node. subdomains correspond to processors. The number of subdomains that contain a node does not always equal the number of processors that share a node. There are two kinds of "sharing" – the "normal" kind, where a node is shared because it is connected to elements that reside on more than one processor, and the "wierd" kind where nodes are considered shared simply because of being in cross-processor constraints. This function describes how many processors share this node in the "normal" sense. Thus, in general, this relationship holds: getNumSubdomains(nodeNum) <= getNumSharingProcs(nodeNum)

Parameters
nodeNumber
Returns
numSubdomains

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int* Lookup::getSubdomainList ( int  nodeNumber)
pure virtual

Given a nodeNumber, return a list of the subdomains that contain this node. This is the list of subdomains that's counted by the method 'getNumSubdomains' above.

Parameters
nodeNumber
Returns
pointer to list of subdomains (processor ranks). NULL if the specified nodeNumber is not found.

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::getNumSharedNodes ( )
pure virtual

Return the number of local nodes that are shared by multiple processors

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual const int* Lookup::getSharedNodeNumbers ( )
pure virtual

Return a pointer to the list of shared nodeNumbers

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual const int* Lookup::getSharedNodeProcs ( int  nodeNumber)
pure virtual

Given a shared nodeNumber, return a pointer to the list of sharing procs.

Parameters
nodeNumberThe subject of the query. Function returns NULL if 'nodeNumber' is not a shared node.

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::getNumSharingProcs ( int  nodeNumber)
pure virtual

Given a shared nodeNumber, return the number of processors that share it.

Parameters
nodeNumberFunction returns -1 if 'nodeNumber' is not a shared node.

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual bool Lookup::isExactlyBlkEqn ( int  ptEqn)
pure virtual

Query whether a pt-eqn corresponds exactly to a blk-eqn. in other words, is pt-eqn the first point equation in a block-equation.

Parameters
ptEqn

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::ptEqnToBlkEqn ( int  ptEqn)
pure virtual

Given a pt-eqn, return the corresponding blk-eqn.

Parameters
ptEqn

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::getOffsetIntoBlkEqn ( int  blkEqn,
int  ptEqn 
)
pure virtual

Given a blk-eqn and a pt-eqn, return the pt-eqn's offset into the blk-eqn (i.e., distance from the 'beginning' of the blk-eqn)

Parameters
blkEqn
ptEqn

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.

virtual int Lookup::getBlkEqnSize ( int  blkEqn)
pure virtual

Given a blk-eqn, return the 'size', or number of pt-eqns corresponding to it.

Parameters
blkEqn

Implemented in SNL_FEI_Structure, and fei::Lookup_Impl.


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