DUNE PDELab (git)

Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time > Class Template Reference

A local operator that scales the result of another local operator. More...

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

Control flags

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.
 
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().
 
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().
 
enum  
 Whether to visit the skeleton methods from both sides.
 

Methods for temporal local operators

typedef Time RealType
 
void setTime (Time t)
 set time for subsequent evaluation More...
 
Time getTime () const
 get current time More...
 
void preStep (Time time, Time dt, int stages)
 to be called once before each time step More...
 
void postStep ()
 to be called once at the end of each time step More...
 
void preStage (Time time, int r)
 to be called once before each stage More...
 
int getStage () const
 get current stage More...
 
void postStage ()
 to be called once at the end of each stage
 
Time suggestTimestep (Time dt) const
 to be called after stage 1 More...
 

Construction and modification

 ScaledLocalOperator (Backend &backend, Factor factor_=0)
 construct a ScaledLocalOperator More...
 
 ScaledLocalOperator (Factor factor_=0)
 construct a ScaledLocalOperator More...
 
void setFactor (Factor factor_)
 set the scaling factor
 
Factor getFactor () const
 get the scaling factor
 
void setBackend (Backend &backend)
 set the backend
 
Backend & getBackend () const
 get a reference to the backend
 

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

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) const
 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) const
 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) const
 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) const
 apply a boundary intersections's jacobian More...
 

Methods to extract the jacobian

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

Detailed Description

template<typename Backend, typename Factor, typename Time = double>
class Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >

A local operator that scales the result of another local operator.

This local operator takes another local operator as a backend and a scaling factor. It forwards calls to the evaluation methods to the backend local operator, but scales the result by the scaling factor before returning it. Calls to the instationary methods setTime() et. al. are forwarded as well, without modification to the behaviour.

If the scaling factor equal zero, calls to the evaluation methods of the backend operator are eliminated at run time. Note that "equals zero" is determined by the expression factor==0, so if you have values that are nearly zero and you want them to be treated as zero, you should force the factor to zero yourself before passing it to this class.

Template Parameters
BackendType of the backend operator.
FactorType of the scaling factor.
TimeType of time values.

Constructor & Destructor Documentation

◆ ScaledLocalOperator() [1/2]

template<typename Backend , typename Factor , typename Time = double>
Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::ScaledLocalOperator ( Backend &  backend,
Factor  factor_ = 0 
)
inline

construct a ScaledLocalOperator

Parameters
backendReference to the backend local operator
factor_Scaling factor

◆ ScaledLocalOperator() [2/2]

template<typename Backend , typename Factor , typename Time = double>
Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::ScaledLocalOperator ( Factor  factor_ = 0)
inline

construct a ScaledLocalOperator

Note
This constructor does not set a backend local operator – you must call setBackend() before you can use this object
Parameters
factor_Scaling factor

Member Function Documentation

◆ alpha_boundary()

template<typename Backend , typename Factor , typename Time = double>
template<typename IG , typename LFSU , typename X , typename LFSV , typename R >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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().
x_s and r_s are of type std::vector.
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 Backend , typename Factor , typename Time = double>
template<typename IG , typename LFSU , typename X , typename LFSV , typename R >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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().
x_s, x_n, r_s and r_n are of type std::vector.
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 Backend , typename Factor , typename Time = double>
template<typename EG , typename LFSU , typename X , typename LFSV , typename R >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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 Backend , typename Factor , typename Time = double>
template<typename EG , typename LFSU , typename X , typename LFSV , typename R >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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().
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 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).

◆ getStage()

template<typename Backend , typename Factor , typename Time = double>
int Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::getStage ( ) const
inline

get current stage

Returns
The current stage number previously set by preStage().

◆ getTime()

template<typename Backend , typename Factor , typename Time = double>
Time Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::getTime ( ) const
inline

get current time

Returns
The time previously set by setTime().

◆ jacobian_apply_boundary()

template<typename Backend , typename Factor , typename Time = double>
template<typename IG , typename LFSU , typename X , typename LFSV , typename Y >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::jacobian_apply_boundary ( const IG &  ig,
const LFSU &  lfsu_s,
const X &  x_s,
const LFSV &  lfsv_s,
Y &  y_s 
) const
inline

apply 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.
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.
x_s and y_s are of type std::vector.
The method should not clear 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().

◆ jacobian_apply_skeleton()

template<typename Backend , typename Factor , typename Time = double>
template<typename IG , typename LFSU , typename X , typename LFSV , typename Y >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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 
) 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.
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 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.
The method should not clear 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().

◆ jacobian_apply_volume()

template<typename Backend , typename Factor , typename Time = double>
template<typename EG , typename LFSU , typename X , typename LFSV , typename Y >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::jacobian_apply_volume ( const EG &  eg,
const LFSU &  lfsu,
const X &  x,
const LFSV &  lfsv,
Y &  y 
) const
inline

apply 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.
yWhere to store the result.
Note
This is different from alpha_volume(), since the result will be linear in x, whereas alpha_volume() may include contributions to the the residual which are constant in x.
x and y are of type std::vector.
The method should not clear 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).

◆ jacobian_apply_volume_post_skeleton()

