dune-fem  2.4.1-rc
Classes | Public Types | Public Member Functions | Protected Member Functions | List of all members
Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits > Struct Template Reference

Interface class for problem definition in the LDG context. More...

#include </local/tomalk/somewhere/tmp/dune-fem/dune/fem/pass/localdg/discretemodel.hh>

Inheritance diagram for Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >:
Inheritance graph

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

DGDiscreteModelTypeasImp ()
 
const DGDiscreteModelTypeasImp () const
 

Detailed Description

template<class DGDiscreteModelTraits>
struct Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >

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).

Note
The definition of f(u) differs from the usual definition for conservation laws in the leading sign.

Member Typedef Documentation

template<class DGDiscreteModelTraits >
typedef Traits::DGDiscreteModelType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::DGDiscreteModelType

Implementation type for Barton-Nackman trick.

template<class DGDiscreteModelTraits >
typedef Traits::DiscreteFunctionSpaceType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::DiscreteFunctionSpaceType

Discrete function space.

template<class DGDiscreteModelTraits >
typedef FunctionSpaceType::DomainType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::DomainType

Coordinate type (world coordinates)

template<class DGDiscreteModelTraits >
typedef DiscreteFunctionSpaceType::EntityType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::EntityType

Element (codim 0 entity) of the grid.

template<class DGDiscreteModelTraits >
typedef DiscreteFunctionSpaceType::FunctionSpaceType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::FunctionSpaceType

Function space type.

template<class DGDiscreteModelTraits >
typedef DiscreteFunctionSpaceType::GridPartType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::GridPartType

Type of GridPart.

template<class DGDiscreteModelTraits >
typedef DiscreteFunctionSpaceType::IntersectionType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::IntersectionType

Intersection type.

template<class DGDiscreteModelTraits >
typedef FunctionSpaceType::JacobianRangeType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::JacobianRangeType

Matrix type of the function space's jacobian.

template<class DGDiscreteModelTraits >
typedef EntityType::Geometry::LocalCoordinate Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::LocalCoordinateType

Type of local coordinate.

template<class DGDiscreteModelTraits >
typedef ExtractMatrix< JacobianRangeType >::Type Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::MassFactorType

Type of mass factor.

template<class DGDiscreteModelTraits >
typedef FunctionSpaceType::RangeFieldType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::RangeFieldType

Vector type of the function space's range field type.

template<class DGDiscreteModelTraits >
typedef FunctionSpaceType::RangeType Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::RangeType

Vector type of the function space's range.

template<class DGDiscreteModelTraits >
typedef DGDiscreteModelTraits Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::Traits

Traits class defined by the user.

Member Enumeration Documentation

template<class DGDiscreteModelTraits >
anonymous enum

dimRange

Enumerator
dimRange 

Member Function Documentation

template<class DGDiscreteModelTraits >
template<class ArgumentTuple >
void Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::analyticalFlux ( const EntityType en,
const double  time,
const LocalCoordinateType x,
const ArgumentTuple &  u,
JacobianRangeType f 
)
inline

Computes the analytical flux of the problem. Analytical flux of the problem.

Parameters
enEntity where the flux is to be evaluated
timeGlobal time.
xReference element coordinate where the flux is evaluated
uTuple of the states of the previous passes and of the global argument.
fThe analytical flux (return value)

References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().

template<class DGDiscreteModelTraits >
DGDiscreteModelType& Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp ( )
inlineprotected
template<class DGDiscreteModelTraits >
const DGDiscreteModelType& Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp ( ) const
inlineprotected
template<class DGDiscreteModelTraits >
template<class ArgumentTuple , class FaceDomainType >
double Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::boundaryFlux ( const IntersectionType intersection,
const double  time,
const FaceDomainType &  x,
const ArgumentTuple &  uLeft,
RangeType gLeft 
)
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.

Parameters
intersectionIntersection in consideration.
timeGlobal time.
xLocal coordinate of the point where the numerical flux gets evaluated.
uLeftTuple of the states on the inner side of the intersection.
gLeftThe numerical flux contribution to the inner element (return value).
Returns
Speed of the fastest wave. This information is used to compute the maximum admissible timestep size.

References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().

template<class DGDiscreteModelTraits >
bool Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::hasFlux ( ) const
inline
template<class DGDiscreteModelTraits >
bool Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::hasMass ( ) const
inline

Returns true if problem has a mass matrix factor disctinct from the identity.

References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().

template<class DGDiscreteModelTraits >
bool Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::hasSource ( ) const
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().

template<class DGDiscreteModelTraits >
template<class ArgumentTuple >
void Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::mass ( const EntityType en,
const double  time,
const LocalCoordinateType x,
const ArgumentTuple &  u,
MassFactorType m 
)
inline

Implements the mass factor term of the problem. Mass term contribution.

Parameters
enEntity where the mass is to be evaluated
timeGlobal time.
xReference element coordinate where the mass is evaluated
uTuple of the states of the previous passes and of the global
mThe mass contribution (return value). default implementation sets this factor to 1.0

References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().

template<class DGDiscreteModelTraits >
template<class ArgumentTuple , class FaceDomainType >
double Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::numericalFlux ( const IntersectionType intersection,
const double  time,
const FaceDomainType &  x,
const ArgumentTuple &  uLeft,
const ArgumentTuple &  uRight,
RangeType gLeft,
RangeType gRight 
)
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.

Parameters
intersectionThe intersection in consideration.
timeGlobal time.
xLocal coordinate of the point where the numerical flux gets evaluated.
uLeftTuple of the states on the inner side of the intersection.
uRightTuple of the states on the outer side of the intersection.
gLeftThe numerical flux contribution to the inner element (return value).
gRightThe numerical flux contribution to the outer element (return value).
Returns
Speed of the fastest wave. This information is used to compute the maximum admissible timestep size.

References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().

template<class DGDiscreteModelTraits >
void Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::setEntity ( EntityType en)
inline

Passes the active entity to the model. This can be used, to set local functions required as data function in the model.

Parameters
enactive Entity

References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().

template<class DGDiscreteModelTraits >
void Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::setNeighbor ( EntityType nb)
inline

Passes the active neigbor entity to the model. This can be used, to set local functions required as data functions in the model.

Parameters
nbactive neighbor Entity

References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().

template<class DGDiscreteModelTraits >
template<class ArgumentTuple , class JacobianTuple >
double Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::source ( const EntityType en,
const double  time,
const LocalCoordinateType x,
const ArgumentTuple &  u,
const JacobianTuple &  jac,
RangeType s 
)
inline

Implements the source term of the problem. Source term contribution. The source term can contain parts of the non-conservative contributions.

Parameters
enEntity where the flux is to be evaluated
timeGlobal time.
xReference element coordinate where the flux is evaluated
uTuple of the states of the previous passes and of the global argument.
jacTuple of the gradients of the precedings passes' results (needed for the non-conservative contributions)
sThe source contribution (return value).

References Dune::Fem::DGDiscreteModelInterface< DGDiscreteModelTraits >::asImp().


The documentation for this struct was generated from the following file: