MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::DGElasticityIntegrator Class Reference

#include <bilininteg.hpp>

Inheritance diagram for mfem::DGElasticityIntegrator:
[legend]
Collaboration diagram for mfem::DGElasticityIntegrator:
[legend]

Public Member Functions

 DGElasticityIntegrator (real_t alpha_, real_t kappa_)
 
 DGElasticityIntegrator (Coefficient &lambda_, Coefficient &mu_, real_t alpha_, real_t kappa_)
 
virtual void AssembleFaceMatrix (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleFaceMatrix (const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly.
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
virtual void AssembleNURBSPA (const FiniteElementSpace &fes)
 Method defining partial assembly on NURBS patches.
 
virtual void AssemblePABoundary (const FiniteElementSpace &fes)
 
virtual void AssemblePAInteriorFaces (const FiniteElementSpace &fes)
 
virtual void AssemblePABoundaryFaces (const FiniteElementSpace &fes)
 
virtual void AssembleDiagonalPA (Vector &diag)
 Assemble diagonal and add it to Vector diag.
 
virtual void AssembleDiagonalPA_ADAt (const Vector &D, Vector &diag)
 Assemble diagonal of \(A D A^T\) ( \(A\) is this integrator) and add it to diag.
 
virtual void AddMultPA (const Vector &x, Vector &y) const
 Method for partially assembled action.
 
virtual void AddMultNURBSPA (const Vector &x, Vector &y) const
 Method for partially assembled action on NURBS patches.
 
virtual void AddMultTransposePA (const Vector &x, Vector &y) const
 Method for partially assembled transposed action.
 
virtual void AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add=true)
 Method defining element assembly.
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining matrix-free assembly.
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
virtual void AddMultTransposeMF (const Vector &x, Vector &y) const
 
virtual void AssembleDiagonalMF (Vector &diag)
 Assemble diagonal and add it to Vector diag.
 
virtual void AssembleEAInteriorFaces (const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add=true)
 
virtual void AssembleEABoundaryFaces (const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add=true)
 
virtual void AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
 Given a particular Finite Element computes the element matrix elmat.
 
virtual void AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
 
virtual void AssemblePatchMatrix (const int patch, const FiniteElementSpace &fes, SparseMatrix *&smat)
 
virtual void AssembleTraceFaceMatrix (int elem, const FiniteElement &trial_face_fe, const FiniteElement &test_fe, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the BilinearFormIntegrator. Note that the default implementation in the base class is general but not efficient.
 
virtual void AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the BilinearFormIntegrator resulting from a face integral term. Note that the default implementation in the base class is general but not efficient.
 
virtual void AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat)
 Assemble the local gradient matrix.
 
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 integral term.
 
virtual void ComputeElementFlux (const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true, const IntegrationRule *ir=NULL)
 Virtual method required for Zienkiewicz-Zhu type error estimators.
 
virtual real_t ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
 Virtual method required for Zienkiewicz-Zhu type error estimators.
 
virtual bool RequiresFaceNormalDerivatives () const
 For bilinear forms on element faces, specifies if the normal derivatives are needed on the faces or just the face restriction.
 
virtual void AddMultPAFaceNormalDerivatives (const Vector &x, const Vector &dxdn, Vector &y, Vector &dydn) const
 Method for partially assembled action.
 
virtual ~BilinearFormIntegrator ()
 
- Public Member Functions inherited from mfem::NonlinearFormIntegrator
virtual void SetIntRule (const IntegrationRule *ir)
 Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL).
 
void SetIntegrationMode (Mode m)
 
void SetNURBSPatchIntRule (NURBSMeshRules *pr)
 For patchwise integration, SetNURBSPatchIntRule must be called.
 
bool HasNURBSPatchIntRule () const
 
bool Patchwise () const
 
void SetIntegrationRule (const IntegrationRule &ir)
 Prescribe a fixed IntegrationRule to use.
 
void SetPAMemoryType (MemoryType mt)
 
const IntegrationRuleGetIntegrationRule () const
 Get the integration rule of the integrator (possibly NULL).
 
virtual real_t GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
 Compute the local energy.
 
virtual void AssembleGradPA (const Vector &x, const FiniteElementSpace &fes)
 Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at the state x.
 
virtual real_t GetLocalStateEnergyPA (const Vector &x) const
 Compute the local (to the MPI rank) energy with partial assembly.
 
virtual void AddMultGradPA (const Vector &x, Vector &y) const
 Method for partially assembled gradient action.
 
virtual void AssembleGradDiagonalPA (Vector &diag) const
 Method for computing the diagonal of the gradient with partial assembly.
 
virtual bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend.
 
ceed::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 

Static Protected Member Functions

static void AssembleBlock (const int dim, const int row_ndofs, const int col_ndofs, const int row_offset, const int col_offset, const real_t jmatcoef, const Vector &col_nL, const Vector &col_nM, const Vector &row_shape, const Vector &col_shape, const Vector &col_dshape_dnM, const DenseMatrix &col_dshape, DenseMatrix &elmat, DenseMatrix &jmat)
 

Protected Attributes

Coefficientlambda
 
Coefficientmu
 
real_t alpha
 
real_t kappa
 
Vector shape1
 
Vector shape2
 
DenseMatrix dshape1
 
DenseMatrix dshape2
 
DenseMatrix adjJ
 
DenseMatrix dshape1_ps
 
DenseMatrix dshape2_ps
 
Vector nor
 
Vector nL1
 
Vector nL2
 
Vector nM1
 
Vector nM2
 
Vector dshape1_dnM
 
Vector dshape2_dnM
 
DenseMatrix jmat
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
const IntegrationRuleIntRule
 
Mode integrationMode = Mode::ELEMENTWISE
 
NURBSMeshRulespatchRules = nullptr
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

Additional Inherited Members

- Public Types inherited from mfem::NonlinearFormIntegrator
enum  Mode { ELEMENTWISE = 0 , PATCHWISE = 1 , PATCHWISE_REDUCED = 2 }
 
- Protected Member Functions inherited from mfem::BilinearFormIntegrator
 BilinearFormIntegrator (const IntegrationRule *ir=NULL)
 
- Protected Member Functions inherited from mfem::NonlinearFormIntegrator
 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 

Detailed Description

Integrator for the DG elasticity form, for the formulations see:

  • PhD Thesis of Jonas De Basabe, High-Order Finite Element Methods for Seismic Wave Propagation, UT Austin, 2009, p. 23, and references therein
  • Peter Hansbo and Mats G. Larson, Discontinuous Galerkin and the Crouzeix-Raviart Element: Application to Elasticity, PREPRINT 2000-09, p.3

    \[ - \left< \{ \tau(u) \}, [v] \right> + \alpha \left< \{ \tau(v) \}, [u] \right> + \kappa \left< h^{-1} \{ \lambda + 2 \mu \} [u], [v] \right> \]

where \( \left<u, v\right> = \int_{F} u \cdot v \), and \( F \) is a face which is either a boundary face \( F_b \) of an element \( K \) or an interior face \( F_i \) separating elements \( K_1 \) and \( K_2 \).

In the bilinear form above \( \tau(u) \) is traction, and it's also \( \tau(u) = \sigma(u) \cdot \vec{n} \), where \( \sigma(u) \) is stress, and \( \vec{n} \) is the unit normal vector w.r.t. to \( F \).

In other words, we have

\[ - \left< \{ \sigma(u) \cdot \vec{n} \}, [v] \right> + \alpha \left< \{ \sigma(v) \cdot \vec{n} \}, [u] \right> + \kappa \left< h^{-1} \{ \lambda + 2 \mu \} [u], [v] \right> \]

For isotropic media

\[ \begin{split} \sigma(u) &= \lambda \nabla \cdot u I + 2 \mu \varepsilon(u) \\ &= \lambda \nabla \cdot u I + 2 \mu \frac{1}{2} (\nabla u + \nabla u^{\mathrm{T}}) \\ &= \lambda \nabla \cdot u I + \mu (\nabla u + \nabla u^{\mathrm{T}}) \end{split} \]

where \( I \) is identity matrix, \( \lambda \) and \( \mu \) are Lame coefficients (see ElasticityIntegrator), \( u, v \) are the trial and test functions, respectively.

The parameters \( \alpha \) and \( \kappa \) determine the DG method to use (when this integrator is added to the "broken" ElasticityIntegrator):

  • IIPG, \(\alpha = 0\), C. Dawson, S. Sun, M. Wheeler, Compatible algorithms for coupled flow and transport, Comp. Meth. Appl. Mech. Eng., 193(23-26), 2565-2580, 2004.
  • SIPG, \(\alpha = -1\), M. Grote, A. Schneebeli, D. Schotzau, Discontinuous Galerkin Finite Element Method for the Wave Equation, SINUM, 44(6), 2408-2431, 2006.
  • NIPG, \(\alpha = 1\), B. Riviere, M. Wheeler, V. Girault, A Priori Error Estimates for Finite Element Methods Based on Discontinuous Approximation Spaces for Elliptic Problems, SINUM, 39(3), 902-931, 2001.

This is a 'Vector' integrator, i.e. defined for FE spaces using multiple copies of a scalar FE space.

Definition at line 3422 of file bilininteg.hpp.

Constructor & Destructor Documentation

◆ DGElasticityIntegrator() [1/2]

