G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsElasticityFunctions.h
Go to the documentation of this file.
1 
16 #pragma once
17 
18 #include <gsCore/gsMultiPatch.h>
20 #include <gsIO/gsOptionList.h>
21 
22 namespace gismo
23 {
24 
28 template <class T>
30 {
31 public:
32 
34  const gsOptionList & options,
35  const gsMultiPatch<T> * geometry,
36  const gsMultiPatch<T> * displacement,
37  const gsMultiPatch<T> * pressure = nullptr,
38  const gsMultiPatch<T> * velocity = nullptr)
39  : m_geometry(geometry),
40  m_displacement(displacement),
41  m_pressure(pressure),
42  m_velocity(velocity),
43  m_patch(patch),
44  m_dim(m_geometry->patch(m_patch).parDim()),
45  m_options(options),
46  m_type(comp)
47  {}
48 
49 
50 
51  virtual short_t domainDim() const
52  {
53  return m_geometry->patch(m_patch).parDim();
54  }
55 
56  virtual short_t targetDim() const
57  {
58  switch(m_type)
59  {
60  case stress_components::von_mises: return 1;
61  case stress_components::all_2D_vector: return 3;
62  case stress_components::all_2D_matrix: return 2;
65  case stress_components::all_3D_matrix: return 3;
66  default: return 0;
67  };
68  }
69 
75  virtual void eval_into(const gsMatrix<T> & u, gsMatrix<T> & result) const
76  {
77  switch (material_law::law(m_options.getInt("MaterialLaw")))
78  {
79  case material_law::hooke : linearElastic(u,result); return;
80  case material_law::saint_venant_kirchhoff : nonLinearElastic(u,result); return;
81  case material_law::neo_hooke_ln : nonLinearElastic(u,result); return;
82  case material_law::neo_hooke_quad : nonLinearElastic(u,result); return;
83  case material_law::mixed_hooke : mixedLinearElastic(u,result); return;
84  case material_law::mixed_neo_hooke_ln : mixedNonLinearElastic(u,result); return;
85  case material_law::muscle : mixedNonLinearElastic(u,result); return;
86  default: return;
87  }
88  }
89 
90 protected:
91 
93  index_t outputCols(index_t inputCols) const
94  {
95  switch (m_type)
96  {
97  case stress_components::von_mises: return inputCols;
98  case stress_components::all_2D_vector: return inputCols;
99  case stress_components::all_2D_matrix: return 2*inputCols;
100  case stress_components::normal_3D_vector: return inputCols;
101  case stress_components::shear_3D_vector: return inputCols;
102  case stress_components::all_3D_matrix: return 3*inputCols;
103  default: return 0;
104  }
105  }
107  void saveStress(const gsMatrix<T> & S, gsMatrix<T> & result, index_t q) const;
108 
110  void linearElastic(const gsMatrix<T> & u, gsMatrix<T> & result) const;
111  void nonLinearElastic(const gsMatrix<T> & u, gsMatrix<T> & result) const;
112  void mixedLinearElastic(const gsMatrix<T> & u, gsMatrix<T> & result) const;
113  void mixedNonLinearElastic(const gsMatrix<T> & u, gsMatrix<T> & result) const;
114 
115 protected:
116  const gsMultiPatch<T> * m_geometry;
117  const gsMultiPatch<T> * m_displacement;
118  const gsMultiPatch<T> * m_pressure;
119  const gsMultiPatch<T> * m_velocity;
120  index_t m_patch;
121  short_t m_dim;
122  const gsOptionList & m_options;
124 
125 }; // class definition ends
126 
127 
128 
132 template <class T>
133 class gsDetFunction : public gsFunction<T>
134 {
135 public:
136 
137  gsDetFunction(const gsMultiPatch<T> & geo, index_t patch)
138  : m_geo(geo),
139  m_patch(patch)
140  {}
141 
142  virtual short_t domainDim() const { return m_geo.domainDim(); }
143 
144  virtual short_t targetDim() const { return 1; }
145 
149  virtual void eval_into(const gsMatrix<T> & u, gsMatrix<T> & result) const;
150 
151 protected:
152 
153  gsMultiPatch<T> const & m_geo;
154  index_t m_patch;
155 
156 }; // class definition ends
157 
162 template <class T>
163 class gsFsiLoad : public gsFunction<T>
164 {
165 public:
166 
167  gsFsiLoad(const gsMultiPatch<T> & geoRef, const gsMultiPatch<T> & ALEdisplacement,
168  index_t patchGeo, boxSide sideGeo,
169  const gsMultiPatch<T> & velocity, const gsMultiPatch<T> & pressure,
170  index_t patchVelPres, T viscosity, T density)
171  : m_geo(geoRef),
172  m_ale(ALEdisplacement),
173  m_patchGeo(patchGeo),
174  m_sideGeo(sideGeo),
175  m_vel(velocity),
176  m_pres(pressure),
177  m_patchVP(patchVelPres),
178  m_viscosity(viscosity),
179  m_density(density)
180  {}
181 
182  virtual short_t domainDim() const { return m_geo.domainDim(); }
183 
184  virtual short_t targetDim() const { return m_geo.domainDim(); }
185 
189  virtual void eval_into(const gsMatrix<T> & u, gsMatrix<T> & result) const;
190 
191 protected:
192 
193  gsMultiPatch<T> const & m_geo;
194  gsMultiPatch<T> const & m_ale;
195  index_t m_patchGeo;
196  boxSide m_sideGeo;
197  gsMultiPatch<T> const & m_vel;
198  gsMultiPatch<T> const & m_pres;
199  index_t m_patchVP;
200  T m_viscosity;
201  T m_density;
202 
203 }; // class definition ends
204 
205 } // namespace ends
206 
207 
208 #ifndef GISMO_BUILD_LIB
209 #include GISMO_HPP_HEADER(gsElasticityFunctions.hpp)
210 #endif
Compute Cauchy stresses for a previously computed/defined displacement field. Can be pushed into gsPi...
Definition: gsElasticityFunctions.h:29
#define short_t
Definition: gsConfig.h:35
S = lambda*ln(J)*C^-1 + mu*(I-C^-1)
Definition: gsBaseUtils.h:134
virtual short_t domainDim() const
Dimension of the (source) domain.
Definition: gsElasticityFunctions.h:51
#define index_t
Definition: gsConfig.h:32
A function from a n-dimensional domain to an m-dimensional image.
Definition: gsFunction.h:59
const index_t & getInt(const std::string &label) const
Reads value for option label from options.
Definition: gsOptionList.cpp:37
sigma = 2*mu*eps + p*I
Definition: gsBaseUtils.h:131
return von Mises stress as a scala
Definition: gsBaseUtils.h:114
sigma = 2*mu*eps + lambda*tr(eps)*I
Definition: gsBaseUtils.h:132
Provides several simple utility and naming classes.
Provides a list of labeled parameters/options that can be set and accessed easily.
Provides declaration of the MultiPatch class.
S = p*C^-1 + mu*(I-C^-1)
Definition: gsBaseUtils.h:136
S = lambda/2*(J^2-1)C^-1 + mu*(I-C^-1)
Definition: gsBaseUtils.h:135
Loading function to transfer fluid action to the solid. Used in Fluid-Structure Interaction simulatio...
Definition: gsElasticityFunctions.h:163
virtual void eval_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Each column of the input matrix (u) corresponds to one evaluation point. Columns of the output matrix...
Definition: gsElasticityFunctions.h:75
void linearElastic(const gsMatrix< T > &u, gsMatrix< T > &result) const
computation routines for different material laws
Definition: gsElasticityFunctions.hpp:26
virtual short_t domainDim() const
Dimension of the (source) domain.
Definition: gsElasticityFunctions.h:182
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition: gsMultiPatch.h:33
Struct which represents a certain side of a box.
Definition: gsBoundary.h:84
law
Definition: gsBaseUtils.h:128
virtual short_t targetDim() const
Dimension of the target space.
Definition: gsElasticityFunctions.h:144
virtual void eval_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Each column of the input matrix (u) corresponds to one evaluation point. Each column of the output ma...
Definition: gsElasticityFunctions.hpp:227
void saveStress(const gsMatrix< T > &S, gsMatrix< T > &result, index_t q) const
save components of the stress tensor to the output matrix according to the m_type ...
Definition: gsElasticityFunctions.hpp:204
virtual void eval_into(const gsMatrix< T > &u, gsMatrix< T > &result) const
Each column of the input matrix (u) corresponds to one evaluation point. Each column of the output ma...
Definition: gsElasticityFunctions.hpp:241
Compute jacobian determinant of the geometry mapping. Can be pushed into gsPiecewiseFunction to const...
Definition: gsElasticityFunctions.h:133
virtual short_t targetDim() const
Dimension of the target space.
Definition: gsElasticityFunctions.h:184
index_t outputCols(index_t inputCols) const
size of the output matrix according to the m_type
Definition: gsElasticityFunctions.h:93
return all components of the 2D stress tensor as a 2x2 matrix
Definition: gsBaseUtils.h:117
virtual short_t domainDim() const
Dimension of the (source) domain.
Definition: gsElasticityFunctions.h:142
Class which holds a list of parameters/options, and provides easy access to them. ...
Definition: gsOptionList.h:32
Definition: gsBaseUtils.h:119
S = 2*mu*E + lambda*tr(E)*I.
Definition: gsBaseUtils.h:133
Definition: gsBaseUtils.h:116
components
Definition: gsBaseUtils.h:111
virtual short_t targetDim() const
Dimension of the target space.
Definition: gsElasticityFunctions.h:56
Definition: gsBaseUtils.h:121