MFEM v4.7.0
Finite element discretization library
|
#include <lininteg.hpp>
Public Member Functions | |
DGElasticityDirichletLFIntegrator (VectorCoefficient &uD_, Coefficient &lambda_, Coefficient &mu_, real_t alpha_, real_t kappa_) | |
virtual void | AssembleRHSElementVect (const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) |
virtual void | AssembleRHSElementVect (const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect) |
virtual void | AssembleRHSElementVect (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, Vector &elvect) |
Public Member Functions inherited from mfem::LinearFormIntegrator | |
virtual bool | SupportsDevice () const |
Method probing for assembly on device. | |
virtual void | AssembleDevice (const FiniteElementSpace &fes, const Array< int > &markers, Vector &b) |
Method defining assembly on device. | |
virtual void | SetIntRule (const IntegrationRule *ir) |
const IntegrationRule * | GetIntRule () |
virtual | ~LinearFormIntegrator () |
Protected Attributes | |
VectorCoefficient & | uD |
Coefficient * | lambda |
Coefficient * | mu |
real_t | alpha |
real_t | kappa |
Vector | shape |
DenseMatrix | dshape |
DenseMatrix | adjJ |
DenseMatrix | dshape_ps |
Vector | nor |
Vector | dshape_dn |
Vector | dshape_du |
Vector | u_dir |
Protected Attributes inherited from mfem::LinearFormIntegrator | |
const IntegrationRule * | IntRule |
Additional Inherited Members | |
Protected Member Functions inherited from mfem::LinearFormIntegrator | |
LinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Boundary linear form integrator for imposing non-zero Dirichlet boundary conditions, in a DG elasticity formulation. Specifically, the linear form is given by
\[ \alpha \langle u_D, (\lambda \mathrm{div}(v) I + \mu (\nabla v + \nabla v^{\mathrm{T}})) \cdot n \rangle + + \kappa \langle h^{-1} (\lambda + 2 \mu) u_D, v \rangle, \]
where u_D is the given Dirichlet data. The parameters \(\alpha\), \(\kappa\), \(\lambda\) and \(\mu\), should match the parameters with the same names used in the bilinear form integrator, DGElasticityIntegrator.
Definition at line 578 of file lininteg.hpp.
|
inline |
Definition at line 597 of file lininteg.hpp.
|
virtual |
Given a particular Finite Element and a transformation (Tr) computes the element vector, elvect.
Implements mfem::LinearFormIntegrator.
Definition at line 867 of file lininteg.cpp.
|
virtual |
Reimplemented from mfem::LinearFormIntegrator.
Definition at line 873 of file lininteg.cpp.
|
virtual |
Reimplemented from mfem::LinearFormIntegrator.
Definition at line 49 of file lininteg.cpp.
|
protected |
Definition at line 588 of file lininteg.hpp.
|
protected |
Definition at line 583 of file lininteg.hpp.
|
protected |
Definition at line 587 of file lininteg.hpp.
|
protected |
Definition at line 591 of file lininteg.hpp.
|
protected |
Definition at line 592 of file lininteg.hpp.
|
protected |
Definition at line 589 of file lininteg.hpp.
|
protected |
Definition at line 583 of file lininteg.hpp.
|
protected |
Definition at line 582 of file lininteg.hpp.
|
protected |
Definition at line 582 of file lininteg.hpp.
|
protected |
Definition at line 590 of file lininteg.hpp.
|
protected |
Definition at line 586 of file lininteg.hpp.
|
protected |
Definition at line 593 of file lininteg.hpp.
|
protected |
Definition at line 581 of file lininteg.hpp.