Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Namespaces
LocalIntegrators Namespace Reference

Library of integrals over cells and faces. More...

Namespaces

namespace  Advection
 Local integrators related to advection along a vector field and its DG formulations.
 
namespace  Divergence
 Local integrators related to the divergence operator and its trace.
 
namespace  Elasticity
 Local integrators related to elasticity problems.
 
namespace  GradDiv
 
namespace  L2
 Local integrators related to L2-inner products.
 
namespace  Laplace
 Local integrators related to the Laplacian and its DG formulations.
 
namespace  Maxwell
 Local integrators related to curl operators and their traces.
 
namespace  Patches
 Integrators writing patches with values in quadrature points.
 

Detailed Description

Library of integrals over cells and faces.

This namespace contains application specific local integrals for bilinear forms, forms and error estimates. It is a collection of functions organized into namespaces devoted to certain applications. For instance, the namespace Laplace contains functions for computing cell matrices and cell residuals for the Laplacian operator, as well as functions for the weak boundary conditions by Nitsche or the interior penalty discontinuous Galerkin method. The namespace Maxwell does the same for curl-curl type problems.

The namespace L2 contains functions for mass matrices and L2-inner products.

Notational conventions

In most cases, the action of a function in this namespace can be described by a single integral. We distinguish between integrals over cells Z and over faces F. If an integral is denoted as

\[ \int_Z u \otimes v \,dx, \]

it will yield the following results, depending on the type of operation

We will use regular cursive symbols \(u\) for scalars and bold symbols \(\mathbf u\) for vectors. Test functions are always v and trial functions are always u. Parameters are Greek and the face normal vectors are \(\mathbf n = \mathbf n_1 = -\mathbf n_2\).

Signature of functions

Functions in this namespace follow a generic signature. In the simplest case, you have two related functions

template <int dim>
void
cell_matrix (
const FEValuesBase<dim>& fe,
const double factor = 1.);
template <int dim>
void
cell_residual (
const FEValuesBase<dim>& fe,
const std::vector<Tensor<1,dim> >& input,
const double factor = 1.);

There is typically a pair of functions for the same operator, the function cell_residual implementing the mapping of the operator from the finite element space into its dual, and the function cell_matrix generating the bilinear form corresponding to the Frechet derivative of cell_residual.

The first argument of these functions is the return type, which is

The next argument is the FEValuesBase object representing the finite element for integration. If the integrated operator maps from one finite element space into the dual of another (for instance an off-diagonal matrix in a block system), then first the FEValuesBase for the trial space and after this the one for the test space are specified.

This list is followed by the set of required data in the order

  1. Data vectors from finite element functions
  2. Data vectors from other objects
  3. Additional data
  4. A factor which is multiplied with the whole result

Usage

The local integrators can be used wherever a local integration loop would have been implemented instead. The following example is from the implementation of a Stokes solver, using MeshWorker::Assembler::LocalBlocksToGlobalBlocks. The matrices are

With these matrices, the function called by MeshWorker::loop() could be written like

using namespace dealii::LocalIntegrators;
template <int dim>
void MatrixIntegrator<dim>::cell(MeshWorker::DoFInfo<dim> &dinfo,
{
info.fe_values(0));
info.fe_values(0),
info.fe_values(1));
L2::cell_matrix (dinfo.matrix(2,false).matrix,
info.fe_values(1));
}
MatrixType matrix
const FEValuesBase< dim, spacedim > & fe_values() const
Access to finite element.
MatrixBlock< FullMatrix< number > > & matrix(const unsigned int i, const bool external=false)
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, double factor=1.)
Definition divergence.h:54
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition laplace.h:52

See step-39 for a worked out example of this code.