FemProcess
The GeMA Fem Process Plugin
Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
GmpFemAssembler Class Reference

A class responsible for assembling the global stiffness matrix, handling multiple physics with possibly diferent degrees of freedom for each node. More...

#include <gmpFemAssembler.h>

Collaboration diagram for GmpFemAssembler:
Collaboration graph
[legend]

Classes

struct  TemplateData
 Struct storing node metadata information. More...
 

Public Types

enum  FixedDofMode { NO_DOF_HANDLING, REMOVE_DOF }
 
enum  GhostDofMode { GHOST_SUPPORT, NO_GHOST_SUPPORT, AUTO_GHOST_SUPPORT }
 
enum  UpdateMode {
  NODES_ADDED = 0x001, DOFS_ADDED = 0x002, DOFS_REMOVED = 0x004, ELEMENTS_ADDED = 0x008,
  ELEMENTS_REMOVED = 0x010, FIXED_BCS_UPDATE = 0x020, FIXED_BCS_CHANGED = 0x040, FIXED_BCS_ALL = 0x080,
  FIXED_BCS_KEEPOLD = 0x100, FULL_UPDATE = 0x4000, SAVE_MAPPING_INFO = 0x8000
}
 

Public Member Functions

 GmpFemAssembler (GmElementMesh *mesh, const QList< GmpFemPhysics * > &physics)
 Constructor. Receives as parameters the mesh and the list of physics objects that will cooperate to build the global stiffness matrix. After object creation, a call to init() must be made before the assembler can be used.
 
 ~GmpFemAssembler ()
 Destructor.
 
bool init (GmpFemMatrixSet *matSet, GmpFemVectorSet *vecSet, FixedDofMode fixedMode, GhostDofMode ghostMode, bool enableFastUpdate, bool createReverseMapping, int numThreads, QString &err)
 Initializes the assembler. Must be called after object construction. More...
 
bool update (int mode, QString &err)
 Updates the assembler structure after a change in the model free or fixed dofs, according to the given mode. Depending on the requested mode and on the enableFastUpdate flag passed to init(), this might result in an update of only the affected dofs or in a full update, roughly equivalent to destroying and rebuilding the assembler (except that recreating the assembler will not create remaping information). More...
 
void setPartitions (int npart)
 Informs the assembler about the number of partitions used when working with multiple threads. Must be called once prior to calling addElementData()
 
int numDof () const
 Returns the number of degrees of freedom of the global assembled matrix.
 
int numTotalDof () const
 Returns the total number of the system degrees of freedom, including fixed ones. This value will be equal to the value returned by numDof() if the mode parameter given in the call to init() was equal to NO_DOF_HANDLING.
 
int numFixedDof () const
 
int numElementDof (const GmElement *e, int physicsIndex) const
 Returns the number of degrees of freedom for the local matrix calculated by the referenced physics for the given element. More...
 
int maxElementDof () const
 Returns the maximum number of degrees of freedom used by a physics/element type.
 
int dofIndex (int node, int dof) const
 Given a global node index and a global degree of freedom index, returns its line/column in the assembled matrix or a negative value if the pair does not belong to it. More...
 
int fillIndexMap (const GmElement *e, int physicsIndex, int *indexMap, const GmElementDof *&dofMap)
 Fills indexMap with the translation map from local to global matrix/vector indices and returns the number of entries filled in the map (which will be equal to the local element matrix/vector sizes). Also makes dofMap point to the physics/element dof map. More...
 
double dofFixedValue (int index)
 Returns the fixed dof value associated to a negative index in an entry in the vector returned by fillIndexMap()
 
GmStateVarstateVar (int globalIndex)
 Given a global index in the force vector (global dof index), returns the associated state variable.
 
int svGroupId (int globalIndex)
 Given a global index in the force vector (global dof index), returns the associated state var group.
 
const QVector< QPair< GmStateVar *, GmValueAccessor * > > & stateVars () const
 Returns the vector storing the set of state variables/accesors to them in use by the assembler.
 
Unit timeUnit () const
 Returns the unit in which time is handled by the physics.
 
bool beginAssembly ()
 Notifies the assembler that the assembly process will begin. Also initializes global active vectors and matrices with zeros. More...
 
bool endAssembly (bool discardData=false)
 Notifies the assembler that the assembly process has ended. See aditional comments on beginAssembly(). When discardData is set to true, the assembly process has not been completed due to an error and the resulting intermediate data should be discarded.
 
const int * addElementData (const GmElement *e, int physicsIndex, const GmpFemAssemblerMatrixCombiner *matCombiner, const GmpFemAssemblerVectorCombiner *vecCombiner)
 Adds data from a set of local element matrix & vectors to their respective global data. Depending on the mode used for creating the matrix and the vector sets (the mode might be different), the assembled data will be added to matching global matrices/vectors and/or to the equivalent matrix/vector. Enabling both of them is usefull for debugging purposes. More...
 
bool addLocalMatrixData (GmMatrixDof &localDofMap, const GmpFemAssemblerMatrixCombiner *matCombiner, const GmpFemAssemblerVectorCombiner *vecCombiner)
 Simmilar to addElementData() but instead of the matrices in the matrix and vector set being tied to an element, they are tied to a local matrix with arbitrary size. The localDofMap parameter gives the matrix layout. It also provides a vector that should be filled with the local to global dof translation map. More...
 
void saveResultData (const GmVector &data, double lambda=1.0)
 Saves the contents of the data vector to the appropriate state vars / nodes according to the mapping between vector lines and the corresponding degree of freedom. The lambda factor will be applied only to fixed values handled by the assembler.
 
void saveResultData (const GmVector &data, int baseDof, GmValueAccessor *ac, double lambda=1.0)
 Saves the contents of the data vector to the given accessor. Only the subset in the source vector associated with the given dof (or dofs for vector based accessors) are saved to the accessor. The lambda factor will be applied only to fixed values handled by the assembler. More...
 