mfem::DGElasticityIntegrator::DGElasticityIntegrator ( real_t alpha_,
real_t kappa_ )
inline

Definition at line 3425 of file bilininteg.hpp.

◆ DGElasticityIntegrator() [2/2]

mfem::DGElasticityIntegrator::DGElasticityIntegrator ( Coefficient & lambda_,
Coefficient & mu_,
real_t alpha_,
real_t kappa_ )
inline

Definition at line 3428 of file bilininteg.hpp.

Member Function Documentation

◆ AssembleBlock()

void mfem::DGElasticityIntegrator::AssembleBlock ( const int dim,
const int row_ndofs,
const int col_ndofs,
const int row_offset,
const int col_offset,
const real_t jmatcoef,
const Vector & col_nL,
const Vector & col_nM,
const Vector & row_shape,
const Vector & col_shape,
const Vector & col_dshape_dnM,
const DenseMatrix & col_dshape,
DenseMatrix & elmat,
DenseMatrix & jmat )
staticprotected

Definition at line 3664 of file bilininteg.cpp.

◆ AssembleFaceMatrix() [1/2]

void mfem::DGElasticityIntegrator::AssembleFaceMatrix ( const FiniteElement & el1,
const FiniteElement & el2,
FaceElementTransformations & Trans,
DenseMatrix & elmat )
virtual

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 3707 of file bilininteg.cpp.

◆ AssembleFaceMatrix() [2/2]

void mfem::BilinearFormIntegrator::AssembleFaceMatrix ( const FiniteElement & trial_face_fe,
const FiniteElement & test_fe1,
const FiniteElement & test_fe2,
FaceElementTransformations & Trans,
DenseMatrix & elmat )
virtual

Abstract method used for assembling TraceFaceIntegrators in a MixedBilinearForm.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 162 of file bilininteg.cpp.

Member Data Documentation

◆ adjJ

DenseMatrix mfem::DGElasticityIntegrator::adjJ
protected

Definition at line 3450 of file bilininteg.hpp.

◆ alpha

real_t mfem::DGElasticityIntegrator::alpha
protected

Definition at line 3440 of file bilininteg.hpp.

◆ dshape1

DenseMatrix mfem::DGElasticityIntegrator::dshape1
protected

Definition at line 3448 of file bilininteg.hpp.

◆ dshape1_dnM

Vector mfem::DGElasticityIntegrator::dshape1_dnM
protected

Definition at line 3458 of file bilininteg.hpp.

◆ dshape1_ps

DenseMatrix mfem::DGElasticityIntegrator::dshape1_ps
protected

Definition at line 3454 of file bilininteg.hpp.

◆ dshape2

DenseMatrix mfem::DGElasticityIntegrator::dshape2
protected

Definition at line 3448 of file bilininteg.hpp.

◆ dshape2_dnM

Vector mfem::DGElasticityIntegrator::dshape2_dnM
protected

Definition at line 3458 of file bilininteg.hpp.

◆ dshape2_ps

DenseMatrix mfem::DGElasticityIntegrator::dshape2_ps
protected

Definition at line 3454 of file bilininteg.hpp.

◆ jmat

DenseMatrix mfem::DGElasticityIntegrator::jmat
protected

Definition at line 3460 of file bilininteg.hpp.

◆ kappa

real_t mfem::DGElasticityIntegrator::kappa
protected

Definition at line 3440 of file bilininteg.hpp.

◆ lambda

Coefficient* mfem::DGElasticityIntegrator::lambda
protected

Definition at line 3439 of file bilininteg.hpp.

◆ mu

Coefficient * mfem::DGElasticityIntegrator::mu
protected

Definition at line 3439 of file bilininteg.hpp.

◆ nL1

Vector mfem::DGElasticityIntegrator::nL1
protected

Definition at line 3456 of file bilininteg.hpp.

◆ nL2

Vector mfem::DGElasticityIntegrator::nL2
protected

Definition at line 3456 of file bilininteg.hpp.

◆ nM1

Vector mfem::DGElasticityIntegrator::nM1
protected

Definition at line 3457 of file bilininteg.hpp.

◆ nM2

Vector mfem::DGElasticityIntegrator::nM2
protected

Definition at line 3457 of file bilininteg.hpp.

◆ nor

Vector mfem::DGElasticityIntegrator::nor
protected

Definition at line 3455 of file bilininteg.hpp.

◆ shape1

Vector mfem::DGElasticityIntegrator::shape1
protected

Definition at line 3445 of file bilininteg.hpp.

◆ shape2

Vector mfem::DGElasticityIntegrator::shape2
protected

Definition at line 3445 of file bilininteg.hpp.


The documentation for this class was generated from the following files: