#include <dune/pdelab/localoperator/laplace.hh>
|
| Laplace (unsigned int quadOrder) |
| Constructor. More...
|
|
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 |
| Compute Laplace matrix times a given vector for one element. More...
|
|
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 &matrix) const |
| Compute the Laplace stiffness matrix for the element given in 'eg'. More...
|
|
template<typename LFSU , typename LFSV , typename LocalPattern > |
void | pattern_volume (const LFSU &lfsu, const LFSV &lfsv, LocalPattern &pattern) const |
|
a local operator for solving the Laplace equation
with conforming finite elements on all types of grids in any dimension.
In other words, it only assembles the Laplace matrix.
Whether to assemble the pattern on the elements, i.e. whether or not pattern_volume() should be called.
Enumerator |
---|
doPatternVolume |
|
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.
Enumerator |
---|
doPatternVolumePostSkeleton |
|
Whether to assemble the pattern on the interior intersections, i.e. whether or not pattern_skeleton() should be called.
Enumerator |
---|
doPatternSkeleton |
|
Whether to assemble the pattern on the boundary intersections, i.e. whether or not pattern_boundary() should be called.
Enumerator |
---|
doPatternBoundary |
|
Whether to call the local operator's alpha_volume(), jacobian_apply_volume() and jacobian_volume().
Whether to call the local operator's alpha_volume_post_skeleton(), jacobian_apply_volume_post_skeleton() and jacobian_volume_post_skeleton().
Enumerator |
---|
doAlphaVolumePostSkeleton |
|
Whether to call the local operator's alpha_skeleton(), jacobian_apply_skeleton() and jacobian_skeleton().
Enumerator |
---|
doAlphaSkeleton |
|
Whether to call the local operator's alpha_boundary(), jacobian_apply_boundary() and jacobian_boundary().
Enumerator |
---|
doAlphaBoundary |
|
Whether to call the local operator's lambda_volume().
Whether to call the local operator's lambda_volume_post_skeleton().
Enumerator |
---|
doLambdaVolumePostSkeleton |
|
Whether to call the local operator's lambda_skeleton().
Enumerator |
---|
doLambdaSkeleton |
|
Whether to call the local operator's lambda_boundary().
Enumerator |
---|
doLambdaBoundary |
|
Whether to visit the skeleton methods from both sides.
Enumerator |
---|
doSkeletonTwoSided |
|
Enumerator |
---|
doPatternVolume |
|
Dune::PDELab::Laplace::Laplace |
( |
unsigned int |
quadOrder | ) |
|
|
inline |
Constructor.
- Parameters
-
quadOrder | Order of the quadrature rule used for integrating over the element |
- Deprecated:
- "Deprecated in Dune-PDELab 2.4, use the local operator ConvectionDiffusionFEM instead!"
template<typename EG , typename LFSU , typename X , typename LFSV , typename R >
void Dune::PDELab::Laplace::alpha_volume |
( |
const EG & |
eg, |
|
|
const LFSU & |
lfsu, |
|
|
const X & |
x, |
|
|
const LFSV & |
lfsv, |
|
|
R & |
r |
|
) |
| const |
|
inline |
Compute Laplace matrix times a given vector for one element.
This is used for matrix-free algorithms for the Laplace equation
- Parameters
-
[in] | eg | The grid element we are assembling on |
[in] | lfsu | Local ansatz function space basis |
[in] | lfsv | Local test function space basis |
[in] | x | Input vector |
[out] | r | The product of the Laplace matrix times x |
References quadOrder_.
template<typename EG , typename LFSU , typename X , typename LFSV , typename M >
void Dune::PDELab::Laplace::jacobian_volume |
( |
const EG & |
eg, |
|
|
const LFSU & |
lfsu, |
|
|
const X & |
x, |
|
|
const LFSV & |
lfsv, |
|
|
M & |
matrix |
|
) |
| const |
|
inline |
Compute the Laplace stiffness matrix for the element given in 'eg'.
- Template Parameters
-
M | Type of the element stiffness matrix |
- Parameters
-
[in] | eg | The grid element we are assembling on |
[in] | lfsu | Local ansatz function space basis |
[in] | lfsv | Local test function space basis |
[in] | x | Current configuration; gets ignored for linear problems like this one |
[out] | matrix | Element stiffness matrix |
References dim, and quadOrder_.
template<typename LFSU , typename LFSV , typename LocalPattern >
void Dune::PDELab::FullVolumePattern::pattern_volume |
( |
const LFSU & |
lfsu, |
|
|
const LFSV & |
lfsv, |
|
|
LocalPattern & |
pattern |
|
) |
| const |
|
inlineinherited |
unsigned int Dune::PDELab::Laplace::quadOrder_ |
|
protected |
The documentation for this class was generated from the following file: