DUNE PDELab (2.8)

Dune::PDELab::LocalOperatorInterface Class Reference

Class to document the stationary local operator interface. More...

#include <dune/pdelab/localoperator/interface.hh>

Methods for selective assembly

template<typename EG >
bool skip_entity (const EG &eg) const
 whether to assembly methods associated with a given entity More...
 
template<typename IG >
bool skip_intersection (const IG &ig) const
 whether to assembly methods associated with a given intersection More...
 

Methods for the sparsity pattern

template<typename LFSU , typename LFSV , typename LocalPattern >
void pattern_volume (const LFSU &lfsu, const LFSV &lfsv, LocalPattern &pattern) const
 get an element's contribution to the sparsity pattern More...
 
template<typename LFSU , typename LFSV , typename LocalPattern >
void pattern_volume_post_skeleton (const LFSU &lfsu, const LFSV &lfsv, LocalPattern &pattern) const
 get an element's contribution to the sparsity pattern after the intersections have been handled More...
 
template<typename LFSU , typename LFSV , typename LocalPattern >
void pattern_skeleton (const LFSU &lfsu_s, const LFSV &lfsv_s, const LFSU &lfsu_n, const LFSV &lfsv_n, LocalPattern &pattern_sn, LocalPattern &pattern_ns) const
 get an internal intersection's contribution to the sparsity pattern More...
 
template<typename LFSU , typename LFSV , typename LocalPattern >
void pattern_boundary (const LFSU &lfsu_s, const LFSV &lfsv_s, LocalPattern &pattern_ss) const
 get a boundary intersection's contribution to the sparsity pattern More...
 

Methods for the residual – non-constant parts

template<typename EG , typename LFSU , typename X , typename LFSV , typename R >
void alpha_volume (const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, R &r) const
 get an element's contribution to alpha More...
 
template<typename EG , typename LFSU , typename X , typename LFSV , typename R >
void alpha_volume_post_skeleton (const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, R &r) const
 get an element's contribution to alpha after the intersections have been handled More...
 
template<typename IG , typename LFSU , typename X , typename LFSV , typename R >
void alpha_skeleton (const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, const LFSU &lfsu_n, const X &x_n, const LFSV &lfsv_n, R &r_s, R &r_n) const
 get an internal intersections's contribution to alpha More...
 
template<typename IG , typename LFSU , typename X , typename LFSV , typename R >
void alpha_boundary (const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, R &r_s) const
 get a boundary intersections's contribution to alpha More...
 

Methods for the residual – constant parts

template<typename EG , typename LFSV , typename R >
void lambda_volume (const EG &eg, const LFSV &lfsv, R &r) const
 get an element's contribution to lambda More...
 
template<typename EG , typename LFSV , typename R >
void lambda_volume_post_skeleton (const EG &eg, const LFSV &lfsv, R &r) const
 get an element's contribution to lambda after the intersections have been handled More...
 
template<typename IG , typename LFSV , typename R >
void lambda_skeleton (const IG &ig, const LFSV &lfsv_s, const LFSV &lfsv_n, R &r_s, R &r_n) const
 get an internal intersections's contribution to lambda More...
 
template<typename IG , typename LFSV , typename R >
void lambda_boundary (const IG &ig, const LFSV &lfsv_s, R &r_s) const
 get a boundary intersections's contribution to lambda More...
 

Methods for the application of the jacobian for linear problems

template<typename EG , typename LFSU , typename X , typename LFSV , typename Y >
void jacobian_apply_volume (const EG &eg, const LFSU &lfsu, const X &z, const LFSV &lfsv, Y &y) const
 Applies an element's jacobian to a vector for a linear problem. More...
 
template<typename EG , typename LFSU , typename X , typename LFSV , typename Y >
void jacobian_apply_volume_post_skeleton (const EG &eg, const LFSU &lfsu, const X &z, const LFSV &lfsv, Y &y) const
 
template<typename IG , typename LFSU , typename X , typename LFSV , typename Y >
void jacobian_apply_skeleton (const IG &ig, const LFSU &lfsu_s, const X &z_s, const LFSV &lfsv_s, const LFSU &lfsu_n, const X &z_n, const LFSV &lfsv_n, Y &y_s, Y &y_n) const
 
template<typename IG , typename LFSU , typename X , typename LFSV , typename Y >
void jacobian_apply_boundary (const IG &ig, const LFSU &lfsu_s, const X &z_s, const LFSV &lfsv_s, Y &y_s) const
 apply a boundary intersections's jacobian for a linear problem. More...
 

Methods for the application of the jacobian for nonlinear problems

template<typename EG , typename LFSU , typename X , typename Z , typename LFSV , typename Y >
void jacobian_apply_volume (const EG &eg, const LFSU &lfsu, const X &x, const Z &z, const LFSV &lfsv, Y &y) const
 Applies an element's jacobian to a vector for a nonlinear problem. More...
 
template<typename EG , typename LFSU , typename X , typename Z , typename LFSV , typename Y >
void jacobian_apply_volume_post_skeleton (const EG &eg, const LFSU &lfsu, const X &x, const Z &z, const LFSV &lfsv, Y &y) const
 
template<typename IG , typename LFSU , typename X , typename Z , typename LFSV , typename Y >
void jacobian_apply_skeleton (const IG &ig, const LFSU &lfsu_s, const X &x_s, const Z &z_s, const LFSV &lfsv_s, const LFSU &lfsu_n, const X &x_n, const Z &z_n, const LFSV &lfsv_n, Y &y_s, Y &y_n) const
 
template<typename IG , typename LFSU , typename X , typename Z , typename LFSV , typename Y >
void jacobian_apply_boundary (const IG &ig, const LFSU &lfsu_s, const X &x_s, const Z &z_s, const LFSV &lfsv_s, Y &y_s) const
 apply a boundary intersections's jacobian for a nonlinear problem. More...
 

Methods to extract the jacobian

template<typename EG , typename LFSU , typename X , typename LFSV , typename LocalMatrix >
void jacobian_volume (const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, LocalMatrix &mat) const
 get an element's jacobian More...
 
template<typename EG , typename LFSU , typename X , typename LFSV , typename LocalMatrix >
void jacobian_volume_post_skeleton (const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, LocalMatrix &mat) const
 get an element's jacobian after the intersections have been handled More...
 
template<typename IG , typename LFSU , typename X , typename LFSV , typename LocalMatrix >
void jacobian_skeleton (const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, const LFSU &lfsu_n, const X &x_n, const LFSV &lfsv_n, LocalMatrix &mat_ss, LocalMatrix &mat_sn, LocalMatrix &mat_ns, LocalMatrix &mat_nn) const
 apply an internal intersections's jacobians More...
 
template<typename IG , typename LFSU , typename X , typename LFSV , typename LocalMatrix >
void jacobian_boundary (const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, LocalMatrix &mat_ss) const
 get a boundary intersections's jacobian More...
 

Flags selective assembly

enum  
 Whether to do selective assembly on the elements, i.e. whether or not skip_entity() should be called.
 
enum  
 Whether to do selective assembly on the intersections, i.e. whether or not skip_intersection() should be called.
 

Flags for the sparsity pattern

enum  
 Whether to assemble the pattern on the elements, i.e. whether or not pattern_volume() should be called.
 
enum  
 Whether to assemble the pattern on the elements after the skeleton has been handled, i.e. whether or not pattern_volume_post_skeleton() should be called.
 
enum  
 Whether to assemble the pattern on the interior intersections, i.e. whether or not pattern_skeleton() should be called.
 
enum  
 Whether to assemble the pattern on the boundary intersections, i.e. whether or not pattern_boundary() should be called.
 

Flags for the non-constant part of the residual and the jacobian

enum  
 Whether to call the local operator's alpha_volume(), jacobian_apply_volume() and jacobian_volume().
 
enum  
 Whether to call the local operator's alpha_volume_post_skeleton(), jacobian_apply_volume_post_skeleton() and jacobian_volume_post_skeleton().
 
enum  
 Whether to call the local operator's alpha_skeleton(), jacobian_apply_skeleton() and jacobian_skeleton().
 
enum  
 Whether to call the local operator's alpha_boundary(), jacobian_apply_boundary() and jacobian_boundary().
 

Flags for the constant part of the residual

enum  
 Whether to call the local operator's lambda_volume().
 
enum  
 Whether to call the local operator's lambda_volume_post_skeleton().
 
enum  
 Whether to call the local operator's lambda_skeleton().
 
enum  
 Whether to call the local operator's lambda_boundary().
 

Special flags

enum  
 Whether to visit the skeleton methods from both sides.
 
enum  
 Wheter the local operator describes a linear problem.
 

Detailed Description

Class to document the stationary local operator interface.

This class is for documentation purposes only. Each method given here is controlled by a flag. If the corresponding flag for a method is false, that method is never called and it is permissible for the method to be missing entirely from the local operator. The flags are those from LocalOperatorDefaultFlags.

There are five categories of methods, denoted by the first part of their name:

  • pattern_*(): Methods for the sparsity pattern, controlled by the doPattern* flags,
  • alpha_*(): Methods for the non-constant parts of the residual, controlled by the doAlpha* flags,
  • lambda_*(): Methods for the constant parts of the residual, controlled by the doLambda* flags,
  • jacobian_apply_*(): Methods for the application of the jacobian, controlled by the doAlpha* flags, and finally
  • jacobian_*(): Methods to extract the jacobian, controlled by the doAlpha* flags.

There are four classes of methods, denoted by the last part of their name:

  • *_volume(): methods called on the entities before iterating over the intersections, controlled by the do*Volume flags,
  • *_volume_post_skeleton(): methods called on the entities after iterating over the intersections, controlled by the do*VolumePostSkeleton flags,
  • *_skeleton(): methods called on the interior intersections, controlled by the do*Skeleton flags, and finally
  • *_boundary(): methods called on the boundary intersections, controlled by the do*Boundary flags.

Not all combinations of categories and methods do actually exist.