void saveResultDataToFixedDofs (const GmVector &data, int baseDof, GmValueAccessor *ac, bool clearFree)
 Saves the contents of the data vector to the fixed dof positions of the given accessor. More...
 
void loadResultData (GmVector &data)
 Loads the contents of the appropriate state vars / nodes to the vector data according to the mapping between vector lines and the corresponding degree of freedom.
 
void loadFixedData (GmVector &data)
 Loads the contents of the appropriate state vars / nodes to the vector data according to the mapping for FIXED degrees of freedom. Vector size should be equal to numTotalDof() - numDof(). Can be used only if the mode parameter given in the call to init() was equal to REMOVE_DOF.
 
const QVector< int > & remapInfo () const
 Returns the remap vector storing information on how to translate values from a vector prior to an update() call to a new vector after the call, when dofs where added / removed. More...
 
void remapData (GmVector &data, double newValue=0.0) const
 Given a vector 'oldData', with size equal to the number of dofs before the last call to update(), resizes and adjusts its contents to the new size equal to the current number of dofs. Places in the updated vector corresponding to new degrees of freedom will be filled with 'newValue'. More...
 
void remapData (const GmVector &oldData, GmVector &newData, double newValue=0.0) const
 Given a vector 'oldData', with size equal to the number of dofs before the last call to update(), transfers its contents to the new vector 'newData' with size equal to the current number of dofs. Places in the 'newData' vector corresponding to new degrees of freedom will be filled with 'newValue'. More...
 
bool fastUpdate () const
 Returns whether the assembler was initialized with the fast update option turned on or off.
 

Private Member Functions

bool addMappedMatrixData (int ndof, const int *indexMap, const GmDofMap *dofMap, const GmpFemAssemblerMatrixCombiner *matCombiner, const GmpFemAssemblerVectorCombiner *vecCombiner)
 Helpper function for addElementData() and addLocalMatrixData()
 
bool initSparseLayout (GmSparseMatrixLayoutBuilder *layoutBuilder)
 Traverses mesh elements & physics to discover sparse matrix layouts.
 
void cleanup ()
 Clear memory resources & used structures.
 
bool initTranslationStructures (bool enableFastUpdate, bool createReverseMapping, QString &err)
 Initializes all data structures used to speed up local to global index translations as described in the class documentation. More...
 
bool initTemplateAssociations (const GmElementDofMap *nodeDofMap, const QHash< int, GmElementDofMap > &fixedNodeDofMap, QString &err)
 Init per node template association for calculating dof offsets & other informations. More...
 
void initDofToNodeIndex (int firstNode)
 Initializes the reverse mapping storing node number per matrix dof.
 
void initUpdateDofPhysicsList ()
 Fills the _nonConstantDofPhysics with those physics that can be used in a 'fast update' mode.
 
GmElementDofMapcreateNodeDofMap (GmElementDofMap *physicsDofMapList, int *maxElemDof, size_t *dofsSum)
 Creates a vector storing for each node its dof map combining information from every element / physics that shares each node. More...
 
bool fillFixedDofMap (QHash< int, GmElementDofMap > &fixedDofMap, const GmElementDofMap *fullNodeDofMap, const QHash< int, GmElementDofMap > *partialNodeDofMap, QString &err)
 Initializes the _fixedValues and the _fixedNodeIndex structures with fixed DOF information obtained from the physics. Also fills the. More...
 
bool fillBcNodeDofMap (const QList< GmpFemPhysics * > &physicsList, QMap< QPair< int, int >, double > &fixedData, bool saveNonConstantPhysics, QString &err)
 Traverses the given set of physics, gathering fixed dof information and merging the result in the fixedData map, whose keys equal (node, dof) and value equals the fixed dof value. More...
 
bool fillUnitMap (const QList< GmpFemPhysics * > &physicsList, GmElementDofMap *physicsDofMapList, QString &err)
 Fills the global list _dofUnitList, indexed by dof number, storing the unit in which each dof associated with the given physics are calculated. Also stores int _timeUnit the unit in which physics handle time. More...
 
int fecthNodeTemplate (const GmElementDofMap &nodeDofMap, const GmElementDofMap &fixedNodeDofMap, int node, QString &err)
 Given a node dof map + the node fixed dof map, returns the index of the corresponding template entry in _templateList. More...
 
int fecthDofStatevar (int dof)
 Given a dof number, returns the index of the corresponding state var entry in _stateVars. More...
 
void fillFixedData (double *data, int size)
 The worker function for loadFixedData()
 
bool updateFixedBcValues (QString &err)
 Updates the assembler internal view of non constant fixed value boundary conditions. More...
 
bool updateTranslationStructures (const QList< GmpFemPhysics * > &physList, const QList< QList< int > > &elementLists, bool nodesAdded, int updateFixedMode, bool updateUnits, QString &err)
 Updates the needed data structures following an update / addition of the given elements, considering only the set of physics given in physList. More...
 
QHash< int, GmElementDofMapcreateUpdateNodeDofMap (const QList< GmpFemPhysics * > &physList, const QList< QList< int > > &elementLists, GmElementDofMap *physicsDofMapList, int *maxElemDof, size_t *dofsSum)
 Creates a QHash with node dof mapings for the nodes in the given element lists, according to their associated physics. More...
 
int updateNodeTemplate (int node, const GmElementDofMap &addedDofMap, const GmElementDofMap &newFixedDofMap, bool useOldFixedDofMap, int *nFreeAdded, int *nFixedAdded, QString &err)
 Returns the updated template index for a node that had the dofs described in addedDofMap added to its current profile and follows the new given fixedDofMap (unless useOldFixedDofMap is set to true when newFixedDofMap is ignored). More...
 
int extendedDofIndex (int node, int dof) const
 Simmilar to dofIndex(), this function differs from it only when it receives a a fixed degree of freedom (and mode passed to init() was REMOVE_DOF). More...
 
bool saveNodeDofInfo (QVector< QPair< int, int > > &map, QString &err) const
 Fills the given map with information about the node and dof number of each index in the current global matrix. Returns false and fills err if the map could not be resized (to the current number of dofs + 1). More...
 
bool buildRemapInfo (const QVector< QPair< int, int > > &oldMap, QString &err)
 Given a map filled by saveNodeDofInfo() before translation structures update, fills the _remapInfo vector. It also sets _hasRemapInfo to true. More...
 
int tgi (int node) const
 Translate Ghost Index - An auxiliary function that given a mesh node index, whose value might be a ghost node index with its high bit set, returns the linear index used to address the _nodexIdex and _nodeTemplateIndex vectors.
 
int nodeIndex (int node) const
 Returns the base node index of the given node. Node can be a ghost node index if ghost support is enabled.
 
void setNodeIndex (int node, int baseIndex)
 Updates the base node index of the given node. Node can be a ghost node index if ghost support is enabled.
 
void incNodeIndex (int node, int offset)
 Updates the base node index of the given node by adding the given value. Node can be a ghost node index if ghost support is enabled.
 
int nodeTemplateIndex (int node) const
 Returns the node template index. Node can be a ghost node index if ghost support is enabled.
 
void setNodeTemplateIndex (int node, int templateIndex)
 Updates the template index of the given node. Node can be a ghost node index if ghost support is enabled.
 
const TemplateDatanodeTemplate (int node) const
 Returns the node template. Node can be a ghost node index if ghost support is enabled.
 

Private Attributes

GmElementMesh_mesh
 The mesh.
 
const QList< GmpFemPhysics * > & _physicsList
 List of physics that will cooperate to generate the global matrix.
 
FixedDofMode _fixedDofMode
 How does the assembler treats fixed degrees of freedom.
 
GmpFemMatrixSet_matSet
 The matrix set used together with the assembler.
 
GmpFemVectorSet_vecSet
 The vector set used together with the assembler.
 
GmpFemLocker_locker
 The locker object used when adding values to the global matrices/vectors in multiple threads.
 
GmpFemAssemblerMatrixAdder_matAdder
 The configured adder object used to assemble matrices.
 
GmpFemAssemblerVectorAdder_vecAdder
 The configured adder object used to assemble vectors.
 
GmpFemAssemblerMatrixAdder_debugMatAdder
 The configured adder object used to assemble matrices in "single debug" mode.
 
GmpFemAssemblerVectorAdder_debugVecAdder
 The configured adder object used to assemble vectors in "single debug" mode.
 
int _nthreads
 The number of threads configured for the assembler. Can be any value >= -1.
 
int _nt
 The number of threads that should be used by the assembler in OMP operations.
 
bool _ghostEnabled
 Should the assembler handle ghost nodes?
 
int _nMatrixDof
 The number of degrees of freedom of the assembled matrix.
 
int _nTotalDof
 The total number of degrees of freedom of the system, including fixed ones.
 
int _maxNodeDofs
 The maximum number of dofs for a node. Might over estimate if number of node dofs are reduced over time.
 
size_t _elemDofsSum
 A sum of the number of local matrix entries for the whole set of elements. When update() is called, this value might over estimate the real value since elements with added dofs are then considered again.
 
int _nNodes
 The number of nodes in _nodeIndex & _nodeTemplateIndex. Will be equal to _mesh->numNodes() or _mesh->totalNumNodes() depending on _ghostEnabled.
 
int * _nodeIndex
 A vector storing for each mesh node its base index in the matrix.
 
QVector< int > _dofNode
 A reverse vector mapping, storing for each dof in the assembled matrix the node it belongs to, with size equal to _nMatrixDof (or 0 if disabled).
 
unsigned char * _nodeTemplateIndex
 A vector storing for each mesh node its dof template index.
 
QVector< double > _fixedValues
 A vector storing the values of the fixed dofs with size equal to _nTotalDof - _nMatrixDof.
 
QVector< double > _oldFixedValues
 A vector storing the values in _fixedValues before the last call to update()
 
QHash< int, int > _fixedNodeIndex
 Map determining the base index of a node in the _fixedValues vector.
 
QVector< TemplateData_templateList
 A list with node metadata information for every node template type.
 
QVector< QPair< GmStateVar *, GmValueAccessor * > > _stateVars
 A list storing pairs of state vars and accessors to them for every different state var that provides degrees of freedom for the assembled matrix.
 
GmTLBuffer< int, true > _elemIndexMap
 An auxiliar vector used by addElementMatrix, with size equal to the biggest local matrix number of dofs possible. Not a simple vector as it can grow when dofs are added to the system.
 
QList< GmpFemPhysics * > _nonConstantFixedValuesPhysics
 A list with the set of physics whose fixed boundary condition values are not constant.
 
QList< GmpFemPhysics * > _nonConstantDofPhysics
 A list with the set of physics that should be queried when updating dofs if _updateDofMode == FAST_UPDATE.
 
QMap< QPair< GmElementDofMap, GmElementDofMap >, int > _templateCache
 A map keyed by both the node dof mapping and the fixed dof mapping, storing the equivalent template index. This map holds values outside the life cycle of initTranslationStructures() only when _updateMode == FAST_UPDATE.
 
QMap< int, int > _dofToSvCache
 A map keyed by dof number storing the corresponding entry in _stateVars. This map holds values outside the life cycle of initTranslationStructures() only when _updateMode == FAST_UPDATE.
 
Unit_dofUnitList
 A vector keyed by dof number, with size equal to GM_MAX_DOF, storing the corresponding dof unit. This vector holds values outside the life cycle of initTranslationStructures() only when _updateMode == FAST_UPDATE.
 
Unit _timeUnit
 The time unit used by all the physics.
 
bool _hasRemapInfo
 Is remap info available (The last call to update included the SAVE_MAPPING_INFO flag)? Needed since a _remapInfo vector with size 0 is conceptually possible if numDof() == 0.
 
QVector< int > _remapInfo
 The remap vector storing information on how to translate values from a vector prior to an update() call to a new vector after the call, when dofs where added / removed. Its contents is valid only if _hasRemapInfo is true. Its size is equal to the current number of dofs. For each dof it stores the index of that dof in the previous vector or -1 if it is a new dof.
 

Detailed Description

A class responsible for assembling the global stiffness matrix, handling multiple physics with possibly diferent degrees of freedom for each node.

To handle multiple degrees of freedom per node, per physics, the assembler starts by defining which dof are active for each node, combining the results of the multiple physics in play. This is done by traversing the mesh and for each elements node, "adding" together the node's dof, as seen by each physics. The result is a vector with a GmElementDofMap entry per mesh node.

Counting the number of ones in each GmElementDofMap entry we reach the global number of dof in the assembled matrix.

The second step involves creating a second array, with an entry per mesh node, keeping the first line/column number of a node in the matrix. As an example, imagine a mesh with 4 nodes having 1, 3, 1 and 2 dof respectively. In this cenario, the node index array will be equal to [0,1,4,5].

The third step is to associate to each pair <element type, physics index> a vector that gives us for each local element matrix line the local node number and the global dof that that line is bound to (notice that we only need to fill this map for the element types that exists in the mesh).

When translating a local element matrix column into its global position, we can combine the previous information to complete the desired translation. Given the local node number of a local element matrix line, we can query the element to get the global node number for that element and use it to index the node index array (created in step two) to reach the base line / column for that node in the global matrix.

The last needed information is the offset from this base number to reach the desired dof associated to this node. The needed information exists in the element dof map associated to each node, encoded in a 'bit stream'. If we want to know the offset of dof 34, we need to count the number of 1 bits before bit 34 (which should be equal to one). This can be done by masking the higher bits and counting the left over 1 bits using a fast count bits implementation (qPopulationCount() offers one possible version). A better option is to recognise that most nodes will follow the same dof pattern and create an array for each existing pattern providing a direct translation. The dof map is then replaced by a vector that ties each node to its template. This option is both faster and uses less memory since the template index vector can have 8 bit entries against 64 bit entries for the dof map. This is the fourth and last step preparing data structures for local to global index translations.

Member Enumeration Documentation

◆ FixedDofMode

Enumerator
NO_DOF_HANDLING 

No fixed dof handling. The full matrix will be generated.

REMOVE_DOF 

Fixed dofs are removed from the matrix.

◆ GhostDofMode

Enumerator
GHOST_SUPPORT 

The assembler will support dofs in ghost nodes.

NO_GHOST_SUPPORT 

The assembler will not support dofs in ghost nodes.

AUTO_GHOST_SUPPORT 

Equal to GHOST_SUPPORT if there is ANY state var in the mesh with ghost support. Equal to NO_GHOST_SUPPORT otherwise.

◆ UpdateMode

Enumerator
NODES_ADDED 

Nodes (or ghost nodes) where added to the mesh. Should be combined with ELEMENTS_ADDED and/or DOFS_ADDED. Fast update possible depending on remaining flags.

DOFS_ADDED 

Dofs where added to existing elements. Fast update possible depending on physics support and remaining flags.

DOFS_REMOVED 

Dofs where removed from existing elements. No fast update.

ELEMENTS_ADDED 

Elements where added to the mesh. Fast update possible depending on remaining flags.

ELEMENTS_REMOVED 

Elements where removed from the mesh. Fast update possible depending on mesh topology support and remaining flags. (NOT IMPLEMENTED YET)

FIXED_BCS_UPDATE 

Fixed boundary conditions marked as non constant should be updated. No other changes where done.

FIXED_BCS_CHANGED 

Fixed boundary conditions where changed but only for added elements and/or elements with added dofs.

FIXED_BCS_ALL 

Fixed boundary conditions should be fully reloaded.

FIXED_BCS_KEEPOLD 

When combined with the previous BCS flags, informs the assembler that old fixed bc data (_oldFixedValues) should not be updated with the current values.

FULL_UPDATE 

The full assembler mapping will be rebuilt, regardless of other flags.

SAVE_MAPPING_INFO 

The information needed to remap vectors after an update will be saved.

Member Function Documentation

◆ addElementData()

const int * GmpFemAssembler::addElementData ( const GmElement e,
int  physicsIndex,
const GmpFemAssemblerMatrixCombiner matCombiner,
const GmpFemAssemblerVectorCombiner vecCombiner 
)

Adds data from a set of local element matrix & vectors to their respective global data. Depending on the mode used for creating the matrix and the vector sets (the mode might be different), the assembled data will be added to matching global matrices/vectors and/or to the equivalent matrix/vector. Enabling both of them is usefull for debugging purposes.

Given a local set of stiffness matrices in '_matSet', created by physics 'physicsIndex' for the element 'e', assembles those matrices into their respective global matrices if their filled attribute is true. Global matrices are stored in the matrix sets and paired to their local sisters by their type. The matrix symmetric attribute should tell if the local matrix is symmetric or not.

The combiner objects are used for defining the composition rule when the matrix and/or the vector mode includes an equivalent matrix. If turned on, fixed DOF transposition is done after calculating the equivalent value. The combiner objects should be a subclass of the GmpFemAssemblerMatrixCombiner and GmpFemAssemblerVectorCombiner types, providing the rules used for elemental matrices / vectors composition. See the comments on those classes for an explanation of the virtual methods that should be reimplemented.

The same operations apply to vectors stored inside '_vecSet'

If the mode parameter given in the call to init() was equal to REMOVE_DOF, each matrix set might have a matching vector set (associated in the matSet). When a local matrix index corresponds to a fixed degree of freedom, that fixed dof "stiffness" will be multiplied by the dof fixed value and the result subtracted from the matching force vector, allowing the dof to be removed from the global matrix. If no association exists, this transposition step won't be done.

This function returns a pointer to a internal vector storing the translation vector mapping local indices onto global indices. Its size is equal to the elements matrix/vector dimension and its contents is kept only until another call to addElementData(). The vector semantics is the same as the vector filled by fillIndexMap(). It returns NULL on errors.

IMPORTANT: All the calls to addElementData() should be placed between a pair of calls to beginAssembly() and endAssembly().

IMPORTANT2: Although begin and endAssembly must be called from the main thread, addElementData() can be called from task threads. Before that, it is necessary to call setPartitions() to setup the needed locks.

◆ addLocalMatrixData()

bool GmpFemAssembler::addLocalMatrixData ( GmMatrixDof localDofMap,
const GmpFemAssemblerMatrixCombiner matCombiner,
const GmpFemAssemblerVectorCombiner vecCombiner 
)

Simmilar to addElementData() but instead of the matrices in the matrix and vector set being tied to an element, they are tied to a local matrix with arbitrary size. The localDofMap parameter gives the matrix layout. It also provides a vector that should be filled with the local to global dof translation map.

IMPORTANT: This function assumes that the dof mappings defined in localDofMap are already known to the assembler, that is, no new node-dof assignments can be created by this dof map. See the documentation for GmpFePhysics::fillContactData() for additional details.

◆ beginAssembly()

bool GmpFemAssembler::beginAssembly ( )

Notifies the assembler that the assembly process will begin. Also initializes global active vectors and matrices with zeros.

This call should be closed by a call to endAssembly(). Between this pair of calls, expects only calls to addElementData().

◆ buildRemapInfo()

bool GmpFemAssembler::buildRemapInfo ( const QVector< QPair< int, int > > &  oldMap,
QString err 
)
private

Given a map filled by saveNodeDofInfo() before translation structures update, fills the _remapInfo vector. It also sets _hasRemapInfo to true.

IMPORTANT: The current algorithm works based on the fact that new nodes are always added to the "end" of the node list.

◆ createNodeDofMap()

GmElementDofMap * GmpFemAssembler::createNodeDofMap ( GmElementDofMap physicsDofMapList,
int *  maxElemDof,
size_t *  dofsSum 
)
private

Creates a vector storing for each node its dof map combining information from every element / physics that shares each node.

This function does NOT take into account the possibility of fixed DOF values. It will fill the physicsDofMapList parameter with a map with all dofs used by each physics. It expects physicsDofMapList to be a vector with size equal to the number of available physics. It also fills the maxElemDof parameter with the maximum number of degrees of freedom for a mesh element / single physics. The dofsSum parameter is filled with the sum of the number of entries for local matrices of each active element.

◆ createUpdateNodeDofMap()

QHash< int, GmElementDofMap > GmpFemAssembler::createUpdateNodeDofMap ( const QList< GmpFemPhysics * > &  physList,
const QList< QList< int > > &  elementLists,
GmElementDofMap physicsDofMapList,
int *  maxElemDof,
size_t *  dofsSum 
)
private

Creates a QHash with node dof mapings for the nodes in the given element lists, according to their associated physics.

This function will include in the hash the dof map combining information from the given elements / physics. Only nodes belonging to an element will be included. For each included node, the returned information will include only those dofs seen by the given elements / physics. That means that the resulting map should be combined with the existing node information for completeness. This function does NOT take into account the possibility of fixed DOF values.

The physList parameter should be a list with the set of physics that will be considered by this function. Its value will probably be equal to either _nonConstantDofPhysics or _physicsList.

The elementLists parameter should contain, per physics, the set of considered elements. In that way it should have a size equal to the size of physList. Optionally, this list can contain a single entry. In that case this single element list will be used for each physics in physList.

If physicsDofMapList is different from NULL, it should be a vector with size equal to the size of physList and it will be filled with a map with all dofs used by each physics (as per this element view).

The maxElemDof parameter will be filled with the maximum number of dofs for this elements/physics while the dofsSum parameter will be filled with the sum of the number of entries for local matrices of each element/physics.

◆ dofIndex()

int GmpFemAssembler::dofIndex ( int  node,
int  dof 
) const

Given a global node index and a global degree of freedom index, returns its line/column in the assembled matrix or a negative value if the pair does not belong to it.

Keep in mind that this function will return -1 if the given node does not include the given dof and -2 for degrees of freedom reported as fixed by some physics when the mode parameter given in the call to init() was equal to REMOVE_DOF.

The node parameter can be a common geometry node index or a ghost node index (either sequentially numbered after common geometry indices or zero based with the high bit set)

◆ extendedDofIndex()

int GmpFemAssembler::extendedDofIndex ( int  node,
int  dof 
) const
private

Simmilar to dofIndex(), this function differs from it only when it receives a a fixed degree of freedom (and mode passed to init() was REMOVE_DOF).

For fixed dofs, instead of returning -2, this function will return a negative value 'v', such that abs(v+1) gives the position in the _fixedValues vector where the fixed dof value can be obtained.

This function should NOT be called with a unrelated node/dof pair.

The node parameter can be a common geometry node index or a ghost node index (either sequentially numbered after common geometry indices or zero based with the high bit set)

◆ fecthDofStatevar()

int GmpFemAssembler::fecthDofStatevar ( int  dof)
private

Given a dof number, returns the index of the corresponding state var entry in _stateVars.

This function will first look-up the dof index by querying _dofToSvCache. If not found, the corresponding state var will be fectched from the mesh, an accessor created and the _stateVars list updated, along with the index cache.

This function expects the global _dofUnitList list to be correctly filled.

