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
Functions
LocalIntegrators::Laplace Namespace Reference

Local integrators related to the Laplacian and its DG formulations. More...

Functions

template<int dim>
void cell_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
 
template<int dim>
void cell_residual (Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< Tensor< 1, dim > > &input, double factor=1.)
 
template<int dim>
void cell_residual (Vector< double > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &input, double factor=1.)
 
template<int dim>
void nitsche_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
 
template<int dim>
void nitsche_tangential_matrix (FullMatrix< double > &M, const FEValuesBase< dim > &fe, double penalty, double factor=1.)
 
template<int dim>
void nitsche_residual (Vector< double > &result, const FEValuesBase< dim > &fe, const std::vector< double > &input, const std::vector< Tensor< 1, dim > > &Dinput, const std::vector< double > &data, double penalty, double factor=1.)
 
template<int dim>
void nitsche_residual (Vector< double > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< double > > &input, const ArrayView< const std::vector< Tensor< 1, dim > > > &Dinput, const ArrayView< const std::vector< double > > &data, double penalty, double factor=1.)
 
template<int dim>
void ip_matrix (FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, double penalty, double factor1=1., double factor2=-1.)
 
template<int dim>
void ip_tangential_matrix (FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, double penalty, double factor1=1., double factor2=-1.)
 
template<int dim>
void ip_residual (Vector< double > &result1, Vector< double > &result2, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const std::vector< double > &input1, const std::vector< Tensor< 1, dim > > &Dinput1, const std::vector< double > &input2, const std::vector< Tensor< 1, dim > > &Dinput2, double pen, double int_factor=1., double ext_factor=-1.)
 
template<int dim>
void ip_residual (Vector< double > &result1, Vector< double > &result2, const FEValuesBase< dim > &fe1, const FEValuesBase< dim > &fe2, const ArrayView< const std::vector< double > > &input1, const ArrayView< const std::vector< Tensor< 1, dim > > > &Dinput1, const ArrayView< const std::vector< double > > &input2, const ArrayView< const std::vector< Tensor< 1, dim > > > &Dinput2, double pen, double int_factor=1., double ext_factor=-1.)
 
template<int dim, int spacedim, typename number >
double compute_penalty (const MeshWorker::DoFInfo< dim, spacedim, number > &dinfo1, const MeshWorker::DoFInfo< dim, spacedim, number > &dinfo2, unsigned int deg1, unsigned int deg2)
 

Detailed Description

Local integrators related to the Laplacian and its DG formulations.

Function Documentation

◆ cell_matrix()

template<int dim>
void LocalIntegrators::Laplace::cell_matrix ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
const double  factor = 1. 
)

Laplacian in weak form, namely on the cell Z the matrix

\[ \int_Z \nu \nabla u \cdot \nabla v \, dx. \]

The FiniteElement in fe may be scalar or vector valued. In the latter case, the Laplacian is applied to each component separately.

Definition at line 52 of file laplace.h.

◆ cell_residual() [1/2]

template<int dim>
void LocalIntegrators::Laplace::cell_residual ( Vector< double > &  result,
const FEValuesBase< dim > &  fe,
const std::vector< Tensor< 1, dim > > &  input,
double  factor = 1. 
)
inline

Laplacian residual operator in weak form

\[ \int_Z \nu \nabla u \cdot \nabla v \, dx. \]

Definition at line 92 of file laplace.h.

◆ cell_residual() [2/2]

template<int dim>
void LocalIntegrators::Laplace::cell_residual ( Vector< double > &  result,
const FEValuesBase< dim > &  fe,
const ArrayView< const std::vector< Tensor< 1, dim > > > &  input,
double  factor = 1. 
)
inline

Vector-valued Laplacian residual operator in weak form

\[ \int_Z \nu \nabla u : \nabla v \, dx. \]

Definition at line 119 of file laplace.h.

◆ nitsche_matrix()

template<int dim>
void LocalIntegrators::Laplace::nitsche_matrix ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
double  penalty,
double  factor = 1. 
)

Weak boundary condition of Nitsche type for the Laplacian, namely on the face F the matrix

\[ \int_F \Bigl(\gamma u v - \partial_n u v - u \partial_n v\Bigr)\;ds. \]

Here, \(\gamma\) is the penalty parameter suitably computed with compute_penalty().

Definition at line 157 of file laplace.h.

◆ nitsche_tangential_matrix()

template<int dim>
void LocalIntegrators::Laplace::nitsche_tangential_matrix ( FullMatrix< double > &  M,
const FEValuesBase< dim > &  fe,
double  penalty,
double  factor = 1. 
)

Weak boundary condition of Nitsche type for the Laplacian applied to the tangential component only, namely on the face F the matrix

\[ \int_F \Bigl(\gamma u_\tau v_\tau - \partial_n u_\tau v_\tau - u_\tau \partial_n v_\tau\Bigr)\;ds. \]

Here, \(\gamma\) is the penalty parameter suitably computed with compute_penalty().

Definition at line 198 of file laplace.h.

◆ nitsche_residual() [1/2]

template<int dim>
void LocalIntegrators::Laplace::nitsche_residual ( Vector< double > &  result,
const FEValuesBase< dim > &  fe,
const std::vector< double > &  input,
const std::vector< Tensor< 1, dim > > &  Dinput,
const std::vector< double > &  data,
double  penalty,
double  factor = 1. 
)

Weak boundary condition for the Laplace operator by Nitsche, scalar version, namely on the face F the vector

\[ \int_F \Bigl(\gamma (u-g) v - \partial_n u v - (u-g) \partial_n v\Bigr)\;ds. \]

Here, u is the finite element function whose values and gradient are given in the arguments input and Dinput, respectively. g is the inhomogeneous boundary value in the argument data. \(\gamma\) is the usual penalty parameter.

Definition at line 261 of file laplace.h.

◆ nitsche_residual() [2/2]

template<int dim>
void LocalIntegrators::Laplace::nitsche_residual ( Vector< double > &  result,
const FEValuesBase< dim > &  fe,
const ArrayView< const std::vector< double > > &  input,
const ArrayView< const std::vector< Tensor< 1, dim > > > &  Dinput,
const ArrayView< const std::vector< double > > &  data,
double  penalty,
double  factor = 1. 
)

Weak boundary condition for the Laplace operator by Nitsche, vector valued version, namely on the face F the vector

\[ \int_F \Bigl(\gamma (\mathbf u- \mathbf g) \cdot \mathbf v - \partial_n \mathbf u \cdot \mathbf v - (\mathbf u-\mathbf g) \cdot \partial_n \mathbf v\Bigr)\;ds. \]

Here, u is the finite element function whose values and gradient are given in the arguments input and Dinput, respectively. g is the inhomogeneous boundary value in the argument data. \(\gamma\) is the usual penalty parameter.

Definition at line 308 of file laplace.h.

◆ ip_matrix()

template<int dim>
void LocalIntegrators::Laplace::ip_matrix ( FullMatrix< double > &  M11,
FullMatrix< double > &  M12,
FullMatrix< double > &  M21,
FullMatrix< double > &  M22,
const FEValuesBase< dim > &  fe1,
const FEValuesBase< dim > &  fe2,
double  penalty,
double  factor1 = 1.,
double  factor2 = -1. 
)

Flux for the interior penalty method for the Laplacian, namely on the face F the matrices associated with the bilinear form

\[ \int_F \Bigl( \gamma [u][v] - \{\nabla u\}[v\mathbf n] - [u\mathbf n]\{\nabla v\} \Bigr) \; ds. \]

The penalty parameter should always be the mean value of the penalties needed for stability on each side. In the case of constant coefficients, it can be computed using compute_penalty().

If factor2 is missing or negative, the factor is assumed the same on both sides. If factors differ, note that the penalty parameter has to be computed accordingly.

Definition at line 359 of file laplace.h.

◆ ip_tangential_matrix()

template<int dim>
void LocalIntegrators::Laplace::ip_tangential_matrix ( FullMatrix< double > &  M11,
FullMatrix< double > &  M12,
FullMatrix< double > &  M21,
FullMatrix< double > &  M22,
const FEValuesBase< dim > &  fe1,
const FEValuesBase< dim > &  fe2,
double  penalty,
double  factor1 = 1.,
double  factor2 = -1. 
)

Flux for the interior penalty method for the Laplacian applied to the tangential components of a vector field, namely on the face F the matrices associated with the bilinear form

\[ \int_F \Bigl( \gamma [u_\tau][v_\tau] - \{\nabla u_\tau\}[v_\tau\mathbf n] - [u_\tau\mathbf n]\{\nabla v_\tau\} \Bigr) \; ds. \]

Warning
This function is still under development!

Definition at line 432 of file laplace.h.

◆ ip_residual() [1/2]

template<int dim>
void LocalIntegrators::Laplace::ip_residual ( Vector< double > &  result1,
Vector< double > &  result2,
const FEValuesBase< dim > &  fe1,
const FEValuesBase< dim > &  fe2,
const std::vector< double > &  input1,
const std::vector< Tensor< 1, dim > > &  Dinput1,
const std::vector< double > &  input2,
const std::vector< Tensor< 1, dim > > &  Dinput2,
double  pen,
double  int_factor = 1.,
double  ext_factor = -1. 
)

Residual term for the symmetric interior penalty method:

\[ \int_F \Bigl( \gamma [u][v] - \{\nabla u\}[v\mathbf n] - [u\mathbf n]\{\nabla v\} \Bigr) \; ds. \]

Definition at line 544 of file laplace.h.

◆ ip_residual() [2/2]

template<int dim>
void LocalIntegrators::Laplace::ip_residual ( Vector< double > &  result1,
Vector< double > &  result2,
const FEValuesBase< dim > &  fe1,
const FEValuesBase< dim > &  fe2,
const ArrayView< const std::vector< double > > &  input1,
const ArrayView< const std::vector< Tensor< 1, dim > > > &  Dinput1,
const ArrayView< const std::vector< double > > &  input2,
const ArrayView< const std::vector< Tensor< 1, dim > > > &  Dinput2,
double  pen,
double  int_factor = 1.,
double  ext_factor = -1. 
)

Vector-valued residual term for the symmetric interior penalty method:

\[ \int_F \Bigl( \gamma [\mathbf u]\cdot[\mathbf v] - \{\nabla \mathbf u\}[\mathbf v\otimes \mathbf n] - [\mathbf u\otimes \mathbf n]\{\nabla \mathbf v\} \Bigr) \; ds. \]

Definition at line 611 of file laplace.h.

◆ compute_penalty()

template<int dim, int spacedim, typename number >
double LocalIntegrators::Laplace::compute_penalty ( const MeshWorker::DoFInfo< dim, spacedim, number > &  dinfo1,
const MeshWorker::DoFInfo< dim, spacedim, number > &  dinfo2,
unsigned int  deg1,
unsigned int  deg2 
)

Auxiliary function computing the penalty parameter for interior penalty methods on rectangles.

Computation is done in two steps: first, we compute on each cell Zi the value Pi = pi(pi+1)/hi, where pi is the polynomial degree on cell Zi and hi is the length of Zi orthogonal to the current face.

Definition at line 685 of file laplace.h.