To assemble the global sparsity pattern, residual or jacobian, the GridOperator iterates over the elements of the grid. For each element, it will call the appropriate *_volume() method. Then it will iterate through the elements intersections and call the appropriate *_skeleton() or *_boundary() methods on the intersection. Finally it will call the appropriate *_volume_post_skeleton() method.

The special flag doSkeletonTwoSided controls whether each interior intersection is visited once or twice. If it is true, each intersection may be given to *_skeleton() twice – the second time with the meaning of inside and outside exchanged. Note the "may": In the parallel case only interior entities are visited, so intersections at a processor boundary will only be visited once per processor in any case.

If doSkeletonTwoSided is false (the default), each intersection will only be visited once – and the orientation in which it is visited is left unspecified, except that the inside entity will always be an interior entity. Note that it will be visited once per process, so intersecions at processor boundaries are still visited twice when all processes are considered.

The alpha and lambda categories are a bit special in that the GridOperator uses them together – each time a method on the local operator should be called, it will first call the alpha_*() method and then call the lambda_*() method.

If the controlling flag for a method is false, the call to the method is omitted in such a way that the method does not even have to be present on the local operator. If both the do*Skeleton and do*Boundary flags are false, the iteration through the intersections is skipped.

Member Function Documentation

◆ alpha_boundary()

template<typename IG , typename LFSU , typename X , typename LFSV , typename R >
void Dune::PDELab::LocalOperatorInterface::alpha_boundary ( const IG &  ig,
const LFSU &  lfsu_s,
const X &  x_s,
const LFSV &  lfsv_s,
R &  r_s 
) const
inline

get a boundary intersections's contribution to alpha

Parameters
igIntersectionGeometry describing the intersection.
lfsu_sLocalFunctionSpace of the trial GridFunctionSpace in the inside entity.
x_sLocal position in the trial GridFunctionSpace in the inside entity.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
r_sLocal part of the residual in the inside entity.
Note
It is permissible to include contributions of the residual which are independent of x_s here (they have to be omitted from lambda_boundary() in that case, of course). This is the difference to jacobian_apply_boundary().
The method should not clear r_s; it should just add its entries to it.

This method is controlled by the flag doAlphaBoundary. For a given element, it's calls happen intermingled with the calls to alpha_skeleton(), but after the call to alpha_volume() and before the call to alpha_volume_post_skeleton().

◆ alpha_skeleton()

template<typename IG , typename LFSU , typename X , typename LFSV , typename R >
void Dune::PDELab::LocalOperatorInterface::alpha_skeleton ( const IG &  ig,
const LFSU &  lfsu_s,
const X &  x_s,
const LFSV &  lfsv_s,
const LFSU &  lfsu_n,
const X &  x_n,
const LFSV &  lfsv_n,
R &  r_s,
R &  r_n 
) const
inline

get an internal intersections's contribution to alpha

Parameters
igIntersectionGeometry describing the intersection.
lfsu_sLocalFunctionSpace of the trial GridFunctionSpace in the inside entity.
x_sLocal position in the trial GridFunctionSpace in the inside entity.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
lfsu_nLocalFunctionSpace of the trial GridFunctionSpace in the outside entity.
x_nLocal position in the trial GridFunctionSpace in the outside entity.
lfsv_nLocalFunctionSpace of the test GridFunctionSpace in the outside entity.
r_sLocal part of the residual in the inside entity.
r_nLocal part of the residual in the outside entity.
Note
It is permissible to include contributions of the residual which are independent of x_s and x_n here (they have to be omitted from lambda_skeleton() in that case, of course). This is the difference to jacobian_apply_skeleton().
The method should not clear r_s and r_n; it should just add its entries to them.

This method is controlled by the flag doAlphaSkeleton. For a given element, it's calls happen intermingled with the calls to alpha_boundary(), but after the call to alpha_volume() and before the call to alpha_volume_post_skeleton().

◆ alpha_volume()

template<typename EG , typename LFSU , typename X , typename LFSV , typename R >
void Dune::PDELab::LocalOperatorInterface::alpha_volume ( const EG &  eg,
const LFSU &  lfsu,
const X &  x,
const LFSV &  lfsv,
R &  r 
) const
inline

get an element's contribution to alpha

Parameters
egElementGeometry describing the entity.
lfsuLocalFunctionSpace of the trial GridFunctionSpace.
xLocal position in the trial GridFunctionSpace.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
rLocal part of the residual.
Note
It is permissible to include contributions of the residual which are independent of x here (they have to be omitted from lambda_volume() in that case, of course). This is the difference to jacobian_apply_volume().
x and r are of type std::vector.
The method should not clear r; it should just add its entries to it.

This method is controlled by the flag doAlphaVolume. For a given element, it is called before the alpha_skeleton() and/or alpha_boundary() methods are called (if they are called at all).

◆ alpha_volume_post_skeleton()

template<typename EG , typename LFSU , typename X , typename LFSV , typename R >
void Dune::PDELab::LocalOperatorInterface::alpha_volume_post_skeleton ( const EG &  eg,
const LFSU &  lfsu,
const X &  x,
const LFSV &  lfsv,
R &  r 
) const
inline

get an element's contribution to alpha after the intersections have been handled

Parameters
egElementGeometry describing the entity.
lfsuLocalFunctionSpace of the trial GridFunctionSpace.
xLocal position in the trial GridFunctionSpace.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
rLocal part of the residual.
Note
It is permissible to include contributions of the residual which are independent of x here (they have to be omitted from lambda_volume_post_skeleton() in that case, of course). This is the difference to jacobian_apply_volume_post_skeleton().
The method should not clear r; it should just add its entries to it.

This method is controlled by the flag doAlphaVolumePostSkeleton. For a given element, it is called after the alpha_skeleton() and/or alpha_boundary() methods are called (if they are called at all).

◆ jacobian_apply_boundary() [1/2]

template<typename IG , typename LFSU , typename X , typename Z , typename LFSV , typename Y >
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_boundary ( const IG &  ig,
const LFSU &  lfsu_s,
const X &  x_s,
const Z &  z_s,
const LFSV &  lfsv_s,
Y &  y_s 
) const
inline

apply a boundary intersections's jacobian for a nonlinear problem.

This method performs the following update:

\[ y = y + (\nabla_u \alpha_{\text{bnd}})\Big|_{x_s} z, \]

where \(u\) is the solution and the Jacobian is evaluated at \(u = x_s\).

Parameters
igIntersectionGeometry describing the intersection.
lfsu_sLocalFunctionSpace of the trial GridFunctionSpace in the inside entity.
x_sLocal position in the trial space for the inside entity at which to evaluate the Jacobian.
z_sLocal position in the trial space for the inside entity to which to apply the Jacobian.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
y_sLocal part of the residual in the inside entity.
Note
This is different from alpha_boundary(), since the result will be linear in x, whereas alpha_boundary() may include contributions to the the residual which are constant in x.
The method should not clear y_s; it should just add its entries to it.
This method assumes that the residual is not linear and the Jacobian is not constant; the Jacobian should be evaluated for the vector x and applied to the vector z.

This method is controlled by the flag doAlphaBoundary. For a given element, it's calls happen intermingled with the calls to jacobian_apply_skeleton(), but after the call to jacobian_apply_volume() and before the call to jacobian_apply_volume_post_skeleton().

◆ jacobian_apply_boundary() [2/2]

template<typename IG , typename LFSU , typename X , typename LFSV , typename Y >
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_boundary ( const IG &  ig,
const LFSU &  lfsu_s,
const X &  z_s,
const LFSV &  lfsv_s,
Y &  y_s 
) const
inline

apply a boundary intersections's jacobian for a linear problem.

This method performs the following update:

\[ y = y + (\nabla_u \alpha_{\text{bnd}}) z, \]

where \(u\) is the solution and the residual is considered linear in \(u\).

Parameters
igIntersectionGeometry describing the intersection.
lfsu_sLocalFunctionSpace of the trial GridFunctionSpace in the inside entity.
z_sLocal position in the trial space for the inside entity to which to apply the Jacobian.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
y_sLocal part of the residual in the inside entity.
Note
This is different from alpha_boundary(), since the result will be linear in x, whereas alpha_boundary() may include contributions to the the residual which are constant in x.
The method should not clear y_s; it should just add its entries to it.
This method assumes that the residual is linear and the Jacobian is constant; as such, the method isn't passed a linearization point at which to evaluate the Jacobian, but only the vector z_s to which the Jacobian should be applied.

This method is controlled by the flag doAlphaBoundary. For a given element, it's calls happen intermingled with the calls to jacobian_apply_skeleton(), but after the call to jacobian_apply_volume() and before the call to jacobian_apply_volume_post_skeleton().

◆ jacobian_apply_skeleton() [1/2]

template<typename IG , typename LFSU , typename X , typename Z , typename LFSV , typename Y >
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_skeleton ( const IG &  ig,
const LFSU &  lfsu_s,
const X &  x_s,
const Z &  z_s,
const LFSV &  lfsv_s,
const LFSU &  lfsu_n,
const X &  x_n,
const Z &  z_n,
const LFSV &  lfsv_n,
Y &  y_s,
Y &  y_n 
) const
inline

Applies an internal intersections's jacobians to the appropriate vectors for a nonlinear problem. This method performs the following updates:

\begin{align*} y_s & = y_s + (\nabla_{u_s} \alpha_{\text{skel},s})\Big|_{(x_s,x_n)} z_s + (\nabla_{u_n} \alpha_{\text{skel},s})\Big|_{(x_s,x_n)} z_n, \\ y_n & = y_n + (\nabla_{u_s} \alpha_{\text{skel},n})\Big|_{(x_s,x_n)} z_s + (\nabla_{u_n} \alpha_{\text{skel},n})\Big|_{(x_s,x_n)} z_n, \end{align*}

where \(u_s\) and \(u_n\) are the solution on the inside and the outside of the intersection and the Jacobian is evaluated at \(u = (x_s,x_n)\)$.