◆ fecthNodeTemplate()

int GmpFemAssembler::fecthNodeTemplate ( const GmElementDofMap nodeDofMap,
const GmElementDofMap fixedNodeDofMap,
int  node,
QString err 
)
private

Given a node dof map + the node fixed dof map, returns the index of the corresponding template entry in _templateList.

This function will first look-up the template index by querying _templateCache. If not found, a new template will be created and the _templateList list updated, along with its index cache.

Since it might call fecthDofStatevar(), this function expects the global _dofUnitList list to be correctly filled.

Returns -1 and fills err in the event of a conflict between a state var data storage type and the binding of the said state var dof to a node (i.e, a state var bound to a geometry node only has its dof tied to a ghost node or vice versa). The node parameter is only used by this check.

◆ fillBcNodeDofMap()

bool GmpFemAssembler::fillBcNodeDofMap ( const QList< GmpFemPhysics * > &  physicsList,
QMap< QPair< int, int >, double > &  fixedData,
bool  saveNonConstantPhysics,
QString err 
)
private

Traverses the given set of physics, gathering fixed dof information and merging the result in the fixedData map, whose keys equal (node, dof) and value equals the fixed dof value.

This function returns false on errors while querying a physics. Even when returning true, err might be filled with warning messages.

If saveNonConstantPhysics is true, this function will append any physics who reports that fixed values can change over time to the _nonConstantFixedValuesPhysics list attribute.

Node values in the map key are 'linearized', i.e, eventual ghost nodes are translated with tgi()

◆ fillFixedDofMap()

bool GmpFemAssembler::fillFixedDofMap ( QHash< int, GmElementDofMap > &  fixedDofMap,
const GmElementDofMap fullNodeDofMap,
const QHash< int, GmElementDofMap > *  partialNodeDofMap,
QString err 
)
private

Initializes the _fixedValues and the _fixedNodeIndex structures with fixed DOF information obtained from the physics. Also fills the.

given hash table defining, for each mesh node that has a fixed dof, which dofs are fixed.

If _fixedValues or _fixedNodeIndex are not empty before this call, the structures will be overwritten.

If there was any error while retrieving fixed dof information from the physics, returns false filling err with an explanation. Err can also be filled with warning messages even when a valid vector is returned by the function.

The fullNodeDofMap and partialNodeDofMap parameters are used only to check the consistency between dof and fixed boundary condition information given by the physics. One of them should be NULL while the other not. This check is essential to prevent a fixed dof to be included in the maps wthout being a real degree of freedom. That can happend when a boundary condition points to a node that do not contain dofs. Imagine a boundary condition defined over a quadratic element side. It will contain nodes with no dofs in a super-parametric setting. It will also happend when elemnts are made inactive. When using a partial map during an assembler update, the missing information will be collected from the current node templates.

◆ fillIndexMap()

int GmpFemAssembler::fillIndexMap ( const GmElement e,
int  physicsIndex,
int *  indexMap,
const GmElementDof *&  dofMap 
)

Fills indexMap with the translation map from local to global matrix/vector indices and returns the number of entries filled in the map (which will be equal to the local element matrix/vector sizes). Also makes dofMap point to the physics/element dof map.

Expects that indexMap was allocated by the user and has enough size to receive the translation map. One can use the return value from maxElementDof() to allocate this vector with a size that is guaranteed to be sufficient.

If one of the elements dof corresponds to a fixed dof and the mode parameter given in the call to init() was equal to REMOVE_DOF, the map will contain a NEGATIVE value indicating that this local node value does NOT have a mapping in the global matrix.

Calling 'v' the negative value returned for a fixed dof, the expression abs(v+1) gives the position in the _fixedValues vector where the fixed dof value can be obtained.

Returns 0 if the element type is not supported by the given physics. In that case, neither indexMap nor dofMap are changed.

◆ fillUnitMap()

bool GmpFemAssembler::fillUnitMap ( const QList< GmpFemPhysics * > &  physicsList,
GmElementDofMap physicsDofMapList,
QString err 
)
private

Fills the global list _dofUnitList, indexed by dof number, storing the unit in which each dof associated with the given physics are calculated. Also stores int _timeUnit the unit in which physics handle time.

If there are inconsistencies between physics, returns false and fills err with a relevant message. The parameter physicsDofMapList is a list storing the dof usage of each physics in physicsList

If the _dofUnitList or _timeUnit already has values, they will be used in the validation process. New dofs that didn't had a previous unit will be added to the list.

◆ init()

bool GmpFemAssembler::init ( GmpFemMatrixSet matSet,
GmpFemVectorSet vecSet,
FixedDofMode  fixedMode,
GhostDofMode  ghostMode,
bool  enableFastUpdate,
bool  createReverseMapping,
int  numThreads,
QString err 
)

Initializes the assembler. Must be called after object construction.

Returns false on errors, filling err with the reason. On warnings, returns true but also fills err with a warning message.

The matSet and vecSet parameters are used to register with the assemble the sets used to handle local and global vectors / matrices. The assembler assumes the task of configuring matrix and vector sizes within the sets by calling adjustSizes() at initialization and whenever that data changes (by a call to update()). Please, keep in mind that the assembler only keeps a reference to this sets and does NOT take ownership of them (will not free their memory).

The fixedMode parameter defines how the assembler will handle fixed degrees of freedom. If equal to NO_DOF_HANDLING, the full matrix will be assembled ignoring that some of the matrix degrees of freedom are indeed fixed. With that mode, the responsability of adjusting the matrix to handle fixed degrees of freedom belongs to the process.

If equal to REMOVE_DOF, the assembler will query the given physics to obtain the list of fixed degrees of freedom in the problem. Those dofs will be removed from the global matrix, which will have its size reduced by the number of fixed dofs.

