1#ifndef DUNE_FEM_ConservationLaw_MODEL_HH
2#define DUNE_FEM_ConservationLaw_MODEL_HH
7#include <dune/common/visibility.hh>
14#include <dune/fempy/quadrature/fempyquadratures.hh>
16#define VirtualConservationLawModelMethods(POINT) \
17 virtual void source ( const POINT &x,\
18 const DRangeType &value,\
19 const DJacobianRangeType &gradient,\
20 RRangeType &flux ) const = 0;\
21 virtual void linSource ( const DRangeType& uBar,\
22 const DJacobianRangeType &gradientBar,\
24 const DRangeType &value,\
25 const DJacobianRangeType &gradient,\
26 RRangeType &flux ) const = 0;\
27 virtual void flux ( const POINT &x,\
28 const DRangeType &value,\
29 const DJacobianRangeType &gradient,\
30 RJacobianRangeType &result ) const = 0;\
32 virtual void diffusiveFlux ( const POINT &x,\
33 const DRangeType &value,\
34 const DJacobianRangeType &gradient,\
35 RJacobianRangeType &result ) const \
36 { flux(x, value, gradient, result ); } \
37 virtual void linFlux ( const DRangeType& uBar,\
38 const DJacobianRangeType& gradientBar,\
40 const DRangeType &value,\
41 const DJacobianRangeType &gradient,\
42 RJacobianRangeType &flux ) const = 0;\
44 virtual void linDiffusiveFlux ( const DRangeType& uBar,\
45 const DJacobianRangeType& gradientBar,\
47 const DRangeType &value,\
48 const DJacobianRangeType &gradient,\
49 RJacobianRangeType &flux ) const \
50 { linFlux(uBar, gradientBar, x, value, gradient, flux ); } \
51 virtual void fluxDivergence( const POINT &x,\
52 const DRangeType &value,\
53 const DJacobianRangeType &jacobian,\
54 const DHessianRangeType &hessian,\
55 RRangeType &flux) const = 0;\
56 virtual void alpha(const POINT &x,\
57 const DRangeType &value,\
58 RRangeType &val) const = 0;\
59 virtual void linAlpha(const DRangeType &uBar,\
61 const DRangeType &value,\
62 RRangeType &val) const = 0;\
63 virtual void dirichlet( int bndId, const POINT &x,\
64 RRangeType &value) const = 0;
66#define WrapperConservationLawModelMethods(POINT) \
67 virtual void source ( const POINT &x,\
68 const DRangeType &value,\
69 const DJacobianRangeType &gradient,\
70 RRangeType &flux ) const \
71 { impl().source(x, value, gradient, flux); } \
72 virtual void linSource ( const DRangeType& uBar,\
73 const DJacobianRangeType &gradientBar,\
75 const DRangeType &value,\
76 const DJacobianRangeType &gradient,\
77 RRangeType &flux ) const \
78 { impl().linSource(uBar, gradientBar, x, value, gradient, flux); } \
79 virtual void flux ( const POINT &x,\
80 const DRangeType &value,\
81 const DJacobianRangeType &gradient,\
82 RJacobianRangeType &flux ) const \
83 { impl().flux(x, value, gradient, flux); } \
84 virtual void linFlux ( const DRangeType& uBar,\
85 const DJacobianRangeType& gradientBar,\
87 const DRangeType &value,\
88 const DJacobianRangeType &gradient,\
89 RJacobianRangeType &flux ) const \
90 { impl().linFlux(uBar, gradientBar, x, value, gradient, flux); } \
91 virtual void fluxDivergence( const POINT &x,\
92 const DRangeType &value,\
93 const DJacobianRangeType &jacobian,\
94 const DHessianRangeType &hessian,\
95 RRangeType &flux) const \
96 { impl().fluxDivergence(x, value, jacobian, hessian, flux); } \
97 virtual void alpha(const POINT &x,\
98 const DRangeType &value,\
99 RRangeType &val) const \
100 { impl().alpha(x, value, val); } \
101 virtual void linAlpha(const DRangeType &uBar,\
103 const DRangeType &value,\
104 RRangeType &val) const \
105 { impl().linAlpha(uBar, x, value, val); } \
106 virtual void dirichlet( int bndId, const POINT &x,\
107 RRangeType &value) const \
108 { impl().dirichlet(bndId, x, value); }
110template<
class Gr
idPart,
int dimDomain,
int dimRange=dimDomain,
class RangeField =
double >
114 static const int dimD = dimDomain;
115 static const int dimR = dimRange;
133 typedef typename GridPartType::template Codim<0>::EntityType
EntityType;
138 template <
class F,
int d>
139 using Traits = Dune::FemPy::FempyQuadratureTraits<F,d>;
169 virtual std::string
name()
const = 0;
176 virtual void setTime(
const double t )
const = 0;
180 virtual double time()
const = 0;
196 typedef std::array<int, dimRange> DirichletComponentType;
205 class CheckTimeMethod
208 static std::true_type testSignature(
double& (T::*)());
211 static decltype(testSignature(&T::time)) test(std::nullptr_t);
214 static std::false_type test(...);
216 using type =
decltype(test<M>(
nullptr));
218 static const bool value = type::value;
222 template <
class Model,
bool>
225 static void setTime( Model&,
const double ) {}
226 static double time(
const Model& ) {
return -1.0; }
229 template <
class Model>
230 struct CallSetTime< Model, true >
232 static void setTime( Model& model,
const double t ) { model.time() = t; }
233 static double time (
const Model& model ) {
return model.time(); }
238template <
class ModelImpl >
241 ModelImpl::dimD, ModelImpl::dimR, typename ModelImpl::RRangeFieldType>
243 typedef ModelImpl Impl;
245 static const int dimD = ModelImpl::dimD;
246 static const int dimR = ModelImpl::dimR;
268 template<
class... Args, std::enable_if_t< std::is_constructible< ModelImpl, Args &&... >::value,
int > = 0 >
270 :
impl_(
std::forward< Args >( args )... )
292 virtual std::string
name()
const
294 return impl().name();
301 detail::CallSetTime< ModelImpl, detail::CheckTimeMethod< ModelImpl >::value >
309 return detail::CallSetTime< ModelImpl, detail::CheckTimeMethod< ModelImpl >::value >
::time(
impl() );
315 return impl().hasDirichletBoundary();
319 return impl().hasNeumanBoundary();
323 return impl().isDirichletIntersection(inter, dirichletComponent);
327 return impl().init(entity);
331 return impl().unbind();
345#define VirtualPenaltyMethods(POINT) \
346 virtual RRangeType penalty ( const POINT &x,\
347 const DRangeType &value ) const = 0;\
348 virtual RRangeType linPenalty ( const POINT &x,\
349 const DRangeType &value ) const = 0;
350#define WrapperPenaltyMethods(POINT) \
351 virtual RRangeType penalty ( const POINT &x,\
352 const DRangeType &value ) const \
353 { return impl().penalty(x, value); } \
354 virtual RRangeType linPenalty ( const POINT &x,\
355 const DRangeType &value ) const \
356 { return impl().linPenalty( x, value); }
358template<
class Gr
idPart,
int dimDomain,
int dimRange=dimDomain,
class RangeField =
double >
362 static const int dimD = dimDomain;
363 static const int dimR = dimRange;
381 typedef typename GridPartType::template Codim<0>::EntityType
EntityType;
386 template <
class F,
int d>
387 using Traits = Dune::FemPy::FempyQuadratureTraits<F,d>;
417 virtual std::string
name()
const = 0;
424 virtual void setTime(
const double t )
const = 0;
428 virtual double time()
const = 0;
444 typedef std::array<int, dimRange> DirichletComponentType;
464template <
class ModelImpl >
466:
public DGConservationLawModel<typename ModelImpl::GridPartType, ModelImpl::dimD, ModelImpl::dimR, typename ModelImpl::RRangeFieldType>
468 typedef ModelImpl Impl;
470 static const int dimD = ModelImpl::dimD;
471 static const int dimR = ModelImpl::dimR;
493 template<
class... Args, std::enable_if_t< std::is_constructible< ModelImpl, Args &&... >::value,
int > = 0 >
495 :
impl_(
std::forward< Args >( args )... )
529 virtual std::string
name()
const
531 return impl().name();
538 detail::CallSetTime< ModelImpl, detail::CheckTimeMethod< ModelImpl >::value >
546 return detail::CallSetTime< ModelImpl, detail::CheckTimeMethod< ModelImpl >::value >
::time(
impl() );
552 return impl().hasDirichletBoundary();
556 return impl().hasNeumanBoundary();
560 return impl().isDirichletIntersection(inter, dirichletComponent);
564 return impl().init(entity);
568 return impl().unbind();
#define VirtualConservationLawModelMethods(POINT)
Definition: conservationlawmodel.hh:16
#define WrapperPenaltyMethods(POINT)
Definition: conservationlawmodel.hh:350
#define VirtualPenaltyMethods(POINT)
Definition: conservationlawmodel.hh:345
Definition: explicitfieldvector.hh:75
quadrature class supporting base function caching
Definition: cachingquadrature.hh:41
quadrature on the codim-0 reference element
Definition: elementquadrature.hh:58
Definition: conservationlawmodel.hh:112
GridPartType::template Codim< 0 >::EntityType EntityType
Definition: conservationlawmodel.hh:133
virtual double time() const =0
DFunctionSpaceType::RangeType DRangeType
Definition: conservationlawmodel.hh:124
Dune::Fem::FunctionSpace< double, RangeFieldType, GridPart::dimensionworld, dimR > RFunctionSpaceType
Definition: conservationlawmodel.hh:122
DFunctionSpaceType::DomainFieldType DDomainFieldType
Definition: conservationlawmodel.hh:127
EntityType::Geometry::LocalCoordinate LocalDomainType
Definition: conservationlawmodel.hh:135
virtual bool isDirichletIntersection(const IntersectionType &inter, DirichletComponentType &dirichletComponent) const =0
virtual std::string name() const =0
DFunctionSpaceType::HessianRangeType DHessianRangeType
Definition: conservationlawmodel.hh:126
virtual void unbind() const =0
ConservationLawModel()
Definition: conservationlawmodel.hh:165
DFunctionSpaceType::JacobianRangeType DJacobianRangeType
Definition: conservationlawmodel.hh:125
virtual VirtualConservationLawModelMethods(Point) VirtualConservationLawModelMethods(ElementPoint) VirtualConservationLawModelMethods(IntersectionPoint) VirtualConservationLawModelMethods(ElementIntersectionPoint) VirtualConservationLawModelMethods(OriginalPoint) VirtualConservationLawModelMethods(OriginalElementPoint) VirtualConservationLawModelMethods(OriginalIntersectionPoint) VirtualConservationLawModelMethods(OriginalElementIntersectionPoint) VirtualConservationLawModelMethods(LocalDomainType) typedef std bool hasDirichletBoundary() const =0
virtual bool hasNeumanBoundary() const =0
virtual ~ConservationLawModel()
Definition: conservationlawmodel.hh:167
DFunctionSpaceType::DomainType DomainType
Definition: conservationlawmodel.hh:123
RangeField RangeFieldType
Definition: conservationlawmodel.hh:117
Dune::Fem::CachingQuadrature< GridPartType, 0, Traits >::QuadraturePointWrapperType Point
Definition: conservationlawmodel.hh:141
virtual bool init(const EntityType &entity) const =0
virtual void setTime(const double t) const =0
RFunctionSpaceType::DomainFieldType rDomainFieldType
Definition: conservationlawmodel.hh:131
static const int dimD
Definition: conservationlawmodel.hh:114
Dune::Fem::CachingQuadrature< GridPartType, 0 >::QuadraturePointWrapperType OriginalPoint
Definition: conservationlawmodel.hh:151
Dune::Fem::FunctionSpace< double, RangeFieldType, GridPart::dimensionworld, dimD > DFunctionSpaceType
Definition: conservationlawmodel.hh:120
GridPartType::IntersectionType IntersectionType
Definition: conservationlawmodel.hh:134
Dune::FemPy::FempyQuadratureTraits< F, d > Traits
Definition: conservationlawmodel.hh:139
GridPart GridPartType
Definition: conservationlawmodel.hh:113
Dune::Fem::ElementQuadrature< GridPartType, 0, Traits >::QuadraturePointWrapperType ElementPoint
Definition: conservationlawmodel.hh:145
RFunctionSpaceType::RangeType RRangeType
Definition: conservationlawmodel.hh:128
Dune::Fem::CachingQuadrature< GridPartType, 1, Traits >::QuadraturePointWrapperType IntersectionPoint
Definition: conservationlawmodel.hh:143
Dune::Fem::CachingQuadrature< GridPartType, 1 >::QuadraturePointWrapperType OriginalIntersectionPoint
Definition: conservationlawmodel.hh:153
ConservationLawModel< GridPartType, dimD, dimR, RangeField > ModelType
Definition: conservationlawmodel.hh:116
Dune::Fem::ElementQuadrature< GridPartType, 1 >::QuadraturePointWrapperType OriginalElementIntersectionPoint
Definition: conservationlawmodel.hh:157
Dune::Fem::ElementQuadrature< GridPartType, 0 >::QuadraturePointWrapperType OriginalElementPoint
Definition: conservationlawmodel.hh:155
static const int dimR
Definition: conservationlawmodel.hh:115
RFunctionSpaceType::HessianRangeType RHessianRangeType
Definition: conservationlawmodel.hh:130
Dune::Fem::ElementQuadrature< GridPartType, 1, Traits >::QuadraturePointWrapperType ElementIntersectionPoint
Definition: conservationlawmodel.hh:147
RFunctionSpaceType::JacobianRangeType RJacobianRangeType
Definition: conservationlawmodel.hh:129
Definition: conservationlawmodel.hh:242
WrapperConservationLawModelMethods(IntersectionPoint)
virtual bool hasNeumanBoundary() const
Definition: conservationlawmodel.hh:317
Base::RJacobianRangeType RJacobianRangeType
Definition: conservationlawmodel.hh:263
ConservationLawModel< GridPartType, dimD, dimR, typename ModelImpl::RRangeFieldType > Base
Definition: conservationlawmodel.hh:247
Base::DomainType DomainType
Definition: conservationlawmodel.hh:258
Base::OriginalPoint OriginalPoint
Definition: conservationlawmodel.hh:253
std::array< int, dimR > DirichletComponentType
Definition: conservationlawmodel.hh:312
virtual bool isDirichletIntersection(const IntersectionType &inter, DirichletComponentType &dirichletComponent) const
Definition: conservationlawmodel.hh:321
Base::IntersectionType IntersectionType
Definition: conservationlawmodel.hh:266
Base::LocalDomainType LocalDomainType
Definition: conservationlawmodel.hh:257
WrapperConservationLawModelMethods(Point)
Base::RHessianRangeType RHessianRangeType
Definition: conservationlawmodel.hh:264
Base::IntersectionPoint IntersectionPoint
Definition: conservationlawmodel.hh:250
virtual double time() const
Definition: conservationlawmodel.hh:307
ModelImpl & impl()
Definition: conservationlawmodel.hh:337
Base::EntityType EntityType
Definition: conservationlawmodel.hh:265
~ConservationLawModelWrapper()
Definition: conservationlawmodel.hh:273
static const int dimD
Definition: conservationlawmodel.hh:245
WrapperConservationLawModelMethods(OriginalPoint)
WrapperConservationLawModelMethods(OriginalElementPoint)
virtual void setTime(const double t) const
Definition: conservationlawmodel.hh:299
Base::OriginalElementIntersectionPoint OriginalElementIntersectionPoint
Definition: conservationlawmodel.hh:256
ConservationLawModelWrapper(Args &&... args)
Definition: conservationlawmodel.hh:269
Base::DRangeType DRangeType
Definition: conservationlawmodel.hh:259
WrapperConservationLawModelMethods(OriginalIntersectionPoint)
Base::ElementPoint ElementPoint
Definition: conservationlawmodel.hh:251
static const int dimR
Definition: conservationlawmodel.hh:246
Base::OriginalElementPoint OriginalElementPoint
Definition: conservationlawmodel.hh:255
WrapperConservationLawModelMethods(OriginalElementIntersectionPoint)
virtual void unbind() const
Definition: conservationlawmodel.hh:329
Base::ElementIntersectionPoint ElementIntersectionPoint
Definition: conservationlawmodel.hh:252
const ModelImpl & impl() const
Definition: conservationlawmodel.hh:333
ModelImpl::GridPartType GridPartType
Definition: conservationlawmodel.hh:244
WrapperConservationLawModelMethods(ElementPoint)
WrapperConservationLawModelMethods(ElementIntersectionPoint)
virtual std::string name() const
Definition: conservationlawmodel.hh:292
Base::OriginalIntersectionPoint OriginalIntersectionPoint
Definition: conservationlawmodel.hh:254
virtual bool init(const EntityType &entity) const
Definition: conservationlawmodel.hh:325
Base::DJacobianRangeType DJacobianRangeType
Definition: conservationlawmodel.hh:260
Base::DHessianRangeType DHessianRangeType
Definition: conservationlawmodel.hh:261
Base::Point Point
Definition: conservationlawmodel.hh:249
Base::RRangeType RRangeType
Definition: conservationlawmodel.hh:262
virtual bool hasDirichletBoundary() const
Definition: conservationlawmodel.hh:313
WrapperConservationLawModelMethods(LocalDomainType)
ModelImpl impl_
Definition: conservationlawmodel.hh:342
Definition: conservationlawmodel.hh:360
RFunctionSpaceType::RangeType RRangeType
Definition: conservationlawmodel.hh:376
virtual bool hasNeumanBoundary() const =0
virtual double time() const =0
Dune::Fem::ElementQuadrature< GridPartType, 1 >::QuadraturePointWrapperType OriginalElementIntersectionPoint
Definition: conservationlawmodel.hh:405
RFunctionSpaceType::JacobianRangeType RJacobianRangeType
Definition: conservationlawmodel.hh:377
Dune::Fem::FunctionSpace< double, RangeFieldType, GridPart::dimensionworld, dimR > RFunctionSpaceType
Definition: conservationlawmodel.hh:370
RFunctionSpaceType::HessianRangeType RHessianRangeType
Definition: conservationlawmodel.hh:378
Dune::Fem::CachingQuadrature< GridPartType, 0, Traits >::QuadraturePointWrapperType Point
Definition: conservationlawmodel.hh:389
Dune::Fem::ElementQuadrature< GridPartType, 0 >::QuadraturePointWrapperType OriginalElementPoint
Definition: conservationlawmodel.hh:403
Dune::Fem::ElementQuadrature< GridPartType, 1, Traits >::QuadraturePointWrapperType ElementIntersectionPoint
Definition: conservationlawmodel.hh:395
virtual bool isDirichletIntersection(const IntersectionType &inter, DirichletComponentType &dirichletComponent) const =0
virtual void setTime(const double t) const =0
GridPartType::template Codim< 0 >::EntityType EntityType
Definition: conservationlawmodel.hh:381
DGConservationLawModel< GridPartType, dimD, dimR, RangeField > ModelType
Definition: conservationlawmodel.hh:364
EntityType::Geometry::LocalCoordinate LocalDomainType
Definition: conservationlawmodel.hh:383
static const int dimD
Definition: conservationlawmodel.hh:362
virtual VirtualConservationLawModelMethods(Point) VirtualConservationLawModelMethods(ElementPoint) VirtualConservationLawModelMethods(IntersectionPoint) VirtualConservationLawModelMethods(ElementIntersectionPoint) VirtualConservationLawModelMethods(OriginalPoint) VirtualConservationLawModelMethods(OriginalElementPoint) VirtualConservationLawModelMethods(OriginalIntersectionPoint) VirtualConservationLawModelMethods(OriginalElementIntersectionPoint) VirtualConservationLawModelMethods(LocalDomainType) typedef std bool hasDirichletBoundary() const =0
virtual std::string name() const =0
static const int dimR
Definition: conservationlawmodel.hh:363
Dune::Fem::CachingQuadrature< GridPartType, 0 >::QuadraturePointWrapperType OriginalPoint
Definition: conservationlawmodel.hh:399
DFunctionSpaceType::JacobianRangeType DJacobianRangeType
Definition: conservationlawmodel.hh:373
virtual bool init(const EntityType &entity) const =0
Dune::Fem::CachingQuadrature< GridPartType, 1, Traits >::QuadraturePointWrapperType IntersectionPoint
Definition: conservationlawmodel.hh:391
RFunctionSpaceType::DomainFieldType rDomainFieldType
Definition: conservationlawmodel.hh:379
DGConservationLawModel()
Definition: conservationlawmodel.hh:413
virtual ~DGConservationLawModel()
Definition: conservationlawmodel.hh:415
DFunctionSpaceType::HessianRangeType DHessianRangeType
Definition: conservationlawmodel.hh:374
GridPartType::IntersectionType IntersectionType
Definition: conservationlawmodel.hh:382
virtual void unbind() const =0
Dune::Fem::CachingQuadrature< GridPartType, 1 >::QuadraturePointWrapperType OriginalIntersectionPoint
Definition: conservationlawmodel.hh:401
Dune::FemPy::FempyQuadratureTraits< F, d > Traits
Definition: conservationlawmodel.hh:387
Dune::Fem::FunctionSpace< double, RangeFieldType, GridPart::dimensionworld, dimD > DFunctionSpaceType
Definition: conservationlawmodel.hh:368
RangeField RangeFieldType
Definition: conservationlawmodel.hh:365
GridPart GridPartType
Definition: conservationlawmodel.hh:361
DFunctionSpaceType::DomainFieldType DDomainFieldType
Definition: conservationlawmodel.hh:375
DFunctionSpaceType::RangeType DRangeType
Definition: conservationlawmodel.hh:372
DFunctionSpaceType::DomainType DomainType
Definition: conservationlawmodel.hh:371
Dune::Fem::ElementQuadrature< GridPartType, 0, Traits >::QuadraturePointWrapperType ElementPoint
Definition: conservationlawmodel.hh:393
Definition: conservationlawmodel.hh:467
Base::DJacobianRangeType DJacobianRangeType
Definition: conservationlawmodel.hh:485
ModelImpl::GridPartType GridPartType
Definition: conservationlawmodel.hh:469
const ModelImpl & impl() const
Definition: conservationlawmodel.hh:570
virtual bool hasDirichletBoundary() const
Definition: conservationlawmodel.hh:550
ModelImpl impl_
Definition: conservationlawmodel.hh:579
Base::DHessianRangeType DHessianRangeType
Definition: conservationlawmodel.hh:486
~DGConservationLawModelWrapper()
Definition: conservationlawmodel.hh:498
DGConservationLawModel< GridPartType, dimD, dimR, typename ModelImpl::RRangeFieldType > Base
Definition: conservationlawmodel.hh:472
static const int dimD
Definition: conservationlawmodel.hh:470
ModelImpl & impl()
Definition: conservationlawmodel.hh:574
Base::Point Point
Definition: conservationlawmodel.hh:474
Base::RHessianRangeType RHessianRangeType
Definition: conservationlawmodel.hh:489
Base::IntersectionType IntersectionType
Definition: conservationlawmodel.hh:491
WrapperConservationLawModelMethods(IntersectionPoint)
WrapperPenaltyMethods(LocalDomainType)
Base::OriginalPoint OriginalPoint
Definition: conservationlawmodel.hh:478
WrapperConservationLawModelMethods(OriginalIntersectionPoint)
virtual std::string name() const
Definition: conservationlawmodel.hh:529
virtual void unbind() const
Definition: conservationlawmodel.hh:566
Base::ElementPoint ElementPoint
Definition: conservationlawmodel.hh:476
virtual bool init(const EntityType &entity) const
Definition: conservationlawmodel.hh:562
WrapperConservationLawModelMethods(Point)
Base::EntityType EntityType
Definition: conservationlawmodel.hh:490
Base::OriginalElementPoint OriginalElementPoint
Definition: conservationlawmodel.hh:480
std::array< int, dimR > DirichletComponentType
Definition: conservationlawmodel.hh:549
Base::OriginalIntersectionPoint OriginalIntersectionPoint
Definition: conservationlawmodel.hh:479
virtual bool hasNeumanBoundary() const
Definition: conservationlawmodel.hh:554
Base::DomainType DomainType
Definition: conservationlawmodel.hh:483
WrapperConservationLawModelMethods(OriginalElementPoint)
WrapperConservationLawModelMethods(ElementPoint)
Base::OriginalElementIntersectionPoint OriginalElementIntersectionPoint
Definition: conservationlawmodel.hh:481
WrapperPenaltyMethods(OriginalElementPoint)
WrapperPenaltyMethods(Point) WrapperPenaltyMethods(ElementPoint) WrapperPenaltyMethods(IntersectionPoint) WrapperPenaltyMethods(ElementIntersectionPoint) WrapperPenaltyMethods(OriginalPoint)
Base::IntersectionPoint IntersectionPoint
Definition: conservationlawmodel.hh:475
WrapperConservationLawModelMethods(LocalDomainType)
DGConservationLawModelWrapper(Args &&... args)
Definition: conservationlawmodel.hh:494
Base::RRangeType RRangeType
Definition: conservationlawmodel.hh:487
virtual void setTime(const double t) const
Definition: conservationlawmodel.hh:536
WrapperConservationLawModelMethods(OriginalPoint)
static const int dimR
Definition: conservationlawmodel.hh:471
Base::RJacobianRangeType RJacobianRangeType
Definition: conservationlawmodel.hh:488
virtual bool isDirichletIntersection(const IntersectionType &inter, DirichletComponentType &dirichletComponent) const
Definition: conservationlawmodel.hh:558
WrapperPenaltyMethods(OriginalElementIntersectionPoint)
Base::ElementIntersectionPoint ElementIntersectionPoint
Definition: conservationlawmodel.hh:477
WrapperConservationLawModelMethods(ElementIntersectionPoint)
WrapperConservationLawModelMethods(OriginalElementIntersectionPoint)
Base::DRangeType DRangeType
Definition: conservationlawmodel.hh:484
virtual double time() const
Definition: conservationlawmodel.hh:544
Base::LocalDomainType LocalDomainType
Definition: conservationlawmodel.hh:482
WrapperPenaltyMethods(OriginalIntersectionPoint)
A vector valued function space.
Definition: functionspace.hh:60
FunctionSpaceTraits::DomainFieldType DomainFieldType
Intrinsic type used for values in the domain field (usually a double)
Definition: functionspaceinterface.hh:60
FunctionSpaceTraits::RangeType RangeType
Type of range vector (using type of range field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:71
FunctionSpaceTraits::LinearMappingType JacobianRangeType
Intrinsic type used for the jacobian values has a Dune::FieldMatrix type interface.
Definition: functionspaceinterface.hh:75
FunctionSpaceTraits::DomainType DomainType
Type of domain vector (using type of domain field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:67