Parameters
igIntersectionGeometry describing the intersection.
lfsu_sLocalFunctionSpace of the trial GridFunctionSpace in the inside entity.
x_sLocal position in the trial space for the inside entity at which to evaluate the Jacobian.
z_sLocal position in the trial space for the inside entity to which to apply the Jacobian.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
lfsu_nLocalFunctionSpace of the trial GridFunctionSpace in the outside entity.
x_nLocal position in the trial space for the outside entity at which to evaluate the Jacobian.
z_nLocal position in the trial space for the outside entity to which to apply the Jacobian.
lfsv_nLocalFunctionSpace of the test GridFunctionSpace in the outside entity.
y_sWhere to store the inside entity's result.
y_nWhere to store the outside entity's result.
Note
This is different from alpha_skeleton(), since the result will be linear in z_s and z_n, whereas alpha_skeleton() may include contributions to the the residual which are constant in z_s and z_n.
The method should not clear y_s and y_n; it should just add its entries to them.
This method assumes that the residual is not linear and the Jacobian is not constant; the Jacobian should be evaluated for the vector x and applied to the vector z.

This method is controlled by the flag doAlphaSkeleton. For a given element, it's calls happen intermingled with the calls to jacobian_apply_boundary(), but after the call to jacobian_apply_volume() and before the call to jacobian_apply_volume_post_skeleton().

◆ jacobian_apply_skeleton() [2/2]

template<typename IG , typename LFSU , typename X , typename LFSV , typename Y >
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_skeleton ( const IG &  ig,
const LFSU &  lfsu_s,
const X &  z_s,
const LFSV &  lfsv_s,
const LFSU &  lfsu_n,
const X &  z_n,
const LFSV &  lfsv_n,
Y &  y_s,
Y &  y_n 
) const
inline

Applies an internal intersections's jacobians to the appropriate vectors for a linear problem. This method performs the following updates:

\begin{align*} y_s & = y_s + (\nabla_{u_s} \alpha_{\text{skel},s}) z_s + (\nabla_{u_n} \alpha_{\text{skel},s}) z_n, \\ y_n & = y_n + (\nabla_{u_s} \alpha_{\text{skel},n}) z_s + (\nabla_{u_n} \alpha_{\text{skel},n}) z_n, \end{align*}

where \(u_s\) and \(u_n\) are the solution on the inside and outside of the intersection and the residual is considered linear in \(u_{\{s,n\}}\).

Parameters
igIntersectionGeometry describing the intersection.
lfsu_sLocalFunctionSpace of the trial GridFunctionSpace in the inside entity.
z_sLocal position in the trial space for the inside entity to which to apply the Jacobian.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
lfsu_nLocalFunctionSpace of the trial GridFunctionSpace in the outside entity.
z_nLocal position in the trial space for the outside entity to which to apply the Jacobian.
lfsv_nLocalFunctionSpace of the test GridFunctionSpace in the outside entity.
y_sWhere to store the inside entity's result.
y_nWhere to store the outside entity's result.
Note
This is different from alpha_skeleton(), since the result will be linear in z_s and z_n, whereas alpha_skeleton() may include contributions to the the residual which are constant in z_s and z_n.
The method should not clear y_s and y_n; it should just add its entries to them.
This method assumes that the residual is linear and the Jacobian is constant; as such, the method isn't passed linearization points at which to evaluate the Jacobian, but only the vectors z_s and z_n to which the Jacobian should be applied.

This method is controlled by the flag doAlphaSkeleton. For a given element, it's calls happen intermingled with the calls to jacobian_apply_boundary(), but after the call to jacobian_apply_volume() and before the call to jacobian_apply_volume_post_skeleton().

◆ jacobian_apply_volume() [1/2]

template<typename EG , typename LFSU , typename X , typename Z , typename LFSV , typename Y >
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_volume ( const EG &  eg,
const LFSU &  lfsu,
const X &  x,
const Z &  z,
const LFSV &  lfsv,
Y &  y 
) const
inline

Applies an element's jacobian to a vector for a nonlinear problem.

This method performs the following update:

\[ y = y + (\nabla_u \alpha_{\text{vol}})\Big|_x z, \]

where \(u\) is the solution and the Jacobian is evaluated at \(u=x\).

Parameters
egElementGeometry describing the entity.
lfsuLocalFunctionSpace of the trial GridFunctionSpace.
xLocal position in the trial space at which to evaluate the Jacobian.
zLocal position in the trial space to which to apply the Jacobian.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
yWhere to store the result.
Note
This is different from alpha_volume(), since the result will be linear in z, whereas alpha_volume() may include contributions to the the residual which are constant in z.
The method should not clear y; it should just add its entries to it.
This method assumes that the residual is not linear and the Jacobian is not constant; the Jacobian should be evaluated for the vector x and applied to the vector z.
Your linear operator should export the flag isLinear=false in this case, see LocalOperatorDefaultFlags.

This method is controlled by the flag doAlphaVolume. For a given element, it is called before the jacobian_apply_skeleton() and/or jacobian_apply_boundary() methods are called (if they are called at all).

◆ jacobian_apply_volume() [2/2]

template<typename EG , typename LFSU , typename X , typename LFSV , typename Y >
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_volume ( const EG &  eg,
const LFSU &  lfsu,
const X &  z,
const LFSV &  lfsv,
Y &  y 
) const
inline

Applies an element's jacobian to a vector for a linear problem.

This method performs the following update:

\[ y = y + (\nabla_u \alpha_{\text{vol}}) z, \]

where \(u\) is the solution and the residual is considered linear in \(u\).

Parameters
egElementGeometry describing the entity.
lfsuLocalFunctionSpace of the trial GridFunctionSpace.
zLocal position in the trial space to which to apply the Jacobian.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
yWhere to store the result.
Note
This is different from alpha_volume(), since the result will be linear in z, whereas alpha_volume() may include contributions to the the residual which are constant in z.
The method should not clear y; it should just add its entries to it.
This method assumes that the residual is linear and the Jacobian is constant; as such, the method isn't passed a linearization point at which to evaluate the Jacobian, but only the vector z to which the Jacobian should be applied.
Your operator should export the flag isLinear=true in this case. This is done by default.

This method is controlled by the flag doAlphaVolume. For a given element, it is called before the jacobian_apply_skeleton() and/or jacobian_apply_boundary() methods are called (if they are called at all).

◆ jacobian_apply_volume_post_skeleton() [1/2]

template<typename EG , typename LFSU , typename X , typename Z , typename LFSV , typename Y >
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_volume_post_skeleton ( const EG &  eg,
const LFSU &  lfsu,
const X &  x,
const Z &  z,
const LFSV &  lfsv,
Y &  y 
) const
inline

Applies an element's jacobian to a vector for a nonlinear problem after the intersections have been handled. This method performs the following update:

\[ y = y + (\nabla_u \alpha_{\text{vol\_post\_skel}})\Big|_x z, \]

where \(u\) is the solution and the Jacobian is evaluated at \(u=x\).

Parameters
egElementGeometry describing the entity.
lfsuLocalFunctionSpace of the trial GridFunctionSpace.
xLocal position in the trial space at which to evaluate the Jacobian.
zLocal position in the trial space to which to apply the Jacobian.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
yWhere to store the result.
Note
This is different from alpha_volume_post_skeleton(), since the result will be linear in z, whereas alpha_volume_post_skeleton() may include contributions to the the residual which are constant in z.
The method should not clear y; it should just add its entries to it.
This method assumes that the residual is not linear and the Jacobian is not constant; the Jacobian should be evaluated for the vector x and applied to the vector z.

This method is controlled by the flag doAlphaVolumePostSkeleton. For a given element, it is called after the jacobian_apply_skeleton() and/or jacobian_apply_boundary() methods are called (if they are called at all).

◆ jacobian_apply_volume_post_skeleton() [2/2]

template<typename EG , typename LFSU , typename X , typename LFSV , typename Y >
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_volume_post_skeleton ( const EG &  eg,
const LFSU &  lfsu,
const X &  z,
const LFSV &  lfsv,
Y &  y 
) const
inline

Applies an element's jacobian to a vector for a linear problem after the intersections have been handled. This method performs the following update:

\[ y = y + (\nabla_u \alpha_{\text{vol\_post\_skel}}) z, \]

where \(u\) is the solution and the residual is considered linear in \(u\).

Parameters
egElementGeometry describing the entity.
lfsuLocalFunctionSpace of the trial GridFunctionSpace.
zLocal position in the trial space to which to apply the Jacobian.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
yWhere to store the result.
Note
This is different from alpha_volume_post_skeleton(), since the result will be linear in z, whereas alpha_volume_post_skeleton() may include contributions to the the residual which are constant in z.
The method should not clear y; it should just add its entries to it.
This method assumes that the residual is linear and the Jacobian is constant; as such, the method isn't passed a linearization point at which to evaluate the Jacobian, but only the vector z to which the Jacobian should be applied.

This method is controlled by the flag doAlphaVolumePostSkeleton. For a given element, it is called after the jacobian_apply_skeleton() and/or jacobian_apply_boundary() methods are called (if they are called at all).

◆ jacobian_boundary()

template<typename IG , typename LFSU , typename X , typename LFSV , typename LocalMatrix >
void Dune::PDELab::LocalOperatorInterface::jacobian_boundary ( const IG &  ig,
const LFSU &  lfsu_s,
const X &  x_s,
const LFSV &  lfsv_s,
LocalMatrix mat_ss 
) const
inline

get a boundary intersections's jacobian

Parameters
igIntersectionGeometry describing the intersection.
lfsu_sLocalFunctionSpace of the trial GridFunctionSpace in the inside entity.
x_sLocal position in the trial GridFunctionSpace in the inside entity.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
mat_ssWhere to store the contribution to the inside entity's jacobian.
Note
The method should not clear mat_ss; it should just add its entries to it.

This method is controlled by the flag doAlphaBoundary. For a given element, it's calls happen intermingled with the calls to jacobian_skeleton(), but after the call to jacobian_volume() and before the call to jacobian_volume_post_skeleton().

◆ jacobian_skeleton()

template<typename IG , typename LFSU , typename X , typename LFSV , typename LocalMatrix >
void Dune::PDELab::LocalOperatorInterface::jacobian_skeleton ( const IG &  ig,
const LFSU &  lfsu_s,
const X &  x_s,
const LFSV &  lfsv_s,
const LFSU &  lfsu_n,
const X &  x_n,
const LFSV &  lfsv_n,
LocalMatrix mat_ss,
LocalMatrix mat_sn,
LocalMatrix mat_ns,
LocalMatrix mat_nn 
) const
inline

apply an internal intersections's jacobians

Parameters
igIntersectionGeometry describing the intersection.
lfsu_sLocalFunctionSpace of the trial GridFunctionSpace in the inside entity.
x_sLocal position in the trial GridFunctionSpace in the inside entity.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
lfsu_nLocalFunctionSpace of the trial GridFunctionSpace in the outside entity.
x_nLocal position in the trial GridFunctionSpace in the outside entity.
lfsv_nLocalFunctionSpace of the test GridFunctionSpace in the outside entity.
mat_ssWhere to store the contribution to the inside entity's jacobian.
mat_snWhere to store the contribution to the interaction jacobian between the inside and the outside entity.
mat_nsWhere to store the contribution to the interaction jacobian between the outside and the inside entity.
mat_nnWhere to store the contribution to the outside entity's jacobian.
Note
The method should not clear mat_ss, mat_sn, mat_ns, mat_nn; it should just add its entries to them.

This method is controlled by the flag doAlphaSkeleton. For a given element, it's calls happen intermingled with the calls to jacobian_boundary(), but after the call to jacobian_volume() and before the call to jacobian_volume_post_skeleton().

◆ jacobian_volume()

template<typename EG , typename LFSU , typename X , typename LFSV , typename LocalMatrix >
void Dune::PDELab::LocalOperatorInterface::jacobian_volume ( const EG &  eg,
const LFSU &  lfsu,
const X &  x,
const LFSV &  lfsv,
LocalMatrix mat 
) const
inline

get an element's jacobian

Parameters
egElementGeometry describing the entity.
lfsuLocalFunctionSpace of the trial GridFunctionSpace.
xLocal position in the trial GridFunctionSpace.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
matWhere to store the contribution to the jacobian.
Note
The method should not clear mat; it should just add its entries to it.

This method is controlled by the flag doAlphaVolume. For a given element, it is called before the jacobian_skeleton() and/or jacobian_boundary() methods are called (if they are called at all).

◆ jacobian_volume_post_skeleton()

template<typename EG , typename LFSU , typename X , typename LFSV , typename LocalMatrix >
void Dune::PDELab::LocalOperatorInterface::jacobian_volume_post_skeleton ( const EG &  eg,
const LFSU &  lfsu,
const X &  x,
const LFSV &  lfsv,
LocalMatrix mat 
) const
inline

get an element's jacobian after the intersections have been handled

Parameters
egElementGeometry describing the entity.
lfsuLocalFunctionSpace of the trial GridFunctionSpace.
xLocal position in the trial GridFunctionSpace.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
matWhere to store the contribution to the jacobian.
Note
The method should not clear mat; it should just add its entries to it.

This method is controlled by the flag doAlphaVolumePostSkeleton. For a given element, it is called after the jacobian_skeleton() and/or jacobian_boundary() methods are called (if they are called at all).

◆ lambda_boundary()

template<typename IG , typename LFSV , typename R >
void Dune::PDELab::LocalOperatorInterface::lambda_boundary ( const IG &  ig,
const LFSV &  lfsv_s,
R &  r_s 
) const
inline

get a boundary intersections's contribution to lambda

Parameters
igIntersectionGeometry describing the intersection. inside entity.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
r_sLocal part of the residual in the inside entity.
Note
The method should not clear r_s; it should just add its entries to it.

This method is controlled by the flag doLambdaBoundary. For a given element, it's calls happen intermingled with the calls to lambda_skeleton(), but after the call to lambda_volume() and before the call to lambda_volume_post_skeleton().

◆ lambda_skeleton()

template<typename IG , typename LFSV , typename R >
void Dune::PDELab::LocalOperatorInterface::lambda_skeleton ( const IG &  ig,
const LFSV &  lfsv_s,
const LFSV &  lfsv_n,
R &  r_s,
R &  r_n 
) const
inline

get an internal intersections's contribution to lambda

Parameters
igIntersectionGeometry describing the intersection. inside entity.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
lfsv_nLocalFunctionSpace of the test GridFunctionSpace in the outside entity.
r_sLocal part of the residual in the inside entity.
r_nLocal part of the residual in the outside entity.
Note
The method should not clear r_s and r_n; it should just add its entries to them.

This method is controlled by the flag doLambdaSkeleton. For a given element, it's calls happen intermingled with the calls to lambda_boundary(), but after the call to lambda_volume() and before the call to lambda_volume_post_skeleton().

◆ lambda_volume()

template<typename EG , typename LFSV , typename R >
void Dune::PDELab::LocalOperatorInterface::lambda_volume ( const EG &  eg,
const LFSV &  lfsv,
R &  r 
) const
inline

get an element's contribution to lambda

Parameters
egElementGeometry describing the entity.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
rLocal part of the residual.
Note
The method should not clear r; it should just add its entries to it.

This method is controlled by the flag doLambdaVolume. For a given element, it is called before the lambda_skeleton() and/or lambda_boundary() methods are called (if they are called at all).

◆ lambda_volume_post_skeleton()

template<typename EG , typename LFSV , typename R >
void Dune::PDELab::LocalOperatorInterface::lambda_volume_post_skeleton ( const EG &  eg,
const LFSV &  lfsv,
R &  r 
) const
inline

get an element's contribution to lambda after the intersections have been handled

Parameters
egElementGeometry describing the entity.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
rLocal part of the residual.
Note
The method should not clear r; it should just add its entries to it.

This method is controlled by the flag doLambdaVolumePostSkeleton. For a given element, it is called after the lambda_skeleton() and/or lambda_boundary() methods are called (if they are called at all).

◆ pattern_boundary()

template<typename LFSU , typename LFSV , typename LocalPattern >
void Dune::PDELab::LocalOperatorInterface::pattern_boundary ( const LFSU &  lfsu_s,
const LFSV &  lfsv_s,
LocalPattern &  pattern_ss 
) const
inline

get a boundary intersection's contribution to the sparsity pattern

Parameters
lfsu_sLocalFunctionSpace of the trial GridFunctionSpace in the inside entity.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
pattern_ssLocal sparsity pattern for the inside entity.
Note
The method should not clear the pattern; it should just add its entries to it.

This method is controlled by the flag doPatternBoundary. For a given element, it's calls happen intermingled with the calls to pattern_skeleton(), but after the call to pattern_volume() and before the call to pattern_volume_post_skeleton().

◆ pattern_skeleton()

template<typename LFSU , typename LFSV , typename LocalPattern >
void Dune::PDELab::LocalOperatorInterface::pattern_skeleton ( const LFSU &  lfsu_s,
const LFSV &  lfsv_s,
const LFSU &  lfsu_n,
const LFSV &  lfsv_n,
LocalPattern &  pattern_sn,
LocalPattern &  pattern_ns 
) const
inline

get an internal intersection's contribution to the sparsity pattern

Parameters
lfsu_sLocalFunctionSpace of the trial GridFunctionSpace in the inside entity.
lfsv_sLocalFunctionSpace of the test GridFunctionSpace in the inside entity.
lfsu_nLocalFunctionSpace of the trial GridFunctionSpace in the outside entity.
lfsv_nLocalFunctionSpace of the test GridFunctionSpace in the outside entity.
pattern_snLocal sparsity pattern.
pattern_nsLocal sparsity pattern.
Note
The method should not clear the patterns; it should just add its entries to them.

This method is controlled by the flag doPatternSkeleton. For a given element, it's calls happen intermingled with the calls to pattern_boundary(), but after the call to pattern_volume() and before the call to pattern_volume_post_skeleton().

◆ pattern_volume()

template<typename LFSU , typename LFSV , typename LocalPattern >
void Dune::PDELab::LocalOperatorInterface::pattern_volume ( const LFSU &  lfsu,
const LFSV &  lfsv,
LocalPattern &  pattern 
) const
inline

get an element's contribution to the sparsity pattern

Parameters
lfsuLocalFunctionSpace of the trial GridFunctionSpace.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
patternLocal sparsity pattern.
Note
The method should not clear the pattern; it should just add its entries to it.

This method is controlled by the flag doPatternVolume. For a given element, it is called before the pattern_skeleton() and/or pattern_boundary() methods are called (if they are called at all).

◆ pattern_volume_post_skeleton()

template<typename LFSU , typename LFSV , typename LocalPattern >
void Dune::PDELab::LocalOperatorInterface::pattern_volume_post_skeleton ( const LFSU &  lfsu,
const LFSV &  lfsv,
LocalPattern &  pattern 
) const
inline

get an element's contribution to the sparsity pattern after the intersections have been handled

Parameters
lfsuLocalFunctionSpace of the trial GridFunctionSpace.
lfsvLocalFunctionSpace of the test GridFunctionSpace.
patternLocal sparsity pattern.
Note
The method should not clear the pattern; it should just add its entries to it.

This method is controlled by the flag doPatternVolume. For a given element, it is called before the pattern_skeleton() and/or pattern_boundary() methods are called (if they are called at all).

◆ skip_entity()

template<typename EG >
bool Dune::PDELab::LocalOperatorInterface::skip_entity ( const EG &  eg) const
inline

whether to assembly methods associated with a given entity

Parameters
egElementGeometry describing the entity.

This method is controlled by the flag doSkipEntity. For a given element, it is called before the any other local method. If return value is true, all volume, skeleton, and boundary methods are skipped.

◆ skip_intersection()

template<typename IG >
bool Dune::PDELab::LocalOperatorInterface::skip_intersection ( const IG &  ig) const
inline

whether to assembly methods associated with a given intersection

Parameters
igIntersectionGeometry describing the intersection.

This method is controlled by the flag doSkipIntersection. For a given intersection, it is called after local volume methods and before any skeleton and boundary methods. If return value is true, all skeleton, and boundary methods are skipped.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)