When merging local matrices in the global matrix, contributions associated with fixed dofs will be subtracted from the matching force vector, if any, as configured in the matrix set by calls to GmpFemMatrixSet::setTransposedDofVector(). If no association exists for a matrix, the matrix and vectors will be reduced but no transposition will be made.

When the matrix is reduced, calls to saveResultData() will automatically fill the mesh state vars associated with the fixed dofs with the given fixed value.

The ghostMode flag controls whether the assembler will enable or not ghost nodes to be degrees of freedom. The automatic option means that support will be enabled if any of the mesh state vars supports ghost nodes.

When the enableFastUpdate flag is set to true, aditional structures will be filled by the initialization procedure to enable faster response when calling update(). This flag should be set to true only if there will be dof changes in the mesh, either by changes in element dofs or by adding/removing elements. When set to false, a call to update() for reevaluating dofs (not including fixed BC dofs) will be roughly equivalent to destroying and rebuilding the assembler. See the update() function documentation for aditional information of when seting enableFastUpdate to true is effective.

The createReverseMapping flag controls the need to create a reverse mapping that given a global dof number (index in the assembled matrices/vectors), returns the node associated to this dof. This index is used only by functions stateVar() and svGroupId() and so can be set to false if those functions aren't needed.

The numThreads parameter controls the number of threads that should be used by paralelized operations. A value of -1 means the maximum number of threads as seen by the thread manager.

◆ initTemplateAssociations()

bool GmpFemAssembler::initTemplateAssociations ( const GmElementDofMap nodeDofMap,
const QHash< int, GmElementDofMap > &  fixedNodeDofMap,
QString err 
)
private

Init per node template association for calculating dof offsets & other informations.

Returns false and fills err in the event of a conflict between a state var data storage type and the binding of the said state var dof to a node (i.e, a state var bound to a geometry node only has its dof tied to a ghost node or vice versa)

◆ initTranslationStructures()

bool GmpFemAssembler::initTranslationStructures ( bool  enableFastUpdate,
bool  createReverseMapping,
QString err 
)
private

Initializes all data structures used to speed up local to global index translations as described in the class documentation.

Returns false on errors, filling err with the reason. On warnings, returns true but also fills err with a warning message.

◆ numElementDof()

int GmpFemAssembler::numElementDof ( const GmElement e,
int  physicsIndex 
) const

Returns the number of degrees of freedom for the local matrix calculated by the referenced physics for the given element.

Since this method is used to define the number of columns in the local element matrix, it does NOT take into consideration fixed dofs.

Parameters
Theelement
Thephysics index in the physics list received in the object constructor
Returns
Returns the number of degrees of freedom or 0 if this physics does not supports the requested element.

◆ remapData() [1/2]

void GmpFemAssembler::remapData ( GmVector data,
double  newValue = 0.0 
) const

Given a vector 'oldData', with size equal to the number of dofs before the last call to update(), resizes and adjusts its contents to the new size equal to the current number of dofs. Places in the updated vector corresponding to new degrees of freedom will be filled with 'newValue'.

Must be called only if SAVE_MAPPING_INFO was enabled in the last call to update().

◆ remapData() [2/2]

void GmpFemAssembler::remapData ( const GmVector oldData,
GmVector newData,
double  newValue = 0.0 
) const

Given a vector 'oldData', with size equal to the number of dofs before the last call to update(), transfers its contents to the new vector 'newData' with size equal to the current number of dofs. Places in the 'newData' vector corresponding to new degrees of freedom will be filled with 'newValue'.

Expects 'newData' to have the desired size. Must be called only if SAVE_MAPPING_INFO was enabled in the last call to update().

◆ remapInfo()

const QVector< int > & GmpFemAssembler::remapInfo ( ) const

Returns the remap vector storing information on how to translate values from a vector prior to an update() call to a new vector after the call, when dofs where added / removed.

The vector size is equal to the current nuber of dofs. For each dof it stores the index of that dof in the previous vector or -1 if it is a new dof. The returned vector is valid only if SAVE_MAPPING_INFO was enabled in the last call to update(). Its contents remains valid only untill the next call to update().

◆ saveNodeDofInfo()

bool GmpFemAssembler::saveNodeDofInfo ( QVector< QPair< int, int > > &  map,
QString err 
) const
private

Fills the given map with information about the node and dof number of each index in the current global matrix. Returns false and fills err if the map could not be resized (to the current number of dofs + 1).

The filled map has one entrie per matrix column + one last enty with a key greater than any possible key to serve as a stop condition.

◆ saveResultData()

void GmpFemAssembler::saveResultData ( const GmVector data,
int  baseDof,
GmValueAccessor ac,
double  lambda = 1.0 
)

Saves the contents of the data vector to the given accessor. Only the subset in the source vector associated with the given dof (or dofs for vector based accessors) are saved to the accessor. The lambda factor will be applied only to fixed values handled by the assembler.

◆ saveResultDataToFixedDofs()

void GmpFemAssembler::saveResultDataToFixedDofs ( const GmVector data,
int  baseDof,
GmValueAccessor ac,
bool  clearFree 
)

Saves the contents of the data vector to the fixed dof positions of the given accessor.

Only the subset in the source vector associated with the given dof (or dofs for vector based accessors) are saved to the accessor. If clearFree is true, free dof positions on the accessor are set to 0.0.

◆ update()

bool GmpFemAssembler::update ( int  mode,
QString err 
)

Updates the assembler structure after a change in the model free or fixed dofs, according to the given mode. Depending on the requested mode and on the enableFastUpdate flag passed to init(), this might result in an update of only the affected dofs or in a full update, roughly equivalent to destroying and rebuilding the assembler (except that recreating the assembler will not create remaping information).

The mode flag is an 'or' of the flags defined in the UpdateMode enumeration. The following rules define its behaviour:

