![]() |
MFEM
v4.2.0
Finite element discretization library
|
#include <bilininteg.hpp>
Public Member Functions | |
DGElasticityIntegrator (double alpha_, double kappa_) | |
DGElasticityIntegrator (Coefficient &lambda_, Coefficient &mu_, double alpha_, double kappa_) | |
virtual void | AssembleFaceMatrix (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat) |
![]() | |
virtual void | AssemblePA (const FiniteElementSpace &fes) |
Method defining partial assembly. More... | |
virtual void | AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_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. More... | |
virtual void | AssembleDiagonalPA_ADAt (const Vector &D, Vector &diag) |
Assemble diagonal of ADA^T (A is this integrator) and add it to diag. More... | |
virtual void | AddMultPA (const Vector &x, Vector &y) const |
Method for partially assembled action. More... | |
virtual void | AddMultTransposePA (const Vector &x, Vector &y) const |
Method for partially assembled transposed action. More... | |
virtual void | AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add=true) |
Method defining element assembly. More... | |
virtual void | AssembleMF (const FiniteElementSpace &fes) |
Method defining matrix-free assembly. More... | |
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. More... | |
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. More... | |
virtual void | AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat) |
virtual void | AssembleFaceMatrix (const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, 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. More... | |
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. More... | |
virtual void | AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat) |
Assemble the local gradient matrix. More... | |
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. More... | |
virtual void | ComputeElementFlux (const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true) |
Virtual method required for Zienkiewicz-Zhu type error estimators. More... | |
virtual double | ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL) |
Virtual method required for Zienkiewicz-Zhu type error estimators. More... | |
virtual | ~BilinearFormIntegrator () |
![]() | |
void | SetIntRule (const IntegrationRule *ir) |
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL). More... | |
void | SetIntegrationRule (const IntegrationRule &irule) |
Prescribe a fixed IntegrationRule to use. More... | |
virtual double | GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun) |
Compute the local energy. More... | |
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 double 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 | |
Coefficient * | lambda |
Coefficient * | mu |
double | alpha |
double | 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 |
![]() | |
const IntegrationRule * | IntRule |
Additional Inherited Members | |
![]() | |
BilinearFormIntegrator (const IntegrationRule *ir=NULL) | |
![]() | |
NonlinearFormIntegrator (const IntegrationRule *ir=NULL) | |
Integrator for the DG elasticity form, for the formulations see:
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^T) \\ &= \lambda \nabla \cdot u I + \mu (\nabla u + \nabla u^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):
This is a 'Vector' integrator, i.e. defined for FE spaces using multiple copies of a scalar FE space.
Definition at line 2815 of file bilininteg.hpp.
|
inline |
Definition at line 2818 of file bilininteg.hpp.
|
inline |
Definition at line 2821 of file bilininteg.hpp.
|
staticprotected |
Definition at line 2975 of file bilininteg.cpp.
|
virtual |
Reimplemented from mfem::BilinearFormIntegrator.
Definition at line 3018 of file bilininteg.cpp.
|
protected |
Definition at line 2843 of file bilininteg.hpp.
|
protected |
Definition at line 2833 of file bilininteg.hpp.
|
protected |
Definition at line 2841 of file bilininteg.hpp.
|
protected |
Definition at line 2851 of file bilininteg.hpp.
|
protected |
Definition at line 2847 of file bilininteg.hpp.
|
protected |
Definition at line 2841 of file bilininteg.hpp.
|
protected |
Definition at line 2851 of file bilininteg.hpp.
|
protected |
Definition at line 2847 of file bilininteg.hpp.
|
protected |
Definition at line 2853 of file bilininteg.hpp.
|
protected |
Definition at line 2833 of file bilininteg.hpp.
|
protected |
Definition at line 2832 of file bilininteg.hpp.
|
protected |
Definition at line 2832 of file bilininteg.hpp.
|
protected |
Definition at line 2849 of file bilininteg.hpp.
|
protected |
Definition at line 2849 of file bilininteg.hpp.
|
protected |
Definition at line 2850 of file bilininteg.hpp.
|
protected |
Definition at line 2850 of file bilininteg.hpp.
|
protected |
Definition at line 2848 of file bilininteg.hpp.
|
protected |
Definition at line 2838 of file bilininteg.hpp.
|
protected |
Definition at line 2838 of file bilininteg.hpp.