FemProcess
The GeMA Fem Process Plugin
Public Types | Public Member Functions | Protected Types | Protected Slots | Protected Member Functions | Protected Attributes | List of all members
GmpFemSolver Class Reference

Basic class for the FEM solving process. More...

#include <gmpFemSolver.h>

Inheritance diagram for GmpFemSolver:
Inheritance graph
[legend]
Collaboration diagram for GmpFemSolver:
Collaboration graph
[legend]

Public Types

enum  PrintOptions {
  PRINT_ELEMENT_MATRICES = 0x001, PRINT_ELEMENT_VECTORS = 0x002, PRINT_ELEMENT_DOF_MAPPING = 0x004, PRINT_EQ_MATRIX = 0x008,
  PRINT_EQ_VECTOR = 0x010, PRINT_GLOBAL_MATRICES = 0x020, PRINT_GLOBAL_VECTORS = 0x040, PRINT_LINEAR_MATRIX = 0x080,
  PRINT_LINEAR_VECTOR = 0x100, PRINT_LINEAR_RESULT = 0x200
}
 Flag options used to compose the printOptions constructor parameter. More...
 

Public Member Functions

 GmpFemSolver (GmElementMesh *mesh, GmSimulationData *simulation, const QList< GmpFemPhysics * > &physics, GmNumSolver *solver, const GmpFemSolverOptions &options, const GmLogCategory &logger)
 Constructor. Expects to receive as parameters the mesh we are acting upon, the list of physics objects that will cooperate to create the global stiffness matrix and the solver used to solve the resulting linear system. More...
 
virtual ~GmpFemSolver (void)
 Destructor.
 
virtual bool init ()
 Prepares the solver for assembling matrices by creating the assembler object and allocating the needed memory. This operations only need to be done once if the solver is used in a loop. More...
 
bool run ()
 Execute the process. Returns true on success.
 
bool calcLinearResidual (double *rnorm, double *maxNodeDiff, double *avgNodeDiff)
 When solving a non linear system by repeated iterations, this function aims to calculate the residual of the actual solution. More...
 
void setMatrixCombinerObject (const GmpFemAssemblerMatrixCombiner *combiner)
 Sets the combiner object that will be used by the assembler to merge elemental data from multiple matrices (eg. C and K) into an equivalent matrix. More...
 
void setVectorCombinerObject (const GmpFemAssemblerVectorCombiner *combiner)
 Sets the combiner objects that will be used by the assembler to merge elemental data from multiple vectors (eg. Fe and Fi) into an equivalent vector. More...
 
bool update (int mode, QString &err)
 Informs the solver of an external change in model dofs, elements or BCs that must be forwarded to the assembler. See the documentation of GmpFemAssembler::update() for more details.
 
virtual bool addStateItemsToGroup (GmStateDump *state, int groupId)
 Adds to 'state' the data items that should be saved for this FEM process. Should probably be overriden by derived classes to add their own requirements. More...
 
virtual bool stateAboutToBeSaved (GmStateDump *state)
 
virtual bool stateSaved (GmStateDump *state)
 
virtual bool stateAboutToBeLoaded (GmStateDump *state)
 
virtual bool stateLoaded (GmStateDump *state)
 Informs physics objects that a load operation was completed giving them an opportunity to adjust its state in case the loaded data is not exactly what the internal state needs to be. It will also call derivedResults() so that previous results can be recovered. In that call, the nonLinearSolver parameter is obtained from the stateLoadedCalcDerivedResultsNonLinear() virtual function.
 
virtual bool fillStateControlMapData (QVariantMap *map)
 Method called for filling the state dump control regustered in addStateItemsToGroup()
 
virtual bool stateControlMapDataLoaded (QVariantMap *map)
 Method called when the state dump control registered in addStateItemsToGroup() has been loaded.
 
virtual bool stateLoadedCalcDerivedResultsNonLinear ()
 Returns the value of the "nonLinearSolver" parameter that will be passed to the calcDerivedResults() call from the standard implementation of stateLoaded()
 
- Public Member Functions inherited from QObject
virtual const QMetaObjectmetaObject () const const
 
virtual void * qt_metacast (const char *)
 
virtual int qt_metacall (QMetaObject::Call, int, void **)
 
 QObject (QObject *parent)
 
