dune-pdelab
2.4-dev
|
Class to document the stationary local operator interface. More...
#include <dune/pdelab/localoperator/interface.hh>
Flags for the sparsity pattern | |
enum | { doPatternVolume } |
Whether to assemble the pattern on the elements, i.e. whether or not pattern_volume() should be called. More... | |
enum | { doPatternVolumePostSkeleton } |
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. More... | |
enum | { doPatternSkeleton } |
Whether to assemble the pattern on the interior intersections, i.e. whether or not pattern_skeleton() should be called. More... | |
enum | { doPatternBoundary } |
Whether to assemble the pattern on the boundary intersections, i.e. whether or not pattern_boundary() should be called. More... | |
Flags for the non-constant part of the residual and the jacobian | |
enum | { doAlphaVolume } |
Whether to call the local operator's alpha_volume(), jacobian_apply_volume() and jacobian_volume(). More... | |
enum | { doAlphaVolumePostSkeleton } |
Whether to call the local operator's alpha_volume_post_skeleton(), jacobian_apply_volume_post_skeleton() and jacobian_volume_post_skeleton(). More... | |
enum | { doAlphaSkeleton } |
Whether to call the local operator's alpha_skeleton(), jacobian_apply_skeleton() and jacobian_skeleton(). More... | |
enum | { doAlphaBoundary } |
Whether to call the local operator's alpha_boundary(), jacobian_apply_boundary() and jacobian_boundary(). More... | |
Flags for the constant part of the residual | |
enum | { doLambdaVolume } |
Whether to call the local operator's lambda_volume(). More... | |
enum | { doLambdaVolumePostSkeleton } |
Whether to call the local operator's lambda_volume_post_skeleton(). More... | |
enum | { doLambdaSkeleton } |
Whether to call the local operator's lambda_skeleton(). More... | |
enum | { doLambdaBoundary } |
Whether to call the local operator's lambda_boundary(). More... | |
Special flags | |
enum | { doSkeletonTwoSided } |
Whether to visit the skeleton methods from both sides. More... | |
Methods for the sparsity pattern | |
template<typename LFSU , typename LFSV , typename LocalPattern > | |
void | pattern_volume (const LFSU &lfsu, const LFSV &lfsv, LocalPattern &pattern) |
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) |
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) |
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) |
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) |
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) |
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) |
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) |
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) |
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) |
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) |
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) |
get a boundary intersections's contribution to lambda More... | |
Methods for the application of the jacobian | |
template<typename EG , typename LFSU , typename X , typename LFSV , typename Y > | |
void | jacobian_apply_volume (const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, Y &y) |
apply an element's jacobian 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 &x, const LFSV &lfsv, Y &y) |
apply an element's jacobian after the intersections have been handled More... | |
template<typename IG , typename LFSU , typename X , typename LFSV , typename Y > | |
void | jacobian_apply_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, Y &y_s, Y &y_n) |
apply an internal intersections's jacobians More... | |
template<typename IG , typename LFSU , typename X , typename LFSV , typename Y > | |
void | jacobian_apply_boundary (const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, Y &y_s) |
apply a boundary intersections's jacobian 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) |
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) |
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) |
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) |
get a boundary intersections's jacobian More... | |
Class to document the stationary local operator interface.
This class is for documentation purposes only. Each method given here is controlled by an 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 bounary intersections, controlled by the do*Boundary
flags.Not all combinations of categories and methods to actually exits.
To assemble the global sparsity pattern, residual or jacobian, the GridOperatorSpace 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 is may be given to *_skeleton
() twice – the second time with the meaning of inside and outside exchanged. Note "may": In the paralell 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 GridOperatorSpace 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 has 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.
void Dune::PDELab::LocalOperatorInterface::alpha_boundary | ( | const IG & | ig, |
const LFSU & | lfsu_s, | ||
const X & | x_s, | ||
const LFSV & | lfsv_s, | ||
R & | r_s | ||
) |
get a boundary intersections's contribution to alpha
ig | IntersectionGeometry describing the intersection. |
lfsu_s | LocalFunctionSpace of the trial GridFunctionSpace in the inside entity. |
x_s | Local position in the trial GridFunctionSpace in the inside entity. |
lfsv_s | LocalFunctionSpace of the test GridFunctionSpace in the inside entity. |
r_s | Local part of the residual in the inside entity. |
x_s
here (they have to be omitted from lambda_boundary() in that case, of course). This is the difference to jacobian_apply_boundary().x_s
and r_s
are of type std::vector.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().
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 | ||
) |
get an internal intersections's contribution to alpha
ig | IntersectionGeometry describing the intersection. |
lfsu_s | LocalFunctionSpace of the trial GridFunctionSpace in the inside entity. |
x_s | Local position in the trial GridFunctionSpace in the inside entity. |
lfsv_s | LocalFunctionSpace of the test GridFunctionSpace in the inside entity. |
lfsu_n | LocalFunctionSpace of the trial GridFunctionSpace in the outside entity. |
x_n | Local position in the trial GridFunctionSpace in the outside entity. |
lfsv_n | LocalFunctionSpace of the test GridFunctionSpace in the outside entity. |
r_s | Local part of the residual in the inside entity. |
r_n | Local part of the residual in the outside entity. |
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().x_s
, x_n
, r_s
and r_n
are of type std::vector.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().
void Dune::PDELab::LocalOperatorInterface::alpha_volume | ( | const EG & | eg, |
const LFSU & | lfsu, | ||
const X & | x, | ||
const LFSV & | lfsv, | ||
R & | r | ||
) |
get an element's contribution to alpha
eg | ElementGeometry describing the entity. |
lfsu | LocalFunctionSpace of the trial GridFunctionSpace. |
x | Local position in the trial GridFunctionSpace. |
lfsv | LocalFunctionSpace of the test GridFunctionSpace. |
r | Local part of the residual. |
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.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).
void Dune::PDELab::LocalOperatorInterface::alpha_volume_post_skeleton | ( | const EG & | eg, |
const LFSU & | lfsu, | ||
const X & | x, | ||
const LFSV & | lfsv, | ||
R & | r | ||
) |
get an element's contribution to alpha after the intersections have been handled
eg | ElementGeometry describing the entity. |
lfsu | LocalFunctionSpace of the trial GridFunctionSpace. |
x | Local position in the trial GridFunctionSpace. |
lfsv | LocalFunctionSpace of the test GridFunctionSpace. |
r | Local part of the residual. |
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().x
and r
are of type std::vector.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).
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_boundary | ( | const IG & | ig, |
const LFSU & | lfsu_s, | ||
const X & | x_s, | ||
const LFSV & | lfsv_s, | ||
Y & | y_s | ||
) |
apply a boundary intersections's jacobian
ig | IntersectionGeometry describing the intersection. |
lfsu_s | LocalFunctionSpace of the trial GridFunctionSpace in the inside entity. |
x_s | Local position in the trial GridFunctionSpace in the inside entity. |
lfsv_s | LocalFunctionSpace of the test GridFunctionSpace in the inside entity. |
y_s | Local part of the residual in the inside entity. |
x
, whereas alpha_boundary() may include contributions to the the residual which are constant in x
.x_s
and y_s
are of type std::vector.y_s
; it should just add its entries to it.x_s
is both the position where the jacobian is evaluated (for non-linear problems) as well as the vector the jacobian is applied to.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().
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_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, | ||
Y & | y_s, | ||
Y & | y_n | ||
) |
apply an internal intersections's jacobians
ig | IntersectionGeometry describing the intersection. |
lfsu_s | LocalFunctionSpace of the trial GridFunctionSpace in the inside entity. |
x_s | Local position in the trial GridFunctionSpace in the inside entity. |
lfsv_s | LocalFunctionSpace of the test GridFunctionSpace in the inside entity. |
lfsu_n | LocalFunctionSpace of the trial GridFunctionSpace in the outside entity. |
x_n | Local position in the trial GridFunctionSpace in the outside entity. |
lfsv_n | LocalFunctionSpace of the test GridFunctionSpace in the outside entity. |
y_s | Where to store the inside entity's result. |
y_n | Where to store the outside entity's result. |
x_s
and x_n
, whereas alpha_skeleton() may include contributions to the the residual which are constant in x_s
and x_n
.x_s
, x_n
, y_s
and y_n
are of type std::vector.y_s
and y_n
; it should just add its entries to them.x_s
and x_n
are both the positions where the jacobian is evaluated (for non-linear problems) as well as the vectors the jacobian is applied to.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().
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_volume | ( | const EG & | eg, |
const LFSU & | lfsu, | ||
const X & | x, | ||
const LFSV & | lfsv, | ||
Y & | y | ||
) |
apply an element's jacobian
eg | ElementGeometry describing the entity. |
lfsu | LocalFunctionSpace of the trial GridFunctionSpace. |
x | Local position in the trial GridFunctionSpace. |
lfsv | LocalFunctionSpace of the test GridFunctionSpace. |
y | Where to store the result. |
x
, whereas alpha_volume() may include contributions to the the residual which are constant in x
.x
and y
are of type std::vector.y
; it should just add its entries to it.x
is both the position where the jacobian is evaluated (for non-linear problems) as well as the vector the jacobian is applied to.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).
void Dune::PDELab::LocalOperatorInterface::jacobian_apply_volume_post_skeleton | ( | const EG & | eg, |
const LFSU & | lfsu, | ||
const X & | x, | ||
const LFSV & | lfsv, | ||
Y & | y | ||
) |
apply an element's jacobian after the intersections have been handled
eg | ElementGeometry describing the entity. |
lfsu | LocalFunctionSpace of the trial GridFunctionSpace. |
x | Local position in the trial GridFunctionSpace. |
lfsv | LocalFunctionSpace of the test GridFunctionSpace. |
y | Where to store the result. |
x
, whereas alpha_volume_post_skeleton() may include contributions to the the residual which are constant in x
.x
and y
are of type std::vector.y
; it should just add its entries to it.x
is both the position where the jacobian is evaluated (for non-linear problems) as well as the vector the jacobian is applied to.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).
void Dune::PDELab::LocalOperatorInterface::jacobian_boundary | ( | const IG & | ig, |
const LFSU & | lfsu_s, | ||
const X & | x_s, | ||
const LFSV & | lfsv_s, | ||
LocalMatrix & | mat_ss | ||
) |
get a boundary intersections's jacobian
ig | IntersectionGeometry describing the intersection. |
lfsu_s | LocalFunctionSpace of the trial GridFunctionSpace in the inside entity. |
x_s | Local position in the trial GridFunctionSpace in the inside entity. |
lfsv_s | LocalFunctionSpace of the test GridFunctionSpace in the inside entity. |
mat_ss | Where to store the contribution to the inside entity's jacobian. |
x_s
is of type std::vector.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().
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 | ||
) |
apply an internal intersections's jacobians
ig | IntersectionGeometry describing the intersection. |
lfsu_s | LocalFunctionSpace of the trial GridFunctionSpace in the inside entity. |
x_s | Local position in the trial GridFunctionSpace in the inside entity. |
lfsv_s | LocalFunctionSpace of the test GridFunctionSpace in the inside entity. |
lfsu_n | LocalFunctionSpace of the trial GridFunctionSpace in the outside entity. |
x_n | Local position in the trial GridFunctionSpace in the outside entity. |
lfsv_n | LocalFunctionSpace of the test GridFunctionSpace in the outside entity. |
mat_ss | Where to store the contribution to the inside entity's jacobian. |
mat_sn | Where to store the contribution to the interaction jacobian between the inside and the outside entity. |
mat_ns | Where to store the contribution to the interaction jacobian between the outside and the inside entity. |
mat_nn | Where to store the contribution to the outside entity's jacobian. |
x_s
and x_n
are of type std::vector.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().
void Dune::PDELab::LocalOperatorInterface::jacobian_volume | ( | const EG & | eg, |
const LFSU & | lfsu, | ||
const X & | x, | ||
const LFSV & | lfsv, | ||
LocalMatrix & | mat | ||
) |
get an element's jacobian
eg | ElementGeometry describing the entity. |
lfsu | LocalFunctionSpace of the trial GridFunctionSpace. |
x | Local position in the trial GridFunctionSpace. |
lfsv | LocalFunctionSpace of the test GridFunctionSpace. |
mat | Where to store the contribution to the jacobian. |
x
is of type std::vector.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).
void Dune::PDELab::LocalOperatorInterface::jacobian_volume_post_skeleton | ( | const EG & | eg, |
const LFSU & | lfsu, | ||
const X & | x, | ||
const LFSV & | lfsv, | ||
LocalMatrix & | mat | ||
) |
get an element's jacobian after the intersections have been handled
eg | ElementGeometry describing the entity. |
lfsu | LocalFunctionSpace of the trial GridFunctionSpace. |
x | Local position in the trial GridFunctionSpace. |
lfsv | LocalFunctionSpace of the test GridFunctionSpace. |
mat | Where to store the contribution to the jacobian. |
x
is of type std::vector.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).
void Dune::PDELab::LocalOperatorInterface::lambda_boundary | ( | const IG & | ig, |
const LFSV & | lfsv_s, | ||
R & | r_s | ||
) |
get a boundary intersections's contribution to lambda
ig | IntersectionGeometry describing the intersection. inside entity. |
lfsv_s | LocalFunctionSpace of the test GridFunctionSpace in the inside entity. |
r_s | Local part of the residual in the inside entity. |
r_s
is of type std::vector.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().
void Dune::PDELab::LocalOperatorInterface::lambda_skeleton | ( | const IG & | ig, |
const LFSV & | lfsv_s, | ||
const LFSV & | lfsv_n, | ||
R & | r_s, | ||
R & | r_n | ||
) |
get an internal intersections's contribution to lambda
ig | IntersectionGeometry describing the intersection. inside entity. |
lfsv_s | LocalFunctionSpace of the test GridFunctionSpace in the inside entity. |
lfsv_n | LocalFunctionSpace of the test GridFunctionSpace in the outside entity. |
r_s | Local part of the residual in the inside entity. |
r_n | Local part of the residual in the outside entity. |
r_s
and r_n
are of type std::vector.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().
void Dune::PDELab::LocalOperatorInterface::lambda_volume | ( | const EG & | eg, |
const LFSV & | lfsv, | ||
R & | r | ||
) |
get an element's contribution to lambda
eg | ElementGeometry describing the entity. |
lfsv | LocalFunctionSpace of the test GridFunctionSpace. |
r | Local part of the residual. |
r
is of type std::vector.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).
void Dune::PDELab::LocalOperatorInterface::lambda_volume_post_skeleton | ( | const EG & | eg, |
const LFSV & | lfsv, | ||
R & | r | ||
) |
get an element's contribution to lambda after the intersections have been handled
eg | ElementGeometry describing the entity. |
lfsv | LocalFunctionSpace of the test GridFunctionSpace. |
r | Local part of the residual. |
r
is of type std::vector.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).
void Dune::PDELab::LocalOperatorInterface::pattern_boundary | ( | const LFSU & | lfsu_s, |
const LFSV & | lfsv_s, | ||
LocalPattern & | pattern_ss | ||
) |
get a boundary intersection's contribution to the sparsity pattern
lfsu_s | LocalFunctionSpace of the trial GridFunctionSpace in the inside entity. |
lfsv_s | LocalFunctionSpace of the test GridFunctionSpace in the inside entity. |
pattern_ss | Local sparsity pattern for the inside entity. |
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().
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 | ||
) |
get an internal intersection's contribution to the sparsity pattern
lfsu_s | LocalFunctionSpace of the trial GridFunctionSpace in the inside entity. |
lfsv_s | LocalFunctionSpace of the test GridFunctionSpace in the inside entity. |
lfsu_n | LocalFunctionSpace of the trial GridFunctionSpace in the outside entity. |
lfsv_n | LocalFunctionSpace of the test GridFunctionSpace in the outside entity. |
pattern_sn | Local sparsity pattern. |
pattern_ns | Local sparsity pattern. |
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().
void Dune::PDELab::LocalOperatorInterface::pattern_volume | ( | const LFSU & | lfsu, |
const LFSV & | lfsv, | ||
LocalPattern & | pattern | ||
) |
get an element's contribution to the sparsity pattern
lfsu | LocalFunctionSpace of the trial GridFunctionSpace. |
lfsv | LocalFunctionSpace of the test GridFunctionSpace. |
pattern | Local sparsity pattern. |
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).
void Dune::PDELab::LocalOperatorInterface::pattern_volume_post_skeleton | ( | const LFSU & | lfsu, |
const LFSV & | lfsv, | ||
LocalPattern & | pattern | ||
) |
get an element's contribution to the sparsity pattern after the intersections have been handled
lfsu | LocalFunctionSpace of the trial GridFunctionSpace. |
lfsv | LocalFunctionSpace of the test GridFunctionSpace. |
pattern | Local sparsity pattern. |
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).