1) If mode includes NODES_ADDED the relevant vectors (_nodeIndex and _nodeTemplateIndex) will be grown and initialized as if the new nodes have no dofs. Their new dofs will be set by new elements or dofs given by DOFS_ADDED or ELEMENTS_ADDED flags.

2) If mode includes FULL_UPDATE, DOFS_REMOVED or the assembler has been initialized with the flag enableFastUpdate set to false, the update will be a full update.

3) If mode includes DOFS_ADDED, the list of elements considered in the update process will be given by the set of physics that report an element based dof behaviour and also that dofs change, are trackable and are only added to the system (see GmpFemPhysics::dofByElement()), through calls to GmpFemPhysics::changedElements().

4) If mode includes ELEMENTS_ADDED, the list of elements considered in the update process will be given by the return of GmCellMesh::addedCells().

5) If mode include both DOFS_ADDED and ELEMENTS_ADDED, the sets described in steps 3) and 4) will be combined.

6) If mode includes ELEMENTS_REMOVED and the mesh does not support topological queries, a full update will be performed, regardless of other flags. If support exists, the set of removed nodes will be obtained by determining, for each node in the removed cells, which of them do not belong to any other cells. Also, nodes that still belong to a cell will have their dofs recomputed. (IMPORTANT: fast mode not implemented yet)

7) If mode includes FIXED_BCS_ALL, the information on fixed dofs will be recomputed for every mesh node, independently of the set of elements considered in steps 3) to 6). Else, if mode includes FIXED_BCS_CHANGED, the information on fixed dofs will be recomputed only for the selected elements (so this option only really makes sense if DOFS_ADDED or ELEMENTS_ADDED where also specified). Note that in this case, fixed dof values will be updated for all fixed dofs.

Only the adition or removal of fixed dofs is constricted to that set of elements. If mode includes FIXED_BCS_UPDATE (and does not include the other BC flags), fixed dof values will be updated for all nodes, assuming that dof aplication points have not changed at all.

8) If mode includes SAVE_MAPPING_INFO, information to enable vector data remaping will be saved, regardless of the update effective mode (fast or full). Including this flag enables the use of the remapInfo() and remapData() calls.

When considering the mode used for updates, one should consider that the fast update mode is effective mainly when the number of mesh nodes/cells that have changed dofs during the simualtion evolution is small in comparisson with the total number of nodes/cells. If this is not true, measurement must be made to check if the fast mode is worthy or if a full update should be done instead.

If, as a result of this function, the number of free, fixed or per element maximum dofs changes, the matrix and vector sets will be resized accordingly.

This function returns false on errors, filling err with the reason. On warnings, it returns true but also fills err with a warning message.

IMPORTANT: This function will clear mesh changes and also physics changes by calling GmCellMesh::clearAddedCells() and GmpFemPhysics::clearChangedElements().

◆ updateFixedBcValues()

bool GmpFemAssembler::updateFixedBcValues ( QString err)
private

Updates the assembler internal view of non constant fixed value boundary conditions.

Since it will reload information only from non constant boundary conditions, as reported by each physics through GmpFemPhysics::fixedNodalDofsBc(), this function can be freelly called with no overhead if all the fixed value boundary conditions are constant.

IMPORTANT: This function expects that bc application points have not changed.

◆ updateNodeTemplate()

int GmpFemAssembler::updateNodeTemplate ( int  node,
const GmElementDofMap addedDofMap,
const GmElementDofMap newFixedDofMap,
bool  useOldFixedDofMap,
int *  nFreeAdded,
int *  nFixedAdded,
QString err 
)
private

Returns the updated template index for a node that had the dofs described in addedDofMap added to its current profile and follows the new given fixedDofMap (unless useOldFixedDofMap is set to true when newFixedDofMap is ignored).

If needed, a new template will be created. This function assumes that _templateCache, _dofToSvCache and _dofUnitList remain filled since the assembler initialization.

Parameters nFreeAdded and nFixedAdded are filled with the number of new dofs really added to the node (which will be different from the dofs in the maps since some of those dofs will probably be already maped to the node).

On the event of a conflict between a state var data storage type and the binding of the said state var dof to the node (i.e, a state var bound to a geometry node only has its dof tied to a ghost node or vice versa), fills err and returns -1;

◆ updateTranslationStructures()

bool GmpFemAssembler::updateTranslationStructures ( const QList< GmpFemPhysics * > &  physList,
const QList< QList< int > > &  elementLists,
bool  nodesAdded,
int  updateFixedMode,
bool  updateUnits,
QString err 
)
private

Updates the needed data structures following an update / addition of the given elements, considering only the set of physics given in physList.

IMPORTANT: This function considers that dofs where added to the system by the addition or by changes in the listed elements. It will NOT work if dofs where removed or added by other elements. Likewise, it expects that dofs that where previously marked as fixed will remain so. New fixed dofs can be handled depending on the given set of flags.

The physList parameter should be a list with the set of physics that will be considered by this function when querying element dofs. Its value will probably be equal to either _nonConstantDofPhysics or _physicsList.

The elementLists parameter should contain, per physics, the set of considered elements. In that way it should have a size equal to the size of physList. Optionally, this list can contain a single entry. In that case this single element list will be used for each physics in physList. If the elementLists is empty, this function will have an effect only if updateFixedMode == FIXED_BCS_ALL.

The supplied flags control whether fixed boundary conditions and their values, or units, need update or not. When updateFixedMode is different from zero and the fixed dof mode is equal to REMOVE_DOF, fixed BC values will be recomputed for every node and fixed BC application points will be reevaluated for the set of elements in elementList, if updateFixedMode is equal to FIXED_BCS_CHANGED, or for all nodes if equal to FIXED_BCS_ALL.

The updateUnits flag tells us that the list of used dof units should be updated. It should be set to false only if there is a gurantee that no different dof types where added to the system.

Returns false on errors, filling err with the reason. On warnings, returns true but also fills err with a warning message.


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