FemProcess
The GeMA Fem Process Plugin
gmpFemPhysics.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2014 by Carlos Augusto Teixera Mendes
4 ** All rights reserved.
5 **
6 ** This file is part of the "GeMA" software. It's use should respect
7 ** the terms in the license agreement that can be found together
8 ** with this source code.
9 ** It is provided AS IS, with NO WARRANTY OF ANY KIND,
10 ** INCLUDING THE WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR
11 ** A PARTICULAR PURPOSE.
12 **
13 ************************************************************************/
14 
24 #ifndef _GEMA_PLUGIN_FEM_PHYSICS_H_
25 #define _GEMA_PLUGIN_FEM_PHYSICS_H_
26 
27 #include <gmPhysics.h>
28 #include <gmElementMesh.h>
29 #include <gmElement.h>
30 
31 #include "gmpFemProcessConfig.h"
32 #include "gmpFemMatrixSet.h"
33 #include "gmpFemVectorSet.h"
34 
35 class GmElementDof;
36 class GmMatrixDof;
39 class GmDiscontinuitySet;
40 class GmGaussAccessor;
41 class GmStateDump;
42 
44 class GMP_FEM_PROCESS_API_EXPORT GmpFemPhysics : public GmPhysics
45 {
46 public:
47 
50  {
54  };
55 
58  {
59  FEM_PARALLEL_FILL_ELEMENT_DATA = 0x01,
60  };
61 
62  GmpFemPhysics(GmSimulationData* simulation, QString id, QString description, const GmLogCategory& logger);
63  virtual ~GmpFemPhysics();
64 
66  GmElementMesh* mesh() const { return _mesh; }
67 
75  const GmCellGroupSet* meshGroupSet() const { return _gs; }
76 
80  int numSupportedExternalLoads() const { return _externalLoadGroupSets.size(); }
81 
90  const GmCellGroupSet* externalLoadGroupSet(int loadId) const { return _externalLoadGroupSets[loadId]; }
91 
94 
96  const QMap<QString, const GmContactBoundaryCondition*>& contactBoundaryConditions() const { return _cbcMap; }
97 
99  GmDiscontinuitySet* discontinuitySet() const { return _discSet; }
100 
108  virtual bool checkPhysicsDependencies(const QList<GmpFemPhysics*>& physicsList, int index)
109  {
110  Q_UNUSED(physicsList);
111  Q_UNUSED(index);
112  return true;
113  }
114 
126  virtual bool dofByElement(bool* fixed = NULL, bool* addOnly = NULL, bool* trackChanges = NULL) const
127  {
128  Q_UNUSED(fixed); Q_UNUSED(addOnly); Q_UNUSED(trackChanges);
129  return false;
130  }
131 
155  virtual const GmElementDof* dofMapping(GmCellType type) const = 0;
156 
176  virtual const GmElementDof* dofMapping(const GmElement* e) const { return dofMapping(e->type()); }
177 
182  virtual Unit dofUnit(int dof) const = 0;
183 
185  virtual Unit timeUnit() const = 0;
186 
195  {
196  static QList<QPair<int, GmValueAccessor*> > dummyList;
197  return dummyList;
198  }
199 
207  virtual bool supportsCellType(GmCellType type) const = 0;
208 
213  virtual bool supportsParallel(FemSupportedParallelMethods methodType) const { Q_UNUSED(methodType); return false; }
214 
235  virtual bool beforeElementStiffnessLoop(const GmpFemMatrixSet& elemMatrices,
236  const GmpFemVectorSet& elemVectors)
237  {
238  Q_UNUSED(elemMatrices); Q_UNUSED(elemVectors);
239  return true;
240  }
241 
245  virtual void afterElementStiffnessLoop() {}
246 
294  virtual FemResultType fillElementData(const GmElement* e, GmpFemMatrixSet& elemMatrices, GmpFemVectorSet& elemVectors) = 0;
295 
300  virtual const QStringList& supportedBcTypes() const
301  {
302  static QStringList dummyList;
303  return dummyList;
304  }
305 
310  virtual bool bcAffectsCalc(const GmBoundaryCondition* bc) const
311  {
312  Q_UNUSED(bc);
313  return false;
314  }
315 
333  virtual FemResultType fillElementDataForBc(const GmElement* e, const GmBoundaryCondition* bc, int bcIndex, int bcListIndex,
334  int border, GmpFemMatrixSet& elemMatrices, GmpFemVectorSet& elemVectors)
335  {
336  Q_UNUSED(e); Q_UNUSED(bc); Q_UNUSED(bcIndex); Q_UNUSED(bcListIndex); Q_UNUSED(border);
337  Q_UNUSED(elemMatrices); Q_UNUSED(elemVectors);
338  return FEM_ERROR;
339  }
340 
358  virtual bool fixedNodalForcesBc(QList<int>& nodes, QList<int>&dof, QList<double>& values) const
359  {
360  Q_UNUSED(nodes); Q_UNUSED(dof); Q_UNUSED(values);
361  return true;
362  }
363 
381  virtual bool fixedNodalDofsBc(QList<int>& nodes, QList<int>&dof, QList<double>& values, bool* constantValues) const
382  {
383  Q_UNUSED(nodes); Q_UNUSED(dof); Q_UNUSED(values);
384  *constantValues = true;
385  return true;
386  }
387 
392  virtual const QStringList& supportedContactBcTypes() const
393  {
394  static QStringList dummyList;
395  return dummyList;
396  }
397 
430  virtual FemResultType fillContactData(const GmContactBoundaryCondition* cbc, int bcIndex1, int bcIndex2,
431  GmMatrixDof& dofMap, GmpFemMatrixSet& matrixSet, GmpFemVectorSet& vecSet)
432  {
433  Q_UNUSED(cbc); Q_UNUSED(bcIndex1); Q_UNUSED(bcIndex2);
434  Q_UNUSED(dofMap); Q_UNUSED(matrixSet); Q_UNUSED(vecSet);
435  return FEM_ERROR;
436  }
437 
444  virtual const QStringList& supportedExternalLoads() const
445  {
446  static QStringList dummyList;
447  return dummyList;
448  }
449 
462  virtual FemResultType fillElementDataForLoads(const GmElement* e, int loadId, GmpFemVectorSet& elemVectors)
463  {
464  Q_UNUSED(e); Q_UNUSED(loadId); Q_UNUSED(elemVectors);
465  return FEM_ERROR;
466  }
467 
480  virtual bool calcDerivedResults(bool nonLinearSolver)
481  {
482  Q_UNUSED(nonLinearSolver);
483  return true;
484  }
485 
491  virtual bool supportsStateDumping() { return false; }
492 
498  virtual bool addStateItemsToGroup(GmStateDump* state, bool fixedHint, int groupId) { Q_UNUSED(state); Q_UNUSED(fixedHint); Q_UNUSED(groupId); return true; }
499 
504  virtual bool fillStateControlMapData(QVariantMap* map) { Q_UNUSED(map); return true; }
505 
512  virtual bool stateControlMapDataLoaded(QVariantMap* map) { Q_UNUSED(map); return true; }
513 
517  virtual bool stateAboutToBeSaved(GmStateDump* state) { Q_UNUSED(state); return true; }
518 
522  virtual bool stateSaved(GmStateDump* state) { Q_UNUSED(state); return true; }
523 
527  virtual bool stateAboutToBeLoaded(GmStateDump* state) { Q_UNUSED(state); return true; }
528 
532  virtual bool stateLoaded(GmStateDump* state) { Q_UNUSED(state); return true; }
533 
534 
542  virtual QList<int> changedElements() const { return QList<int>(); }
543 
547  virtual void clearChangedElements() {}
548 
549  const GmIntegrationRule* elementIntegrationRule(GmCellType type) const;
550  const GmBorderIntegrationRule* borderIntegrationRule (GmCellType type, int borderIndex) const;
551  const GmBorderIntegrationRule* faceIntegrationRule (GmCellType type, int faceIndex) const;
552  const GmBorderIntegrationRule* edgeIntegrationRule (GmCellType type) const;
553 
555  int integrationRuleSet() const { return _intRuleSet; }
556 
557  bool setIntegrationRuleSet(int ruleSet);
558 
559  virtual bool loadPrivateData(LuaTable& table);
560  virtual void printParameters(const GmLogCategory& logger);
561 
562 protected:
563 
572  virtual bool checkMeshDimension(int ndim) = 0;
573 
580  virtual bool checkAndLoadPrivateData(LuaTable& table) = 0;
581 
593  virtual bool checkAndLoadDofMapping(LuaTable& table) = 0;
594 
609  virtual bool checkAndLoadPropertyAccessors(LuaTable& table) = 0;
610 
625  virtual bool checkAndLoadAttributeAccessors(LuaTable& nodeTable, LuaTable& gaussTable) = 0;
626 
640  virtual bool checkAndLoadExternalLoadAccessors(LuaTable& loadsTable) = 0;
641 
655  virtual bool checkAndLoadBcAccessors() = 0;
656 
670  virtual bool checkAndLoadContactBcAccessors() = 0;
671 
687  virtual bool checkAndLoadDiscontinuitySetAccessors(LuaTable& table) = 0;
688 
689  GmGaussAccessor* gaussAttributeAccessor(GmElementMesh* mesh, QString stdAttributeName, QString msgDescription,
690  LuaTable& table, Unit desiredUnit, bool required, int type = -1,
691  int nlin = -1, int ncol = -1, int history = -1, bool canBeFunction = true,
692  bool create = false, QString createFormat = "", int createHistory = -1,
693  bool ignoreWarnings = false);
694 
699  virtual GmGaussAccessor* createGaussAttributeAccessor(GmElementMesh* mesh, QString id, int snum, bool locked,
700  Unit desiredUnit, const GmLogCategory& logger) const
701  {
702  return mesh->gaussAttributeAccessor(id, snum, locked, desiredUnit, logger);
703  }
704 
705 
712 
718 };
719 
720 
721 #endif
722 
virtual bool stateControlMapDataLoaded(QVariantMap *map)
Virtual method called just after the solver control map was loaded from the state....
Definition: gmpFemPhysics.h:512
const GmCellGroupSet * externalLoadGroupSet(int loadId) const
Returns the element group set used to traverse mesh elements tied to the given external loads id.
Definition: gmpFemPhysics.h:90
virtual bool stateSaved(GmStateDump *state)
Virtual method called just after succesfully completing a save operation on the given state....
Definition: gmpFemPhysics.h:522
virtual bool fillStateControlMapData(QVariantMap *map)
Virtual method called just before the solver control map is saved to the state. Should be used to fil...
Definition: gmpFemPhysics.h:504
int numSupportedExternalLoads() const
Returns the number of supported external load types by this physics. Based on the list returned by su...
Definition: gmpFemPhysics.h:80
virtual bool loadPrivateData(LuaTable &table)=0
virtual FemResultType fillElementDataForLoads(const GmElement *e, int loadId, GmpFemVectorSet &elemVectors)
Method called by the FEM process in a similar way as the call to fillElementData() for the physics to...
Definition: gmpFemPhysics.h:462
virtual const QStringList & supportedBcTypes() const
A list with names for the known boundary condition types supported by this physics.
Definition: gmpFemPhysics.h:300
The method finished correctly.
Definition: gmpFemPhysics.h:51
FemResultType
Result type for local matrix calculation methods.
Definition: gmpFemPhysics.h:49
virtual GmGaussAccessor * createGaussAttributeAccessor(GmElementMesh *mesh, QString id, int snum, bool locked, Unit desiredUnit, const GmLogCategory &logger) const
A virtual function called by gaussAttributeAccessor() to get the returned Gauss accessor....
Definition: gmpFemPhysics.h:699
virtual bool addStateItemsToGroup(GmStateDump *state, bool fixedHint, int groupId)
Initialization method, called once, allowing the object to add its state items to the given group of ...
Definition: gmpFemPhysics.h:498
virtual const GmElementDof * dofMapping(const GmElement *e) const
Returns the mapping of node degrees of freedom for the given element or NULL if this physics object d...
Definition: gmpFemPhysics.h:176
virtual void printParameters(const GmLogCategory &logger)
virtual FemResultType fillElementDataForBc(const GmElement *e, const GmBoundaryCondition *bc, int bcIndex, int bcListIndex, int border, GmpFemMatrixSet &elemMatrices, GmpFemVectorSet &elemVectors)
Method called by the FEM process in a similar way as the call to fillElementData() for the physics to...
Definition: gmpFemPhysics.h:333
A convenience class that builds a GmVectorSet with the types given by GmpFemVectorTypes and also asso...
Definition: gmpFemVectorSet.h:47
A convenience class that builds a GmMatrixSet with the types given by GmpFemMatrixTypes and also asso...
Definition: gmpFemMatrixSet.h:56
Declaration of the GmpFemVectorSet class.
const QMap< QString, const GmBoundaryCondition * > & boundaryConditions() const
Returns a reference to the map storing contact boundary conditions associated to this physics,...
Definition: gmpFemPhysics.h:93
virtual bool fixedNodalDofsBc(QList< int > &nodes, QList< int > &dof, QList< double > &values, bool *constantValues) const
Method called by the FEM process asking the physics to provide fixed dof (state variable) values used...
Definition: gmpFemPhysics.h:381
virtual bool supportsParallel(FemSupportedParallelMethods methodType) const
Returns true if this physics object supports multiple calls in parallel (by multiple threads) to the ...
Definition: gmpFemPhysics.h:213
QMap< QString, const GmBoundaryCondition * > _bcMap
A map storing boundary conditions associated to this physics, keyed by bc type.
Definition: gmpFemPhysics.h:710
virtual bool stateAboutToBeLoaded(GmStateDump *state)
Virtual method called just before starting a load operation on the given state. Returning false abort...
Definition: gmpFemPhysics.h:527
FemSupportedParallelMethods
Enum used to specify which physics plugin methods can be called in parallel. Values can be ored toget...
Definition: gmpFemPhysics.h:57
virtual const QList< QPair< int, GmValueAccessor * > > & dofForceAttributes() const
Returns a list that associates degrees of freedom handled by this physics with their associated force...
Definition: gmpFemPhysics.h:194
GmDiscontinuitySet * discontinuitySet() const
Returns the associated discontinuity set, NULL if there is none.
Definition: gmpFemPhysics.h:99
virtual GmGaussAccessor * gaussAttributeAccessor(QString id, Unit desiredUnit, const GmLogCategory &logger) const
virtual const QStringList & supportedExternalLoads() const
A list with names for the known external load conditions supported by this physics.
Definition: gmpFemPhysics.h:444
virtual FemResultType fillContactData(const GmContactBoundaryCondition *cbc, int bcIndex1, int bcIndex2, GmMatrixDof &dofMap, GmpFemMatrixSet &matrixSet, GmpFemVectorSet &vecSet)
Method called by the FEM process in a similar way as the call to fillElementData() for the physics to...
Definition: gmpFemPhysics.h:430
The method finished due to an unexpected error.
Definition: gmpFemPhysics.h:53
virtual bool dofByElement(bool *fixed=NULL, bool *addOnly=NULL, bool *trackChanges=NULL) const
Returns false if degrees of freedom are mapped to nodes depending only on the element type,...
Definition: gmpFemPhysics.h:126
virtual bool stateAboutToBeSaved(GmStateDump *state)
Virtual method called just before starting a save operation on the given state. Returning false abort...
Definition: gmpFemPhysics.h:517
virtual QList< int > changedElements() const
Method called by the assembler to obtain the list of elements that had their dofs changed since the p...
Definition: gmpFemPhysics.h:542
QMap< QString, const GmContactBoundaryCondition * > _cbcMap
A map storing contact boundary conditions associated to this physics, keyed by bc type.
Definition: gmpFemPhysics.h:711
int integrationRuleSet() const
Returns the integration rule set in use by this physics object.
Definition: gmpFemPhysics.h:555
virtual GmCellType type() const=0
virtual bool checkPhysicsDependencies(const QList< GmpFemPhysics * > &physicsList, int index)
A function called by the Fem solver, at its initialization, to let a physics check that it is compati...
Definition: gmpFemPhysics.h:108
GmCellType
QVector< const GmCellGroupSet * > _externalLoadGroupSets
Element group sets for each active external loads or NULL for inactive ones. Indexed in the same orde...
Definition: gmpFemPhysics.h:717
virtual bool fixedNodalForcesBc(QList< int > &nodes, QList< int > &dof, QList< double > &values) const
Method called by the FEM process asking the physics to provide prescribed nodal force values that wil...
Definition: gmpFemPhysics.h:358
const GmCellGroupSet * meshGroupSet() const
Returns the element group set used to traverse mesh elements belonging to the element groups that tog...
Definition: gmpFemPhysics.h:75
virtual void clearChangedElements()
Method called by the assembler to clear the internal list used to track dof changes....
Definition: gmpFemPhysics.h:547
Declaration of the GmpFemMatrixSet class.
const GmCellGroupSet * _gs
The element group set used to traverse the mesh elements that form the physics spatial domain.
Definition: gmpFemPhysics.h:707
virtual void afterElementStiffnessLoop()
Method called by the Fem process to notify the physical object that the assembly process has ended.
Definition: gmpFemPhysics.h:245
virtual bool supportsStateDumping()
Method called by the solver to check if this physics supports state dumping and restoring....
Definition: gmpFemPhysics.h:491
The method finished by unexpected error assembling Fi.
Definition: gmpFemPhysics.h:52
virtual const QStringList & supportedContactBcTypes() const
A list with names for the known contact boundary condition types supported by this physics.
Definition: gmpFemPhysics.h:392
Base interface class for FEM Physics type plugins.
Definition: gmpFemPhysics.h:44
virtual bool stateLoaded(GmStateDump *state)
Virtual method called just after succesfully completing a load operation on the given state....
Definition: gmpFemPhysics.h:532
GmDiscontinuitySet * _discSet
The discontinuity set, NULL if none was associated.
Definition: gmpFemPhysics.h:709
virtual bool calcDerivedResults(bool nonLinearSolver)
Method called after the the system solution asking each physics to calculate derived results based on...
Definition: gmpFemPhysics.h:480
virtual bool bcAffectsCalc(const GmBoundaryCondition *bc) const
Method called by the FEM process asking the physics whether the supplied boundary condition affects t...
Definition: gmpFemPhysics.h:310
Declaration of usefull configuration definitions for the plugin library.
GmElementMesh * _mesh
The mesh that this physics is bound to.
Definition: gmpFemPhysics.h:706
int _intRuleSet
The integration rule set used by this physics.
Definition: gmpFemPhysics.h:708
virtual bool beforeElementStiffnessLoop(const GmpFemMatrixSet &elemMatrices, const GmpFemVectorSet &elemVectors)
Method called by the Fem process to notify the physical object that the assembly process will start.
Definition: gmpFemPhysics.h:235
GmElementMesh * mesh() const
Returns the element mesh that this physics object is tied to.
Definition: gmpFemPhysics.h:66