virtual bool event (QEvent *e)
 
virtual bool eventFilter (QObject *watched, QEvent *event)
 
QString objectName () const const
 
void setObjectName (const QString &name)
 
bool isWidgetType () const const
 
bool isWindowType () const const
 
bool signalsBlocked () const const
 
bool blockSignals (bool block)
 
QThreadthread () const const
 
void moveToThread (QThread *targetThread)
 
int startTimer (int interval, Qt::TimerType timerType)
 
int startTimer (std::chrono::milliseconds time, Qt::TimerType timerType)
 
void killTimer (int id)
 
findChild (const QString &name, Qt::FindChildOptions options) const const
 
QList< T > findChildren (const QString &name, Qt::FindChildOptions options) const const
 
QList< T > findChildren (const QRegExp &regExp, Qt::FindChildOptions options) const const
 
QList< T > findChildren (const QRegularExpression &re, Qt::FindChildOptions options) const const
 
const QObjectList & children () const const
 
void setParent (QObject *parent)
 
void installEventFilter (QObject *filterObj)
 
void removeEventFilter (QObject *obj)
 
QMetaObject::Connection connect (const QObject *sender, const char *signal, const char *method, Qt::ConnectionType type) const const
 
bool disconnect (const char *signal, const QObject *receiver, const char *method) const const
 
bool disconnect (const QObject *receiver, const char *method) const const
 
void dumpObjectTree ()
 
void dumpObjectInfo ()
 
void dumpObjectTree () const const
 
void dumpObjectInfo () const const
 
bool setProperty (const char *name, const QVariant &value)
 
QVariant property (const char *name) const const
 
QList< QByteArraydynamicPropertyNames () const const
 
void destroyed (QObject *obj)
 
void objectNameChanged (const QString &objectName)
 
QObjectparent () const const
 
bool inherits (const char *className) const const
 
void deleteLater ()
 

Protected Types

enum  DisabledWarnings {
  PrescribedForceForInvalidDof, PrescribedForceForFixedDof, FixedBcForInvalidDof, ConflictingFixedBcValue,
  NumDisabledWarnings
}
 Enums describing the set of warnings that can be disabled by simulation options. More...
 

Protected Slots

virtual void meshChanged ()
 Slot called when the underlying mesh has been changed. More...
 

Protected Member Functions

bool initSolver (GmpFemAssembler::FixedDofMode assemblerMode, bool assemblerReverseMapping, bool enableFastUpdate)
 Basic implementation of the init function receiving as parameters the configuration options that are sent to the assembler. More...
 
virtual bool initElementSets (GmNumSolver *solver)
 Helpper function used to initialize the matrix / vector sets with the definition of the needed types and also the relationships of those sets with the assembler. More...
 
virtual bool initResultAttributes (QString prefix)
 Helpper function used to register the set of result attributes managed by the fem solver. To avoid name clashes, all of them should prepend the given prefix (usually "fem") to the attribute name. More...
 
virtual bool cleanup ()
 Dealocates memory and sets allocated resources to NULL. As a convenience, returns false.
 
GmpFemPhysics::FemResultType prepareMatrices (bool skipFixedBcs=false)
 Auxilliary function used to fill the K, and f matrices / vectors along with any other matrices/vectors stored in _matSet and _vecSet. More...
 
GmpFemPhysics::FemResultType saveElementData (int iter)
 Auxilliary function used to save to the configured file data from the set of element matrices and vectors. In this process, physics are called again to fill element matrices and vectors (according to the active matrices and vectors in _vecSet and _matSet) but the results are NOT added to the global matrices.
 
bool solveLinearSystem (bool xFilled)
 Solves the linear system K.x = f taking K and f from the single equivalent matrix/vector or from matrix K / vector f, depending on the assembler configuration. Results are saved in _x. Prints the used matrices / result as configured in _printOptions. More...
 
double timeConvert (double val)
 Auxiliary function created for being used by derived classes that converts a time value given in the current simulation time unit and returns that value converted to the physics expected unit.
 
GmpFemPhysics::FemResultType traverseElements ()
 Fills the stiffness matrix _K and the force vector _f by traversing elements asking physics for the local element matrix / force vector. More...
 
GmpFemPhysics::FemResultType traverseBoundaryElements ()
 Simillar to traverseElements(), this function steps through each element that belongs to an edge or face boundary condition, asking physics for additional contributions to the stiffness matrix _K and the force vector _f.
 
GmpFemPhysics::FemResultType traverseExternalLoads ()
 Simillar to traverseElements(), this function steps through each element associated to an external load force, asking physics for additional contributions to the force vector _f.
 
GmpFemPhysics::FemResultType traverseContactBoundaries ()
 Simillar to traverseElements(), this function steps through each possible contact from the given contact boundary conditions and asks the physics for additional contributions to the solver matrices and vectors. This is different from the other element based loops since the returned matrices are not restricted to a single element, but can combine nodes from several ones.
 
GmpFemPhysics::FemResultType traverseElementsForSaving (FILE *f, int iter)
 Similar to traverseElements() but instead of adding the elements to the assembler, saves the element data to the given file. More...
 
GmpFemPhysics::FemResultType fillElementData (const GmElement *e, int physIndex)
 Helper function used by traverseElements() to get data for a single element calling phys->fillElementData()
 
GmpFemPhysics::FemResultType fillElementBoundaryData (const GmElement *e, int physIndex, const GmBoundaryCondition *bc, int bcIndex, int bcListIndex, int border, const GmCellBoundary *b)
 Helper function used by traverseBoundaryElements() to get data for a single element calling phys->fillElementDataForBc()
 
GmpFemPhysics::FemResultType fillContactData (GmMatrixDof &localDofMap, const GmContactBoundaryCondition *cbc, int physIndex, int index1, int index2)
 Helper function used by traverseContactBoundaries() to get data for a single contact pair calling phys->fillContactData()
 
bool addFixedForces ()
 Adds fixed nodal forces, as seen by each physics, to the global vector _f.
 
virtual bool applyFixedBoundaryConditions (int *numFixed)
 Apply Dirichlet (fixed) boundary conditions to the global equation system. More...
 
bool collectFixedBcs (QVector< bool > &fixedRows, QVector< double > &fixedValues, QList< int > &dofIndex)
 Fills the set of vectors received as parameters with the complete set of fixed boundary conditions, as seen by the full set of physics in use, checking for possible conflicts in conditions. More...
 
virtual void collectGlobalContactPairs (const GmContactBoundaryCondition *cbc, QList< QPair< int, int > > &contactPairs)
 Global search rule used to define contact surfaces. Should be implemented to allow global contact conditions, filling contactPairs with the indices in cbc of the in-contact surfaces.
 
void printElementData (const GmElement *e, int ndof, const int *dofMapping, const GmpFemPhysics *p, const GmBoundaryCondition *bc=NULL, const GmContactBoundaryCondition *cbc=NULL, int loadId1=-1, int id2=-1)
 Print element dof mapping, local matrices and vectors depending on the current solver print options. More...
 
void printGlobalData (bool linear)
 Print global matrices and vectors, either after assembling or before solving the linear system.
 
- Protected Member Functions inherited from QObject
QObjectsender () const const
 
int senderSignalIndex () const const
 
int receivers (const char *signal) const const
 
bool isSignalConnected (const QMetaMethod &signal) const const
 
virtual void timerEvent (QTimerEvent *event)
 
virtual void childEvent (QChildEvent *event)
 
virtual void customEvent (QEvent *event)
 
virtual void connectNotify (const QMetaMethod &signal)
 
virtual void disconnectNotify (const QMetaMethod &signal)
 

Protected Attributes

int _solverId
 A unique index for this solver object used by the result attribute auto prefix feature.
 
GmElementMesh_mesh
 Mesh with the nodes and cells data.
 
GmSimulationData_simulation
 Simulation Data.
 
QList< GmpFemPhysics * > _physicsList
 List of physics that will cooperate to generate the system response.
 
GmNumSolver_solver
 Numeric solver that will be used to solve the set of linear equations.
 
const GmLogCategory_logger
 Logger for eventual messages.
 
QMutex _printMutex
 Mutex for serializing calls to printElementData when traversing elements in parallel.
 
bool _enabledWarn [NumDisabledWarnings]
 Vector storing which warnings are enabled.
 
