dune-fem
2.4.1-rc
|
Interface class for problem definition in the LDG context. More...
#include </local/tomalk/somewhere/tmp/dune-fem/dune/fem/pass/localdg/discretemodel.hh>
Classes | |
struct | ExtractMatrix |
Mass factor type. More... | |
struct | ExtractMatrix< MatrixType< ctype, dimD, dimR > > |
Public Types | |
enum | { dimRange = DiscreteFunctionSpaceType::dimRange } |
dimRange More... | |
typedef DGDiscreteModelTraits | Traits |
Traits class defined by the user. More... | |
typedef Traits::DGDiscreteModelType | DGDiscreteModelType |
Implementation type for Barton-Nackman trick. More... | |
typedef Traits::DiscreteFunctionSpaceType | DiscreteFunctionSpaceType |
Discrete function space. More... | |
typedef DiscreteFunctionSpaceType::FunctionSpaceType | FunctionSpaceType |
Function space type. More... | |
typedef FunctionSpaceType::DomainType | DomainType |
Coordinate type (world coordinates) More... | |
typedef FunctionSpaceType::RangeType | RangeType |
Vector type of the function space's range. More... | |
typedef FunctionSpaceType::RangeFieldType | RangeFieldType |
Vector type of the function space's range field type. More... | |
typedef FunctionSpaceType::JacobianRangeType | JacobianRangeType |
Matrix type of the function space's jacobian. More... | |
typedef DiscreteFunctionSpaceType::GridPartType | GridPartType |
Type of GridPart. More... | |
typedef DiscreteFunctionSpaceType::EntityType | EntityType |
Element (codim 0 entity) of the grid. More... | |
typedef DiscreteFunctionSpaceType::IntersectionType | IntersectionType |
Intersection type. More... | |
typedef EntityType::Geometry::LocalCoordinate | LocalCoordinateType |
Type of local coordinate. More... | |
typedef ExtractMatrix< JacobianRangeType >::Type | MassFactorType |
Type of mass factor. More... | |
Public Member Functions | |
bool | hasFlux () const |
Returns true if problem has a flux contribution. If you program this method to return true, make sure to implement numericalFlux and analyticalFlux as well. More... | |
bool | hasSource () const |
Returns true if problem has a source term. If you program this method to return true, make sure to implement method source well. More... | |
bool | hasMass () const |
Returns true if problem has a mass matrix factor disctinct from the identity. More... | |
template<class ArgumentTuple , class FaceDomainType > | |
double | numericalFlux (const IntersectionType &intersection, const double time, const FaceDomainType &x, const ArgumentTuple &uLeft, const ArgumentTuple &uRight, RangeType &gLeft, RangeType &gRight) |
Computes the numerical flux at a cell interface. More... | |
template<class ArgumentTuple , class FaceDomainType > | |
double | boundaryFlux (const IntersectionType &intersection, const double time, const FaceDomainType &x, const ArgumentTuple &uLeft, RangeType &gLeft) |
Computes the flux at the boundary Special kind of numerical flux. The intersection iterator provides the necessary information to identify the type of boundary. More... | |
template<class ArgumentTuple > | |
void | analyticalFlux (const EntityType &en, const double time, const LocalCoordinateType &x, const ArgumentTuple &u, JacobianRangeType &f) |
Computes the analytical flux of the problem. Analytical flux of the problem. More... | |
template<class ArgumentTuple , class JacobianTuple > | |
double | source (const EntityType &en, const double time, const LocalCoordinateType &x, const ArgumentTuple &u, const JacobianTuple &jac, RangeType &s) |
Implements the source term of the problem. Source term contribution. The source term can contain parts of the non-conservative contributions. More... | |
template<class ArgumentTuple > | |
void | mass (const EntityType &en, const double time, const LocalCoordinateType &x, const ArgumentTuple &u, MassFactorType &m) |
Implements the mass factor term of the problem. Mass term contribution. More... | |
void | setEntity (EntityType &en) |
Passes the active entity to the model. This can be used, to set local functions required as data function in the model. More... | |
void | setNeighbor (EntityType &nb) |
Passes the active neigbor entity to the model. This can be used, to set local functions required as data functions in the model. More... | |
Protected Member Functions | |
DGDiscreteModelType & | asImp () |
const DGDiscreteModelType & | asImp () const |
Interface class for problem definition in the LDG context.
The problem interface prescribes the methods needed by the implementation of local DG passes. Users need to derive from either this class or the DGDiscreteModelDefault class. The methods provided by this class are used by LocalDGPass (and possibly other classes) to solve for the space discretization of the equation d_t u + div f(u) = s(u), where f represents a flux function and s a source term. Even non-conservative fluxes can be treated (for details see paper by Dedner et al).
typedef Traits::DGDiscreteModelType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::DGDiscreteModelType |
Implementation type for Barton-Nackman trick.
typedef Traits::DiscreteFunctionSpaceType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::DiscreteFunctionSpaceType |
Discrete function space.
typedef FunctionSpaceType::DomainType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::DomainType |
Coordinate type (world coordinates)
typedef DiscreteFunctionSpaceType::EntityType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::EntityType |
Element (codim 0 entity) of the grid.
typedef DiscreteFunctionSpaceType::FunctionSpaceType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::FunctionSpaceType |
Function space type.
typedef DiscreteFunctionSpaceType::GridPartType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::GridPartType |
Type of GridPart.
typedef DiscreteFunctionSpaceType::IntersectionType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::IntersectionType |
Intersection type.
typedef FunctionSpaceType::JacobianRangeType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::JacobianRangeType |
Matrix type of the function space's jacobian.
typedef EntityType::Geometry::LocalCoordinate Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::LocalCoordinateType |
Type of local coordinate.
typedef ExtractMatrix< JacobianRangeType >::Type Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::MassFactorType |
Type of mass factor.
typedef FunctionSpaceType::RangeFieldType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::RangeFieldType |
Vector type of the function space's range field type.
typedef FunctionSpaceType::RangeType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::RangeType |
Vector type of the function space's range.
typedef DGDiscreteModelTraits Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::Traits |
Traits class defined by the user.
|
inline |
Computes the analytical flux of the problem. Analytical flux of the problem.
en | Entity where the flux is to be evaluated |
time | Global time. |
x | Reference element coordinate where the flux is evaluated |
u | Tuple of the states of the previous passes and of the global argument. |
f | The analytical flux (return value) |
References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().
|
inlineprotected |
Referenced by Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::analyticalFlux(), Dune::Fem::DGDiscreteModelDefault< DGDiscreteModelTraits, N1, N2, N3, N4, N5, N6, N7, N8, N9 >::analyticalFlux(), Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::boundaryFlux(), Dune::Fem::DGDiscreteModelDefault< DGDiscreteModelTraits, N1, N2, N3, N4, N5, N6, N7, N8, N9 >::boundaryFlux(), Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::hasFlux(), Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::hasMass(), Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::hasSource(), Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::mass(), Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::numericalFlux(), Dune::Fem::DGDiscreteModelDefault< DGDiscreteModelTraits, N1, N2, N3, N4, N5, N6, N7, N8, N9 >::numericalFlux(), Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::setEntity(), Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::setNeighbor(), Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::source(), and Dune::Fem::DGDiscreteModelDefault< DGDiscreteModelTraits, N1, N2, N3, N4, N5, N6, N7, N8, N9 >::source().
|
inlineprotected |
|
inline |
Computes the flux at the boundary Special kind of numerical flux. The intersection iterator provides the necessary information to identify the type of boundary.
intersection | Intersection in consideration. |
time | Global time. |
x | Local coordinate of the point where the numerical flux gets evaluated. |
uLeft | Tuple of the states on the inner side of the intersection. |
gLeft | The numerical flux contribution to the inner element (return value). |
References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().
|
inline |
Returns true if problem has a flux contribution. If you program this method to return true, make sure to implement numericalFlux and analyticalFlux as well.
References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().
Referenced by Dune::Fem::DGDiscreteModelDefault< DGDiscreteModelTraits, N1, N2, N3, N4, N5, N6, N7, N8, N9 >::analyticalFlux(), Dune::Fem::DGDiscreteModelDefault< DGDiscreteModelTraits, N1, N2, N3, N4, N5, N6, N7, N8, N9 >::boundaryFlux(), and Dune::Fem::DGDiscreteModelDefault< DGDiscreteModelTraits, N1, N2, N3, N4, N5, N6, N7, N8, N9 >::numericalFlux().
|
inline |
Returns true if problem has a mass matrix factor disctinct from the identity.
References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().
|
inline |
Returns true if problem has a source term. If you program this method to return true, make sure to implement method source well.
References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().
Referenced by Dune::Fem::DGDiscreteModelDefault< DGDiscreteModelTraits, N1, N2, N3, N4, N5, N6, N7, N8, N9 >::source().
|
inline |
Implements the mass factor term of the problem. Mass term contribution.
en | Entity where the mass is to be evaluated |
time | Global time. |
x | Reference element coordinate where the mass is evaluated |
u | Tuple of the states of the previous passes and of the global |
m | The mass contribution (return value). default implementation sets this factor to 1.0 |
References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().
|
inline |
Computes the numerical flux at a cell interface.
Interface method for calculating the numerical flux at cell interfaces. The method has distinct return values for the contribution to the inner and outer element, since the framework allows for non-conservative terms. The non-conservative terms don't get evaluated over separate methods, but need to be implemented using the flux and source methods.
intersection | The intersection in consideration. |
time | Global time. |
x | Local coordinate of the point where the numerical flux gets evaluated. |
uLeft | Tuple of the states on the inner side of the intersection. |
uRight | Tuple of the states on the outer side of the intersection. |
gLeft | The numerical flux contribution to the inner element (return value). |
gRight | The numerical flux contribution to the outer element (return value). |
References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().
|
inline |
Passes the active entity to the model. This can be used, to set local functions required as data function in the model.
en | active Entity |
References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().
|
inline |
Passes the active neigbor entity to the model. This can be used, to set local functions required as data functions in the model.
nb | active neighbor Entity |
References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().
|
inline |
Implements the source term of the problem. Source term contribution. The source term can contain parts of the non-conservative contributions.
en | Entity where the flux is to be evaluated |
time | Global time. |
x | Reference element coordinate where the flux is evaluated |
u | Tuple of the states of the previous passes and of the global argument. |
jac | Tuple of the gradients of the precedings passes' results (needed for the non-conservative contributions) |
s | The source contribution (return value). |
References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().