KASKADE 7 development version
Public Types | Public Member Functions | List of all members
Kaskade::DirichletNitscheBoundary< GridView, components, ScalarType > Class Template Reference

Dirichlet boundary conditions by Nitsche's method. More...

#include <boundaryConditions.hh>

Detailed Description

template<class GridView, int components, class ScalarType = typename GridView::ctype>
class Kaskade::DirichletNitscheBoundary< GridView, components, ScalarType >

Dirichlet boundary conditions by Nitsche's method.

For the Dirichlet boundary condition \( u|_{\partial \Omega} = u_D \), this defines the term

\[ \frac{\gamma}{2} h_F^{-1} |u-u_D|^2 - u^T u_x \kappa n\]

(and its derivatives), which can be added to the variational functional. The scaling by the face diameter \( h_F \) is supposed to provide nearly optimal convergence rates for linear finite elements.

This can lead to indefinite stiffness matrices if \( \gamma \) is not chosen sufficiently large, but does not cause ill-conditioning on mesh refinement, in contrast to DirichletPenaltyBoundary.

This follows R. Stenberg: On some techniques for approximating boundary conditions in the finite element method, J Comp. Appl. Math. 63:139-148, 1995.

Template Parameters
GridViewthe grid view on which the problem is defined (used for ctype and face type)
componentsthe number of vectorial components of the FE function
Scalarthe type of scalar entries of the FE function
Todo:
Consider extending to M. Juntunen, R. Stenberg: NITSCHE’S METHOD FOR GENERAL BOUNDARY CONDITIONS, Math. Comp. 78(267):1353-1374, 2009.

Definition at line 172 of file boundaryConditions.hh.

Public Types

using Vector = Dune::FieldVector< Scalar, components >
 

Public Member Functions

void moveTo (typename GridView::IntersectionIterator const &fi)
 Moves the boundary condition to a new face. More...
 
void setBoundaryData (Dune::FieldVector< typename GridView::ctype, dim-1 > const &xi, Scalar gamma_, Vector const &u, Vector const &ud, Dune::FieldMatrix< Scalar, components, dim > ux_, Dune::FieldMatrix< Scalar, dim, dim > const &kappa)
 Defines the data for the boundary condition. More...
 
Scalar d0 () const
 
Vector d1 (VariationalArg< Scalar, dim > const &arg) const
 
Dune::FieldMatrix< Scalar, components, components > d2 (VariationalArg< Scalar, dim > const &argTest, VariationalArg< Scalar, dim > const &argAnsatz) const
 

Member Typedef Documentation

◆ Vector

template<class GridView , int components, class ScalarType = typename GridView::ctype>
using Kaskade::DirichletNitscheBoundary< GridView, components, ScalarType >::Vector = Dune::FieldVector<Scalar,components>

Definition at line 180 of file boundaryConditions.hh.

Member Function Documentation

◆ d0()

template<class GridView , int components, class ScalarType = typename GridView::ctype>
Scalar Kaskade::DirichletNitscheBoundary< GridView, components, ScalarType >::d0 ( ) const
inline

Definition at line 220 of file boundaryConditions.hh.

◆ d1()

template<class GridView , int components, class ScalarType = typename GridView::ctype>
Vector Kaskade::DirichletNitscheBoundary< GridView, components, ScalarType >::d1 ( VariationalArg< Scalar, dim > const &  arg) const
inline

Definition at line 222 of file boundaryConditions.hh.

◆ d2()

template<class GridView , int components, class ScalarType = typename GridView::ctype>
Dune::FieldMatrix< Scalar, components, components > Kaskade::DirichletNitscheBoundary< GridView, components, ScalarType >::d2 ( VariationalArg< Scalar, dim > const &  argTest,
VariationalArg< Scalar, dim > const &  argAnsatz 
) const
inline

Definition at line 228 of file boundaryConditions.hh.

◆ moveTo()

template<class GridView , int components, class ScalarType = typename GridView::ctype>
void Kaskade::DirichletNitscheBoundary< GridView, components, ScalarType >::moveTo ( typename GridView::IntersectionIterator const &  fi)
inline

Moves the boundary condition to a new face.

This is called before calling d0,d1,d2 on a new face, and here the weight \( h_F^{-1} \) is computed once.

Definition at line 188 of file boundaryConditions.hh.

◆ setBoundaryData()

template<class GridView , int components, class ScalarType = typename GridView::ctype>
void Kaskade::DirichletNitscheBoundary< GridView, components, ScalarType >::setBoundaryData ( Dune::FieldVector< typename GridView::ctype, dim-1 > const &  xi,
Scalar  gamma_,
Vector const &  u,
Vector const &  ud,
Dune::FieldMatrix< Scalar, components, dim >  ux_,
Dune::FieldMatrix< Scalar, dim, dim > const &  kappa 
)
inline

Defines the data for the boundary condition.

Parameters
xiThe local position \( \xi \) within the face.
gammaThe penalty factor.
uThe current boundary value around which to linearize.
udThe Dirichlet value.
uxThe current solution's derivative at the boundary.
kappaThe diffusivity tensor.

Definition at line 209 of file boundaryConditions.hh.


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