GmpFemSolverOptions _solverOptions
 The set of solver options.
 
UnitConverter_timeConv
 Unit converter from the currentTimeUnit to the physics unit. Used only by derived solvers.
 
GmpFemAssembler_assembler
 Assembler used for adding element contributions to vectors and matrices.
 
const GmpFemAssemblerMatrixCombiner_matCombiner
 The combiner object used with the assembler to merge data from multiple elemental matrices.
 
const GmpFemAssemblerVectorCombiner_vecCombiner
 The combiner object used with the assembler to merge data from multiple elemental vectors.
 
GmpFemVectorSet _vecSet
 Set of global/element vectors filled by physics.
 
GmpFemMatrixSet _matSet
 Set of global/element matrices filled by physics.
 
GmVector _x
 State vector.
 
GmVector _r
 Residual vector.
 
unsigned _preCount
 Number of times that prepareMatrices() was called.
 
unsigned _runCount
 Number of times that run() was called.
 
QVariantMap _controlMap
 Variant map used for storing dump control data for the solver.
 
int _iterResAttr
 The registered id for the "fem.iter" result attribute.
 
int _iterResErrAttr
 The registered id for the "fem.iterErr" result attribute.
 

Additional Inherited Members

- Static Public Member Functions inherited from QObject
QString tr (const char *sourceText, const char *disambiguation, int n)
 
QString trUtf8 (const char *sourceText, const char *disambiguation, int n)
 
QMetaObject::Connection connect (const QObject *sender, const char *signal, const QObject *receiver, const char *method, Qt::ConnectionType type)
 
QMetaObject::Connection connect (const QObject *sender, const QMetaMethod &signal, const QObject *receiver, const QMetaMethod &method, Qt::ConnectionType type)
 
QMetaObject::Connection connect (const QObject *sender, PointerToMemberFunction signal, const QObject *receiver, PointerToMemberFunction method, Qt::ConnectionType type)
 
QMetaObject::Connection connect (const QObject *sender, PointerToMemberFunction signal, Functor functor)
 
QMetaObject::Connection connect (const QObject *sender, PointerToMemberFunction signal, const QObject *context, Functor functor, Qt::ConnectionType type)
 
bool disconnect (const QObject *sender, const char *signal, const QObject *receiver, const char *method)
 
bool disconnect (const QObject *sender, const QMetaMethod &signal, const QObject *receiver, const QMetaMethod &method)
 
bool disconnect (const QMetaObject::Connection &connection)
 
bool disconnect (const QObject *sender, PointerToMemberFunction signal, const QObject *receiver, PointerToMemberFunction method)
 
- Properties inherited from QObject
 objectName
 

Detailed Description

Basic class for the FEM solving process.

This class is tailored for solving a linear FEM problem expressed as an equation K.x = f, where the stiffness matrix K and the force vector f are filled by using the supplied physics objects.

Mesh elements are traversed and physics objects are used to collect their local contribution for the global matrix and force vector. In a second pass, boundary objects are traversed for their boundary contributions.

This class is internally organized as a set of functions to allow for the reuse of its structure in derived classes that implement other, more complicated, models.

Member Enumeration Documentation

◆ DisabledWarnings

Enums describing the set of warnings that can be disabled by simulation options.

Enumerator
PrescribedForceForInvalidDof 

Prescribed force over an invalid Dof warning at addFixedForces()

PrescribedForceForFixedDof 

Prescribed force over a fixed BC Dof warning at addFixedForces()

FixedBcForInvalidDof 

Fixed BC specified over an invalid Dof warning at collectFixedBcs()

ConflictingFixedBcValue 

Different fixed BC values specified for a single Dof at collectFixedBcs()

NumDisabledWarnings 

The number of entrie above.

◆ PrintOptions

Flag options used to compose the printOptions constructor parameter.

Enumerator
PRINT_ELEMENT_MATRICES 

All element matrices will be logged.

PRINT_ELEMENT_VECTORS 

All element vectors will be logged.

PRINT_ELEMENT_DOF_MAPPING 

Element DOF mapping will be logged.

PRINT_EQ_MATRIX 

The equivalent matrix built by the assembler in "single" mode or the K matrix in "match" mode will be logged.

PRINT_EQ_VECTOR 

The equivalent vector built by the assembler in "single" mode or the Fe vector in "match" mode will be logged.

PRINT_GLOBAL_MATRICES 

The set of all global matrices will be logged. In "single" mode they will be created for debugging purposes only.

PRINT_GLOBAL_VECTORS 

The set of all global vectors will be logged. In "single" mode they will be created for debugging purposes only.

PRINT_LINEAR_MATRIX 

The final matrix for the linear system will be logged.

PRINT_LINEAR_VECTOR 

The final vector for the linear system will be logged.

PRINT_LINEAR_RESULT 

The calculated linear system result (dof vector) will be logged.

Constructor & Destructor Documentation

◆ GmpFemSolver()

GmpFemSolver::GmpFemSolver ( GmElementMesh mesh,
GmSimulationData simulation,
const QList< GmpFemPhysics * > &  physics,
GmNumSolver solver,
const GmpFemSolverOptions options,
const GmLogCategory logger 
)

Constructor. Expects to receive as parameters the mesh we are acting upon, the list of physics objects that will cooperate to create the global stiffness matrix and the solver used to solve the resulting linear system.

The options parameter contains the set of parameters sent to the solver.

Member Function Documentation

◆ addStateItemsToGroup()

bool GmpFemSolver::addStateItemsToGroup ( GmStateDump state,
int  groupId 
)
virtual

Adds to 'state' the data items that should be saved for this FEM process. Should probably be overriden by derived classes to add their own requirements.

Auxiliar variant map dump item that notifies us for loading and saving data to the map

< Our "father" solver object

Implements GmGroupDumpItem.

Reimplemented in GmpFemNonLinearSolver, GmpFemNLSolver, and GmpFemTransientSolver.

◆ applyFixedBoundaryConditions()

bool GmpFemSolver::applyFixedBoundaryConditions ( int *  numFixed)
protectedvirtual

Apply Dirichlet (fixed) boundary conditions to the global equation system.

Apply boundary conditions, clearing rows and columns in the stiffness matrix and adjusting the global force vector as needed.

Follows the procedures described by Fellipa in sections 3.6.1 and 3.6.2 (pages 3-12 and 3-13). Our algorithm does the same thing as the one presented in figure 21.5 (page 21-7), but in a different way (hopefully more clear).

Returns in numFixed the number of degrees of freedom "removed" from the system.

Reimplemented in GmpFemNLSolver.

◆ calcLinearResidual()

bool GmpFemSolver::calcLinearResidual ( double *  rnorm,
double *  maxNodeDiff,
double *  avgNodeDiff 
)

When solving a non linear system by repeated iterations, this function aims to calculate the residual of the actual solution.

Consider a model K(x).x = f and an initial state x0, the model can be iterated so that we have K1 = K(x0), x1 = inv(K1).f, K2 = K(x1), ... More generally, we have Ki = K(xi-1), xi = inv(Ki).f.

With that model, our residual vector ri will be equal to K(xi).xi - f and the returned values: rnorm: The level two norm of the vector ri, ||ri|| maxNodeDiff: The maximum value of each entry in ri (makes sense only if all dof are comparable). Calculated only if maxNodeDiff != NULL. avgNodeDiff: The average value of each entry in ri (makes sense only if all dof are comparable) Calculated only if avgNodeDiff != NULL.

IMPORTANT: Keep in mind that this function expects to find in _x the result of the last call to run() and will by itself fill a new version of _K and _f, calculated based on this new value of _x.

The function returns false on errors (while evaluating K(x)).

◆ collectFixedBcs()

bool GmpFemSolver::collectFixedBcs ( QVector< bool > &  fixedRows,
QVector< double > &  fixedValues,
QList< int > &  dofIndex 
)
protected

Fills the set of vectors received as parameters with the complete set of fixed boundary conditions, as seen by the full set of physics in use, checking for possible conflicts in conditions.

Parameters
fixedRowsVector filled with entries set to true if a matrix row has a predefined value. Expects that the received vector has size equal to the system degrees of freedom, initialized with false.
fixedValuesVector filled with fixed values for entries where fixedRows[i] == true. Expects that the received vector has size equal to the system degrees of freedom.
dofIndexCompact list filled with indexes in fixedRows where fixedRows[i] == true

◆ init()