template<typename Backend , typename Factor , typename Time = double>
template<typename EG , typename LFSU , typename X , typename LFSV , typename Y >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::jacobian_apply_volume_post_skeleton ( const EG &  eg,
const LFSU &  lfsu,
const X &  x,
const LFSV &  lfsv,
Y &  y 
) const
inline

apply 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.
yWhere to store the result.
Note
This is different from alpha_volume_post_skeleton(), since the result will be linear in 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.
The method should not clear 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).

◆ jacobian_boundary()

template<typename Backend , typename Factor , typename Time = double>
template<typename IG , typename LFSU , typename X , typename LFSV , typename M >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::jacobian_boundary ( const IG &  ig,
const LFSU &  lfsu_s,
const X &  x_s,
const LFSV &  lfsv_s,
M &  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
x_s is of type std::vector.
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 Backend , typename Factor , typename Time = double>
template<typename IG , typename LFSU , typename X , typename LFSV , typename M >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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,
M &  mat_ss,
M &  mat_sn,
M &  mat_ns,
M &  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
x_s and x_n are of type std::vector.
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 Backend , typename Factor , typename Time = double>
template<typename EG , typename LFSU , typename X , typename LFSV , typename M >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::jacobian_volume ( const EG &  eg,
const LFSU &  lfsu,
const X &  x,
const LFSV &  lfsv,
M &  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
x is of type std::vector.
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 Backend , typename Factor , typename Time = double>
template<typename EG , typename LFSU , typename X , typename LFSV , typename M >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::jacobian_volume_post_skeleton ( const EG &  eg,
const LFSU &  lfsu,
const X &  x,
const LFSV &  lfsv,
M &  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
x is of type std::vector.
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 Backend , typename Factor , typename Time = double>
template<typename IG , typename LFSV , typename R >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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
r_s is of type std::vector.
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 Backend , typename Factor , typename Time = double>
template<typename IG , typename LFSV , typename R >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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
r_s and r_n are of type std::vector.
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 Backend , typename Factor , typename Time = double>
template<typename EG , typename LFSV , typename R >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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
r is 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 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 Backend , typename Factor , typename Time = double>
template<typename EG , typename LFSV , typename R >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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
r is 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 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 Backend , typename Factor , typename Time = double>
template<typename LFSU , typename LFSV , typename LocalPattern >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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 Backend , typename Factor , typename Time = double>
template<typename LFSU , typename LFSV , typename LocalPattern >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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 Backend , typename Factor , typename Time = double>
template<typename LFSU , typename LFSV , typename LocalPattern >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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 Backend , typename Factor , typename Time = double>
template<typename LFSU , typename LFSV , typename LocalPattern >
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::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).

◆ postStep()

template<typename Backend , typename Factor , typename Time = double>
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::postStep ( )
inline

to be called once at the end of each time step

Note
With the OneStepMethod and the ExplicitOneStepMetod, for reasons unknown this is only called for temporal but not for spatial local operators. With the MultiStepMethod this is called for all local operators.

◆ preStage()

template<typename Backend , typename Factor , typename Time = double>
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::preStage ( Time  time,
int  r 
)
inline

to be called once before each stage

Parameters
timeTime of the stage
rNumber of the stage, r ∈ [1, nstages] inclusive, where nstages is the number of stage in the step given in the previous call to preStep()
Note
For ExplicitOneStepMethod the time given here for stage 1 may be incorrect, since the time step size is only finally determined after the first stage has been assembled.
For the MultiStepMethod, this is called once after preStep() with r=1.

◆ preStep()

template<typename Backend , typename Factor , typename Time = double>
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::preStep ( Time  time,
Time  dt,
int  stages 
)
inline

to be called once before each time step

Parameters
timeTime at beginning of the step.
dtSize of time step.
stagesNumber of stages to do in the step. For the MultiStepMethod this is always 1.
Note
For ExplicitOneStepMethod the dt given here may be incorrect, since the time step size is only finally determined after the first stage has been assembled.
For the MultiStepMethod the number of stages is given as 1. There are no since there are no times of evaluation in the middle of the step, a multi-step method is similar to a one step method with one stage.

◆ setTime()

template<typename Backend , typename Factor , typename Time = double>
void Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::setTime ( Time  t)
inline

set time for subsequent evaluation

This method set the time for subsequent calls to the alpha_*(), lambda_*(), jacobian_*() and jacobian_apply_*() methods.

Note
For ExplicitOneStepMethod the time given here in the first stage may be incorrect, since the time step size is only finally determined after the first stage has been assembled.

◆ suggestTimestep()

template<typename Backend , typename Factor , typename Time = double>
Time Dune::PDELab::ScaledLocalOperator< Backend, Factor, Time >::suggestTimestep ( Time  dt) const
inline

to be called after stage 1

Note
Only used by the ExplicitOneStepMethod.

This may be called on the spatial local operator in the case of an explicit one step scheme. It is called after stage 1 has been assembled (so the time given to preStep() may not apply anymore in this case). All the alpha_*() and lambda_*() methods should have been called, so they are a good place to generate the information returned here.


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 (Nov 23, 23:29, 2024)