![]() |
MFEM v4.8.0
Finite element discretization library
|
#include <nonlininteg.hpp>
Public Member Functions | |
NeoHookeanModel (real_t mu_, real_t K_, real_t g_=1.0) | |
NeoHookeanModel (Coefficient &mu_, Coefficient &K_, Coefficient *g_=NULL) | |
real_t | EvalW (const DenseMatrix &J) const override |
Evaluate the strain energy density function, W = W(Jpt). | |
void | EvalP (const DenseMatrix &J, DenseMatrix &P) const override |
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt). | |
void | AssembleH (const DenseMatrix &J, const DenseMatrix &DS, const real_t weight, DenseMatrix &A) const override |
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the local gradient matrix 'A'. | |
![]() | |
HyperelasticModel () | |
virtual | ~HyperelasticModel () |
void | SetTransformation (ElementTransformation &Ttr_) |
Protected Member Functions | |
void | EvalCoeffs () const |
Protected Attributes | |
real_t | mu |
real_t | K |
real_t | g |
Coefficient * | c_mu |
Coefficient * | c_K |
Coefficient * | c_g |
bool | have_coeffs |
DenseMatrix | Z |
DenseMatrix | G |
DenseMatrix | C |
![]() | |
ElementTransformation * | Ttr |
Neo-Hookean hyperelastic model with a strain energy density function given by the formula:
. The parameters
Definition at line 269 of file nonlininteg.hpp.
Definition at line 282 of file nonlininteg.hpp.
|
inline |
Definition at line 285 of file nonlininteg.hpp.
|
overridevirtual |
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the local gradient matrix 'A'.
[in] | Jpt | Represents the target->physical transformation Jacobian matrix. |
[in] | DS | Gradient of the basis matrix (dof x dim). |
[in] | weight | Quadrature weight coefficient for the point. |
[in,out] | A | Local gradient matrix where the contribution from this point will be added. |
Computes weight * d(dW_dxi)_d(xj) at the current point, for all i and j, where x1 ... xn are the FE dofs. This function is usually defined using the matrix invariants and their derivatives.
Implements mfem::HyperelasticModel.
Definition at line 318 of file nonlininteg.cpp.
|
inlineprotected |
Definition at line 271 of file nonlininteg.cpp.
|
overridevirtual |
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
[in] | Jpt | Represents the target->physical transformation Jacobian matrix. |
[out] | P | The evaluated 1st Piola-Kirchhoff stress tensor. |
Implements mfem::HyperelasticModel.
Definition at line 297 of file nonlininteg.cpp.
|
overridevirtual |
Evaluate the strain energy density function, W = W(Jpt).
[in] | Jpt | Represents the target->physical transformation Jacobian matrix. |
Implements mfem::HyperelasticModel.
Definition at line 281 of file nonlininteg.cpp.
|
protected |
Definition at line 277 of file nonlininteg.hpp.
|
protected |
Definition at line 273 of file nonlininteg.hpp.
|
protected |
Definition at line 273 of file nonlininteg.hpp.
|
protected |
Definition at line 273 of file nonlininteg.hpp.
|
mutableprotected |
Definition at line 277 of file nonlininteg.hpp.
|
protected |
Definition at line 272 of file nonlininteg.hpp.
|
protected |
Definition at line 274 of file nonlininteg.hpp.
|
protected |
Definition at line 272 of file nonlininteg.hpp.
|
mutableprotected |
Definition at line 272 of file nonlininteg.hpp.
|
mutableprotected |
Definition at line 276 of file nonlininteg.hpp.