MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::DGDirichletLFIntegrator Class Reference

#include <lininteg.hpp>

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

Public Member Functions

 DGDirichletLFIntegrator (Coefficient &u, const real_t s, const real_t k)
 
 DGDirichletLFIntegrator (Coefficient &u, Coefficient &q, const real_t s, const real_t k)
 
 DGDirichletLFIntegrator (Coefficient &u, MatrixCoefficient &q, const real_t s, const real_t k)
 
void AssembleRHSElementVect (const FiniteElement &el, ElementTransformation &Tr, Vector &elvect) override
 
void AssembleRHSElementVect (const FiniteElement &el, FaceElementTransformations &Tr, Vector &elvect) override
 
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 ~LinearFormIntegrator ()
 
- Public Member Functions inherited from mfem::Integrator
 Integrator (const IntegrationRule *ir=NULL)
 Create a new Integrator, optionally providing a prescribed quadrature rule to use in assembly.
 
virtual void SetIntRule (const IntegrationRule *ir)
 Prescribe a fixed IntegrationRule to use, or set to null to let the integrator choose an appropriate rule.
 
void SetIntegrationRule (const IntegrationRule &ir)
 Prescribe a fixed IntegrationRule to use. Sets the NURBS patch integration rule to null.
 
void SetNURBSPatchIntRule (NURBSMeshRules *pr)
 Sets an integration rule for use on NURBS patches.
 
bool HasNURBSPatchIntRule () const
 Check if a NURBS patch integration rule has been set.
 
const IntegrationRuleGetIntRule () const
 Directly return the IntRule pointer (possibly null) without checking for NURBS patch rules or falling back on a default.
 
const IntegrationRuleGetIntegrationRule () const
 Equivalent to GetIntRule, but retained for backward compatibility with applications.
 

Protected Attributes

CoefficientuD
 
CoefficientQ
 
MatrixCoefficientMQ
 
real_t sigma
 
real_t kappa
 
Vector shape
 
Vector dshape_dn
 
Vector nor
 
Vector nh
 
Vector ni
 
DenseMatrix dshape
 
DenseMatrix mq
 
DenseMatrix adjJ
 
- Protected Attributes inherited from mfem::Integrator
const IntegrationRuleIntRule
 
NURBSMeshRulespatchRules = nullptr
 

Additional Inherited Members

- Protected Member Functions inherited from mfem::LinearFormIntegrator
 LinearFormIntegrator (const IntegrationRule *ir=NULL)
 
- Protected Member Functions inherited from mfem::Integrator
const IntegrationRuleGetIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const
 Returns an integration rule based on the arguments and internal state of the Integrator object.
 
const IntegrationRuleGetIntegrationRule (const FiniteElement &el, const ElementTransformation &trans) const
 Returns an integration rule based on the arguments and internal state. (Version for identical trial_fe and test_fe)
 
virtual const IntegrationRuleGetDefaultIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const
 Subclasses should override to choose a default integration rule.
 

Detailed Description

Boundary linear integrator for imposing non-zero Dirichlet boundary conditions, to be used in conjunction with DGDiffusionIntegrator. Specifically, given the Dirichlet data \(u_D\), the linear form assembles the following integrals on the boundary:

\[ \sigma \langle u_D, (Q \nabla v)) \cdot n \rangle + \kappa \langle {h^{-1} Q} u_D, v \rangle, \]

where Q is a scalar or matrix diffusion coefficient and v is the test function. The parameters \(\sigma\) and \(\kappa\) should be the same as the ones used in the DGDiffusionIntegrator.

Definition at line 552 of file lininteg.hpp.

Constructor & Destructor Documentation

◆ DGDirichletLFIntegrator() [1/3]

mfem::DGDirichletLFIntegrator::DGDirichletLFIntegrator ( Coefficient & u,
const real_t s,
const real_t k )
inline

Definition at line 564 of file lininteg.hpp.

◆ DGDirichletLFIntegrator() [2/3]

mfem::DGDirichletLFIntegrator::DGDirichletLFIntegrator ( Coefficient & u,
Coefficient & q,
const real_t s,
const real_t k )
inline

Definition at line 566 of file lininteg.hpp.

◆ DGDirichletLFIntegrator() [3/3]

mfem::DGDirichletLFIntegrator::DGDirichletLFIntegrator ( Coefficient & u,
MatrixCoefficient & q,
const real_t s,
const real_t k )
inline

Definition at line 569 of file lininteg.hpp.

Member Function Documentation

◆ AssembleRHSElementVect() [1/3]

void mfem::DGDirichletLFIntegrator::AssembleRHSElementVect ( const FiniteElement & el,
ElementTransformation & Tr,
Vector & elvect )
overridevirtual

Given a particular Finite Element and a transformation (Tr) computes the element vector, elvect.

Implements mfem::LinearFormIntegrator.

Definition at line 810 of file lininteg.cpp.

◆ AssembleRHSElementVect() [2/3]

void mfem::DGDirichletLFIntegrator::AssembleRHSElementVect ( const FiniteElement & el,
FaceElementTransformations & Tr,
Vector & elvect )
overridevirtual

Reimplemented from mfem::LinearFormIntegrator.

Definition at line 816 of file lininteg.cpp.

◆ AssembleRHSElementVect() [3/3]

void mfem::LinearFormIntegrator::AssembleRHSElementVect ( const FiniteElement & el1,
const FiniteElement & el2,
FaceElementTransformations & Tr,
Vector & elvect )
virtual

Reimplemented from mfem::LinearFormIntegrator.

Definition at line 49 of file lininteg.cpp.

Member Data Documentation

◆ adjJ

DenseMatrix mfem::DGDirichletLFIntegrator::adjJ
protected

Definition at line 561 of file lininteg.hpp.

◆ dshape

DenseMatrix mfem::DGDirichletLFIntegrator::dshape
protected

Definition at line 561 of file lininteg.hpp.

◆ dshape_dn

Vector mfem::DGDirichletLFIntegrator::dshape_dn
protected

Definition at line 560 of file lininteg.hpp.

◆ kappa

real_t mfem::DGDirichletLFIntegrator::kappa
protected

Definition at line 557 of file lininteg.hpp.

◆ MQ

MatrixCoefficient* mfem::DGDirichletLFIntegrator::MQ
protected

Definition at line 556 of file lininteg.hpp.

◆ mq

DenseMatrix mfem::DGDirichletLFIntegrator::mq
protected

Definition at line 561 of file lininteg.hpp.

◆ nh

Vector mfem::DGDirichletLFIntegrator::nh
protected

Definition at line 560 of file lininteg.hpp.

◆ ni

Vector mfem::DGDirichletLFIntegrator::ni
protected

Definition at line 560 of file lininteg.hpp.

◆ nor

Vector mfem::DGDirichletLFIntegrator::nor
protected

Definition at line 560 of file lininteg.hpp.

◆ Q

Coefficient * mfem::DGDirichletLFIntegrator::Q
protected

Definition at line 555 of file lininteg.hpp.

◆ shape

Vector mfem::DGDirichletLFIntegrator::shape
protected

Definition at line 560 of file lininteg.hpp.

◆ sigma

real_t mfem::DGDirichletLFIntegrator::sigma
protected

Definition at line 557 of file lininteg.hpp.

◆ uD

Coefficient* mfem::DGDirichletLFIntegrator::uD
protected

Definition at line 555 of file lininteg.hpp.


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