MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
nonlininteg.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_NONLININTEG
13 #define MFEM_NONLININTEG
14 
15 #include "../config/config.hpp"
16 #include "fe.hpp"
17 #include "coefficient.hpp"
18 #include "fespace.hpp"
19 
20 namespace mfem
21 {
22 
23 /** @brief This class is used to express the local action of a general nonlinear
24  finite element operator. In addition it may provide the capability to
25  assemble the local gradient operator and to compute the local energy. */
27 {
28 protected:
30 
32  : IntRule(ir) { }
33 
34 public:
35  /** @brief Prescribe a fixed IntegrationRule to use (when @a ir != NULL) or
36  let the integrator choose (when @a ir == NULL). */
37  void SetIntRule(const IntegrationRule *ir) { IntRule = ir; }
38 
39  /// Prescribe a fixed IntegrationRule to use.
40  void SetIntegrationRule(const IntegrationRule &irule) { IntRule = &irule; }
41 
42  /// Perform the local action of the NonlinearFormIntegrator
43  virtual void AssembleElementVector(const FiniteElement &el,
45  const Vector &elfun, Vector &elvect);
46 
47  /// @brief Perform the local action of the NonlinearFormIntegrator resulting
48  /// from a face integral term.
49  virtual void AssembleFaceVector(const FiniteElement &el1,
50  const FiniteElement &el2,
52  const Vector &elfun, Vector &elvect);
53 
54  /// Assemble the local gradient matrix
55  virtual void AssembleElementGrad(const FiniteElement &el,
57  const Vector &elfun, DenseMatrix &elmat);
58 
59  /// @brief Assemble the local action of the gradient of the
60  /// NonlinearFormIntegrator resulting from a face integral term.
61  virtual void AssembleFaceGrad(const FiniteElement &el1,
62  const FiniteElement &el2,
64  const Vector &elfun, DenseMatrix &elmat);
65 
66  /// Compute the local energy
67  virtual double GetElementEnergy(const FiniteElement &el,
69  const Vector &elfun);
70 
71  /// Method defining partial assembly.
72  /** The result of the partial assembly is stored internally so that it can be
73  used later in the methods AddMultPA(). */
74  virtual void AssemblePA(const FiniteElementSpace &fes);
75 
76  /** The result of the partial assembly is stored internally so that it can be
77  used later in the methods AddMultPA().
78  Used with BilinearFormIntegrators that have different spaces. */
79  virtual void AssemblePA(const FiniteElementSpace &trial_fes,
80  const FiniteElementSpace &test_fes);
81 
82  /// Method for partially assembled action.
83  /** Perform the action of integrator on the input @a x and add the result to
84  the output @a y. Both @a x and @a y are E-vectors, i.e. they represent
85  the element-wise discontinuous version of the FE space.
86 
87  This method can be called only after the method AssemblePA() has been
88  called. */
89  virtual void AddMultPA(const Vector &x, Vector &y) const;
90 
92 };
93 
94 /** The abstract base class BlockNonlinearFormIntegrator is
95  a generalization of the NonlinearFormIntegrator class suitable
96  for block state vectors. */
98 {
99 public:
100  /// Compute the local energy
101  virtual double GetElementEnergy(const Array<const FiniteElement *>&el,
103  const Array<const Vector *>&elfun);
104 
105  /// Perform the local action of the BlockNonlinearFormIntegrator
108  const Array<const Vector *> &elfun,
109  const Array<Vector *> &elvec);
110 
111  virtual void AssembleFaceVector(const Array<const FiniteElement *> &el1,
112  const Array<const FiniteElement *> &el2,
114  const Array<const Vector *> &elfun,
115  const Array<Vector *> &elvect);
116 
117  /// Assemble the local gradient matrix
118  virtual void AssembleElementGrad(const Array<const FiniteElement*> &el,
120  const Array<const Vector *> &elfun,
121  const Array2D<DenseMatrix *> &elmats);
122 
123  virtual void AssembleFaceGrad(const Array<const FiniteElement *>&el1,
126  const Array<const Vector *> &elfun,
127  const Array2D<DenseMatrix *> &elmats);
128 
130 };
131 
132 
133 /// Abstract class for hyperelastic models
135 {
136 protected:
137  ElementTransformation *Ttr; /**< Reference-element to target-element
138  transformation. */
139 
140 public:
141  HyperelasticModel() : Ttr(NULL) { }
142  virtual ~HyperelasticModel() { }
143 
144  /// A reference-element to target-element transformation that can be used to
145  /// evaluate Coefficient%s.
146  /** @note It is assumed that _Ttr.SetIntPoint() is already called for the
147  point of interest. */
148  void SetTransformation(ElementTransformation &_Ttr) { Ttr = &_Ttr; }
149 
150  /** @brief Evaluate the strain energy density function, W = W(Jpt).
151  @param[in] Jpt Represents the target->physical transformation
152  Jacobian matrix. */
153  virtual double EvalW(const DenseMatrix &Jpt) const = 0;
154 
155  /** @brief Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
156  @param[in] Jpt Represents the target->physical transformation
157  Jacobian matrix.
158  @param[out] P The evaluated 1st Piola-Kirchhoff stress tensor. */
159  virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const = 0;
160 
161  /** @brief Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor
162  and assemble its contribution to the local gradient matrix 'A'.
163  @param[in] Jpt Represents the target->physical transformation
164  Jacobian matrix.
165  @param[in] DS Gradient of the basis matrix (dof x dim).
166  @param[in] weight Quadrature weight coefficient for the point.
167  @param[in,out] A Local gradient matrix where the contribution from this
168  point will be added.
169 
170  Computes weight * d(dW_dxi)_d(xj) at the current point, for all i and j,
171  where x1 ... xn are the FE dofs. This function is usually defined using
172  the matrix invariants and their derivatives.
173  */
174  virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS,
175  const double weight, DenseMatrix &A) const = 0;
176 };
177 
178 
179 /** Inverse-harmonic hyperelastic model with a strain energy density function
180  given by the formula: W(J) = (1/2) det(J) Tr((J J^t)^{-1}) where J is the
181  deformation gradient. */
183 {
184 protected:
185  mutable DenseMatrix Z, S; // dim x dim
186  mutable DenseMatrix G, C; // dof x dim
187 
188 public:
189  virtual double EvalW(const DenseMatrix &J) const;
190 
191  virtual void EvalP(const DenseMatrix &J, DenseMatrix &P) const;
192 
193  virtual void AssembleH(const DenseMatrix &J, const DenseMatrix &DS,
194  const double weight, DenseMatrix &A) const;
195 };
196 
197 
198 /** Neo-Hookean hyperelastic model with a strain energy density function given
199  by the formula: \f$(\mu/2)(\bar{I}_1 - dim) + (K/2)(det(J)/g - 1)^2\f$ where
200  J is the deformation gradient and \f$\bar{I}_1 = (det(J))^{-2/dim} Tr(J
201  J^t)\f$. The parameters \f$\mu\f$ and K are the shear and bulk moduli,
202  respectively, and g is a reference volumetric scaling. */
204 {
205 protected:
206  mutable double mu, K, g;
209 
210  mutable DenseMatrix Z; // dim x dim
211  mutable DenseMatrix G, C; // dof x dim
212 
213  inline void EvalCoeffs() const;
214 
215 public:
216  NeoHookeanModel(double _mu, double _K, double _g = 1.0)
217  : mu(_mu), K(_K), g(_g), have_coeffs(false) { c_mu = c_K = c_g = NULL; }
218 
220  : mu(0.0), K(0.0), g(1.0), c_mu(&_mu), c_K(&_K), c_g(_g),
221  have_coeffs(true) { }
222 
223  virtual double EvalW(const DenseMatrix &J) const;
224 
225  virtual void EvalP(const DenseMatrix &J, DenseMatrix &P) const;
226 
227  virtual void AssembleH(const DenseMatrix &J, const DenseMatrix &DS,
228  const double weight, DenseMatrix &A) const;
229 };
230 
231 
232 /** Hyperelastic integrator for any given HyperelasticModel.
233 
234  Represents @f$ \int W(Jpt) dx @f$ over a target zone, where W is the
235  @a model's strain energy density function, and Jpt is the Jacobian of the
236  target->physical coordinates transformation. The target configuration is
237  given by the current mesh at the time of the evaluation of the integrator.
238 */
240 {
241 private:
242  HyperelasticModel *model;
243 
244  // Jrt: the Jacobian of the target-to-reference-element transformation.
245  // Jpr: the Jacobian of the reference-to-physical-element transformation.
246  // Jpt: the Jacobian of the target-to-physical-element transformation.
247  // P: represents dW_d(Jtp) (dim x dim).
248  // DSh: gradients of reference shape functions (dof x dim).
249  // DS: gradients of the shape functions in the target (stress-free)
250  // configuration (dof x dim).
251  // PMatI: coordinates of the deformed configuration (dof x dim).
252  // PMatO: reshaped view into the local element contribution to the operator
253  // output - the result of AssembleElementVector() (dof x dim).
254  DenseMatrix DSh, DS, Jrt, Jpr, Jpt, P, PMatI, PMatO;
255 
256 public:
257  /** @param[in] m HyperelasticModel that will be integrated. */
259 
260  /** @brief Computes the integral of W(Jacobian(Trt)) over a target zone
261  @param[in] el Type of FiniteElement.
262  @param[in] Ttr Represents ref->target coordinates transformation.
263  @param[in] elfun Physical coordinates of the zone. */
264  virtual double GetElementEnergy(const FiniteElement &el,
266  const Vector &elfun);
267 
268  virtual void AssembleElementVector(const FiniteElement &el,
270  const Vector &elfun, Vector &elvect);
271 
272  virtual void AssembleElementGrad(const FiniteElement &el,
274  const Vector &elfun, DenseMatrix &elmat);
275 };
276 
277 /** Hyperelastic incompressible Neo-Hookean integrator with the PK1 stress
278  \f$P = \mu F - p F^{-T}\f$ where \f$\mu\f$ is the shear modulus,
279  \f$p\f$ is the pressure, and \f$F\f$ is the deformation gradient */
281 {
282 private:
283  Coefficient *c_mu;
284  DenseMatrix DSh_u, DS_u, J0i, J, J1, Finv, P, F, FinvT;
285  DenseMatrix PMatI_u, PMatO_u, PMatI_p, PMatO_p, Z, G, C;
286  Vector Sh_p;
287 
288 public:
290 
291  virtual double GetElementEnergy(const Array<const FiniteElement *>&el,
293  const Array<const Vector *> &elfun);
294 
295  /// Perform the local action of the NonlinearFormIntegrator
298  const Array<const Vector *> &elfun,
299  const Array<Vector *> &elvec);
300 
301  /// Assemble the local gradient matrix
302  virtual void AssembleElementGrad(const Array<const FiniteElement*> &el,
304  const Array<const Vector *> &elfun,
305  const Array2D<DenseMatrix *> &elmats);
306 };
307 
309 {
310 private:
311  Coefficient *Q{};
312  DenseMatrix dshape, dshapex, EF, gradEF, ELV, elmat_comp;
313  Vector shape;
314  // PA extension
315  Vector pa_data;
316  const DofToQuad *maps; ///< Not owned
317  const GeometricFactors *geom; ///< Not owned
318  int dim, ne, nq;
319 public:
321 
322  VectorConvectionNLFIntegrator() = default;
323 
324  static const IntegrationRule &GetRule(const FiniteElement &fe,
326 
327  virtual void AssembleElementVector(const FiniteElement &el,
329  const Vector &elfun,
330  Vector &elvect);
331 
332  virtual void AssembleElementGrad(const FiniteElement &el,
334  const Vector &elfun,
335  DenseMatrix &elmat);
336 
338 
339  virtual void AssemblePA(const FiniteElementSpace &fes);
340 
341  virtual void AddMultPA(const Vector &x, Vector &y) const;
342 };
343 
344 }
345 
346 #endif
Abstract class for all finite elements.
Definition: fe.hpp:235
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
Definition: nonlininteg.cpp:31
virtual void EvalP(const DenseMatrix &J, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual double EvalW(const DenseMatrix &J) const
Evaluate the strain energy density function, W = W(Jpt).
void trans(const Vector &u, Vector &x)
Definition: ex27.cpp:421
virtual void AssembleElementGrad(const Array< const FiniteElement * > &el, ElementTransformation &Tr, const Array< const Vector * > &elfun, const Array2D< DenseMatrix * > &elmats)
Assemble the local gradient matrix.
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
virtual void AssembleFaceVector(const Array< const FiniteElement * > &el1, const Array< const FiniteElement * > &el2, FaceElementTransformations &Tr, const Array< const Vector * > &elfun, const Array< Vector * > &elvect)
Definition: nonlininteg.cpp:89
NeoHookeanModel(Coefficient &_mu, Coefficient &_K, Coefficient *_g=NULL)
virtual void AssembleH(const DenseMatrix &J, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:467
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
virtual double EvalW(const DenseMatrix &Jpt) const =0
Evaluate the strain energy density function, W = W(Jpt).
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition: mesh.hpp:1394
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
virtual void EvalP(const DenseMatrix &J, DenseMatrix &P) const
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual void AssembleFaceGrad(const Array< const FiniteElement * > &el1, const Array< const FiniteElement * > &el2, FaceElementTransformations &Tr, const Array< const Vector * > &elfun, const Array2D< DenseMatrix * > &elmats)
virtual void EvalP(const DenseMatrix &Jpt, DenseMatrix &P) const =0
Evaluate the 1st Piola-Kirchhoff stress tensor, P = P(Jpt).
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
virtual void AssembleElementGrad(const Array< const FiniteElement * > &el, ElementTransformation &Tr, const Array< const Vector * > &elfun, const Array2D< DenseMatrix * > &elmats)
Assemble the local gradient matrix.
void EvalCoeffs() const
virtual void AssemblePA(const FiniteElementSpace &fes)
Method defining partial assembly.
Definition: nonlininteg.cpp:18
virtual void AssembleFaceGrad(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat)
Assemble the local action of the gradient of the NonlinearFormIntegrator resulting from a face integr...
Definition: nonlininteg.cpp:61
virtual double GetElementEnergy(const Array< const FiniteElement * > &el, ElementTransformation &Tr, const Array< const Vector * > &elfun)
Compute the local energy.
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &Ttr, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
Definition: nonlininteg.cpp:37
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &Ttr, const Vector &elfun)
Computes the integral of W(Jacobian(Trt)) over a target zone.
Dynamic 2D array using row-major layout.
Definition: array.hpp:333
void SetTransformation(ElementTransformation &_Ttr)
virtual void AssembleElementGrad(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat)
Assemble the local gradient matrix.
Definition: nonlininteg.cpp:53
NeoHookeanModel(double _mu, double _K, double _g=1.0)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
virtual void AssembleFaceVector(const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator resulting from a face integral term...
Definition: nonlininteg.cpp:45
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
Definition: fe.hpp:131
ElementTransformation * Ttr
const IntegrationRule * IntRule
Definition: nonlininteg.hpp:29
VectorConvectionNLFIntegrator(Coefficient &q)
virtual void AssembleH(const DenseMatrix &J, const DenseMatrix &DS, const double weight, DenseMatrix &A) const
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
HyperelasticNLFIntegrator(HyperelasticModel *m)
virtual double EvalW(const DenseMatrix &J) const
Evaluate the strain energy density function, W = W(Jpt).
This class is used to express the local action of a general nonlinear finite element operator...
Definition: nonlininteg.hpp:26
void SetIntegrationRule(const IntegrationRule &irule)
Prescribe a fixed IntegrationRule to use.
Definition: nonlininteg.hpp:40
virtual void AssembleH(const DenseMatrix &Jpt, const DenseMatrix &DS, const double weight, DenseMatrix &A) const =0
Evaluate the derivative of the 1st Piola-Kirchhoff stress tensor and assemble its contribution to the...
virtual void AssembleElementVector(const Array< const FiniteElement * > &el, ElementTransformation &Tr, const Array< const Vector * > &elfun, const Array< Vector * > &elvec)
Perform the local action of the NonlinearFormIntegrator.
virtual void AddMultPA(const Vector &x, Vector &y) const
Method for partially assembled action.
void SetIntRule(const IntegrationRule *ir)
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == N...
Definition: nonlininteg.hpp:37
Abstract class for hyperelastic models.
IncompressibleNeoHookeanIntegrator(Coefficient &_mu)
virtual void AssembleElementVector(const FiniteElement &el, ElementTransformation &Ttr, const Vector &elfun, Vector &elvect)
Perform the local action of the NonlinearFormIntegrator.
virtual double GetElementEnergy(const Array< const FiniteElement * > &el, ElementTransformation &Tr, const Array< const Vector * > &elfun)
Compute the local energy.
Vector data type.
Definition: vector.hpp:51
NonlinearFormIntegrator(const IntegrationRule *ir=NULL)
Definition: nonlininteg.hpp:31
virtual double GetElementEnergy(const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
Compute the local energy.
Definition: nonlininteg.cpp:70
static const IntegrationRule & GetRule(const FiniteElement &fe, ElementTransformation &T)
virtual void AssembleElementVector(const Array< const FiniteElement * > &el, ElementTransformation &Tr, const Array< const Vector * > &elfun, const Array< Vector * > &elvec)
Perform the local action of the BlockNonlinearFormIntegrator.
Definition: nonlininteg.cpp:79