bool GmpFemSolver::init ( )
virtual

Prepares the solver for assembling matrices by creating the assembler object and allocating the needed memory. This operations only need to be done once if the solver is used in a loop.

When overriding this function, if you need to pass to the assmebler different configurations from the default (REMOVE_DOF & false), consider calling initSolver() in your reimplementation instead of init().

Reimplemented in GmpFemNonLinearSolver, GmpFemNLSolver, and GmpFemTransientSolver.

◆ initElementSets()

bool GmpFemSolver::initElementSets ( GmNumSolver solver)
protectedvirtual

Helpper function used to initialize the matrix / vector sets with the definition of the needed types and also the relationships of those sets with the assembler.

In a nutshell, this function is responsible for callint the initTypes() methods from matrix and vector sets. It is also its responsibility to call other configuration functions such as GmpFemMatrixSet::setTransposedDofVector(), GmpFemVectorSet::addTransposedDofVector() and GmpFemVectorSet::addSavedDofVector() as needed to drive the method that will be used by the assembler.

This is the function to be reimplemented when a different solver driver needs to ask elements for additional matrices / vectors.

Returns false on errors.

Reimplemented in GmpFemNLSolver, GmpFemNonLinearSolver, and GmpFemTransientSolver.

◆ initResultAttributes()

bool GmpFemSolver::initResultAttributes ( QString  prefix)
protectedvirtual

Helpper function used to register the set of result attributes managed by the fem solver. To avoid name clashes, all of them should prepend the given prefix (usually "fem") to the attribute name.

This is the function to be reimplemented when a different solver driver needs to add other attributes.

Returns false on errors.

◆ initSolver()

bool GmpFemSolver::initSolver ( GmpFemAssembler::FixedDofMode  assemblerMode,
bool  assemblerReverseMapping,
bool  enableFastUpdate 
)
protected

Basic implementation of the init function receiving as parameters the configuration options that are sent to the assembler.

If _solverOptions._assemblerDofMode is equal to automatic (the default), this function will update that option to reflect the assemblerMode received as parameter. If not, the mode parameter passed to the function will be ignored. Either way, _solverOptions._assemblerDofMode will contain the value passed to the assembler.

Parameters
assemblerModeThe dof handling mode used by the assembler. This option is used only if _solverOptions._assemblerDofMode == automatic (the default). If the user has given us a specific assembler option, it will be honored by this function. If you really need to override the user setting, consider changing the value in _solverOptions._assemblerDofMode before calling this function.

\parm assemblerReverseMapping Do we need to enable reverse mapping in the assembler?

Parameters
enableFastUpdateDo we need to enable the fast update option in the assembler? Set it to true only if there will be dof changes in the mesh while solving the model.

◆ meshChanged

void GmpFemSolver::meshChanged ( )
protectedvirtualslot

Slot called when the underlying mesh has been changed.

This default implementation just calls the virtual functions cleanup() and init().

◆ prepareMatrices()

GmpFemPhysics::FemResultType GmpFemSolver::prepareMatrices ( bool  skipFixedBcs = false)
protected

Auxilliary function used to fill the K, and f matrices / vectors along with any other matrices/vectors stored in _matSet and _vecSet.

Returns FEM_RESULT_OK on success, FEM_ERROR if the solver process should be aborted or any other enumeration result if the current prepare matrix was unsuccessfull but the solver process might be able to continue with changes (in the time step, for example).

◆ printElementData()

void GmpFemSolver::printElementData ( const GmElement e,
int  ndof,
const int *  dofMapping,
const GmpFemPhysics p,
const GmBoundaryCondition bc = NULL,
const GmContactBoundaryCondition cbc = NULL,
int  loadId1 = -1,
int  id2 = -1 
)
protected

Print element dof mapping, local matrices and vectors depending on the current solver print options.

Assumes that local matrices and vectors are filled in _matSet and _vecSet. If printing external load data, ignores _matSet and prints only _vecSet values.

  • If bc and cbc are NULL, loadId1 and id2 are -1, we are printing basic element data.
  • If bc and cbc are NULL, but loadId1 is >= 0, we are printing external loads data and id2 should be -1.
  • If bc is not NULL, cbc should be NULL, loadId1 and id2 should be -1 and we are printing boundary condition element data.
  • If cbc is not NULL, bc should be NULL, loadId1 and id2 should be >= 0 and we are printing contact boundary condition data.

◆ setMatrixCombinerObject()

void GmpFemSolver::setMatrixCombinerObject ( const GmpFemAssemblerMatrixCombiner combiner)
inline

Sets the combiner object that will be used by the assembler to merge elemental data from multiple matrices (eg. C and K) into an equivalent matrix.

A combiner object MUST be set before preparing matrices whenever the mode for the matrix set includes an "equivalent" matrix. The solver does NOT take ownership of the combiner, so it can be a reference to a local object, and it must be destroyed by te derived class making use of this function.

◆ setVectorCombinerObject()

void GmpFemSolver::setVectorCombinerObject ( const GmpFemAssemblerVectorCombiner combiner)
inline

Sets the combiner objects that will be used by the assembler to merge elemental data from multiple vectors (eg. Fe and Fi) into an equivalent vector.

A combiner object MUST be set before preparing vectors whenever the mode for the vector set include an "equivalent" vector. The solver does NOT take ownership of the combiner, so it can be a reference to a local object, and it must be destroyed by te derived class making use of this function.

◆ solveLinearSystem()

bool GmpFemSolver::solveLinearSystem ( bool  xFilled)
protected

Solves the linear system K.x = f taking K and f from the single equivalent matrix/vector or from matrix K / vector f, depending on the assembler configuration. Results are saved in _x. Prints the used matrices / result as configured in _printOptions.

If xFilled is true, tells the numeric solver that the current values of _x can be used as an initial guess to solve the linear system.

◆ traverseElements()

GmpFemPhysics::FemResultType GmpFemSolver::traverseElements ( )
protected

Fills the stiffness matrix _K and the force vector _f by traversing elements asking physics for the local element matrix / force vector.

Construtor

The function that will be called for each parallel task

Returns the reported error code in case of an aborted thread

< The Fem solver object, our distinct father

< The physics index inside _self

< Element processing result;

◆ traverseElementsForSaving()

GmpFemPhysics::FemResultType GmpFemSolver::traverseElementsForSaving ( FILE *  f,
int  iter 
)
protected

Similar to traverseElements() but instead of adding the elements to the assembler, saves the element data to the given file.

File format: Iteration d (the argument is the iter parameter value) Physics d s (the arguments are the physics index and its name) ElemId, ElemType, NumDof, DofMap, Sv, Fe, Fi, Ke, Ce, Me Data values for physics 1 (active elements only) ... Physics d s ElemId, ElemType, NumDof, DofMap, Sv, Fe, Fi, Ke, Ce, Me Data values for physics n

  • The iteration, physics and the header lines can be disabled with control options. By default, only the header line is enabled.
  • The set of columns can be changed with control options. Only the ElemId, ElemType and NumDof columns are fixed. The default option includes the Sv and the Fi columns.
  • Even if a column is enabled in the options, it will be skipped if the contents is NOT listed as part of the solver vector set / matrix set.
  • The DofMap column is a string surrounded by "" with the node local number followed by the dof state var name (ex: "1T 2T 3T 4T"). It has NumDof entries. If the state var is multi dimensional, the dimension is given inside [] (like "1u[1] 1u[2] 2u[1] 2u[2] 3u[1] 3u[2]"). If a dof is fixed (fixed B.C), its name will be surrounded by () (like "(1T) (2T) 3T 4T", where the first two dofs are fixed and their values come from the boundary conditions).
  • The number of entries for the Sv, Fe and Fi columns is equal to NumDof. The number of entries for the K, C and M columns is equal to NumDof x NumDof. Values are given in COLUMN MAJOR format.
  • Remember that the mesh can be HETEROGENEOUS, so each line can have a DIFFERENT number of columns. That can also happend if the physics uses a per element dof mapping (like the XFEM)
  • If a specified column exists in the matrix/vector set but was NOT FILLED by the physics, its entries will be replaced by "NA" strings (this maintains the number of line columns)
  • The Sv column (node Dof) values are currently either taken from the _x vector (and not from the state variables) or from the assembler list of fixed dof values if that is the case.
  • Values are formatted according to the print options format definitions (should we have an independent format?)

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