MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMaterialFluid.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 
25 #ifndef _GEMA_PLUGIN_MECHANICALMATERIAL_FLUID_H_
26 #define _GEMA_PLUGIN_MECHANICALMATERIAL_FLUID_H_
27 
28 #include "gmpMechanicalConfig.h"
29 #include <gmpFemPhysics.h>
30 #include <gmMathUtils.h>
31 
32 #include "gmpMechanicalMaterial.h"
33 #include "gmpMechanicPoint.h"
34 
35 class GMP_MECHANICAL_PHYSICS_API_EXPORT GmpMaterialFluid : public GmpMechanicalMaterial
36 {
37 protected:
40  {
42 
43  // --- NO ADDING BELOW THIS LINE
45  };
46 
47 public:
48 
50  GmpMaterialFluid(int typeIndex, QString typeName, const GmLogCategory& logger)
51  : GmpMechanicalMaterial(typeIndex, typeName, logger) {}
52 
54  virtual ~GmpMaterialFluid() {}
55 
56  // Returns a map with material associated properties.
57  virtual const QVariantMap* materialMetaDataMap()
58  {
59  S_TRACE();
60  static QVariantMap m;
61  if (m.isEmpty())
62  {
63  m["properties"] = GmpFemPhysicsCommon::ValueList()
64  << GmpFemPhysicsCommon::ScalarValue(Bf_ID, "Bf", QObject::tr("Fluid compressibility"), "1/kPa", true);
65  }
66  return &m;
67  }
68 
70  static GmpFemPhysicsCommonMaterial* instance(GmSimulationData* simulation, int typeIndex, QString typeName, const GmLogCategory& logger)
71  {
72  Q_UNUSED(simulation);
73  return new GmpMaterialFluid(typeIndex, typeName, logger);
74  }
75 
76  // Returns the stresses according to the material behavior adopted
77  virtual bool mechanicalConstitutiveModel(const GmElement* e, GmMatrix& De, const GmpMechanicPoint* mp, const GmVector* coord, const GmVector& Time, unsigned sc, bool ips) const
78  {
79  S_TRACE();
80  Q_UNUSED(Time);
81  assert(e);
82 
83  //stress components
84  const unsigned ns = GmMathUtils::bit_enum_count(sc, (unsigned)GmpMechanicPoint::stress::last);
85  // Integration point
86  int ip = mp->_index;
87  // Defines scalar varibles
88  GmVector s(ns), dStrain(ns);
89 
90  //old strain
91  const GmVector eu_old(mp->_oldStrain->valueAt(e, ip, coord), ns);
92  // old stress
93  const GmVector s_old(mp->_oldStress->valueAt(e, ip, coord), ns);
94  //residual stress
95  const GmVector sr(mp->_residualStress->valueAt(e, ip, coord), ns);
96 
97  //new strain
98  const GmVector eu(mp->_newStrain->valueAt(e,ip, coord), ns);
99 
100  // Fluid Compressibility
101  double Cf = fluidCompressibility(e, coord, ip);
102  // Bulk modulus
103  double K = 1/Cf;
104 
105  // Fills the elastic matrix;
106  GmpMechanicUtils::elasticStiffnessKG(De, K, 0.0, sc, ips);
107 
108  // Strain increment
109  dStrain = eu - eu_old;
110 
111  //Updates of stresses
112  s = s_old + De*dStrain + sr;
113 
114  //Updates state variables
115  mp->_newStress->setValue(e, ip, s.mem);
116  mp->_newYSR->setScalarValue(e, ip, 1.0);
117  return true;
118  }
119 
121  virtual bool isIsotropic() const { return true; }
122 
124  virtual double fluidCompressibility(const GmElement* e, const GmVector* coord, int ip) const
125  {
126  S_TRACE();
127  double Cf = propertyAc(Bf_ID)->scalarValueAt(e, coord, ip);
128 
129  // Computes Bulk modulus
130  if (Cf == 0.0)
131  {
132  return 1.0e-25;
133  }
134 
135  return Cf;
136  }
137 
138 };
139 
140 #endif
void elasticStiffnessKG(GmMatrix &De, double K, double G, unsigned sc, unsigned ips)
Computes the elastic stiffness matrix using K and G.
Definition: gmpMechanicalMaterial.cpp:983
ElementPropertyIds
IDs for material element properties.
Definition: gmpMaterialFluid.h:39
virtual bool mechanicalConstitutiveModel(const GmElement *e, GmMatrix &De, const GmpMechanicPoint *mp, const GmVector *coord, const GmVector &Time, unsigned sc, bool ips) const
Evaluates stress and tangent matrix according to the material behavior adopted.
Definition: gmpMaterialFluid.h:77
#define S_TRACE()
static GmpFemPhysicsCommonMaterial * instance(GmSimulationData *simulation, int typeIndex, QString typeName, const GmLogCategory &logger)
A "factory" function used to register the material with the physics material factory.
Definition: gmpMaterialFluid.h:70
Declaration of the GmpMechanicalMaterial class.
QString tr(const char *sourceText, const char *disambiguation, int n)
Basic class providing the interface for a mechanical material.
Definition: gmpMechanicalMaterial.h:41
virtual ~GmpMaterialFluid()
Virtual destructor.
Definition: gmpMaterialFluid.h:54
The number of property ids above.
Definition: gmpMaterialFluid.h:44
Definition: gmpMechanicPoint.h:32
GmpMaterialFluid(int typeIndex, QString typeName, const GmLogCategory &logger)
Constructor. Gets as parameters the material index and its name.
Definition: gmpMaterialFluid.h:50
char bit_enum_count(unsigned list, unsigned last)
virtual bool isIsotropic() const
Returns true if the material is isotropic, false otherwise.
Definition: gmpMaterialFluid.h:121
Id for retrieving the Fluid compressibiolity accessor.
Definition: gmpMaterialFluid.h:41
Declaration of the GmpMechanicPoint class.
Definition: gmpMaterialFluid.h:35
Declaration of usefull configuration definitions for the plugin library.
virtual double fluidCompressibility(const GmElement *e, const GmVector *coord, int ip) const
Returns the material elastic modulus.
Definition: gmpMaterialFluid.h:124
arma::vec GmVector
arma::mat GmMatrix