1#ifndef ELLIPTC_MODEL_HH
2#define ELLIPTC_MODEL_HH
7#include <dune/common/visibility.hh>
14#include <dune/fempy/quadrature/fempyquadratures.hh>
16#define VirtualDiffusionModelMethods(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 diffusiveFlux ( const POINT &x,\
28 const DRangeType &value,\
29 const DJacobianRangeType &gradient,\
30 RJacobianRangeType &flux ) const = 0;\
31 virtual void linDiffusiveFlux ( const DRangeType& uBar,\
32 const DJacobianRangeType& gradientBar,\
34 const DRangeType &value,\
35 const DJacobianRangeType &gradient,\
36 RJacobianRangeType &flux ) const = 0;\
37 virtual void fluxDivergence( const POINT &x,\
38 const DRangeType &value,\
39 const DJacobianRangeType &jacobian,\
40 const DHessianRangeType &hessian,\
41 RRangeType &flux) const = 0;\
42 virtual void alpha(const POINT &x,\
43 const DRangeType &value,\
44 RRangeType &val) const = 0;\
45 virtual void linAlpha(const DRangeType &uBar,\
47 const DRangeType &value,\
48 RRangeType &val) const = 0;\
49 virtual void dirichlet( int bndId, const POINT &x,\
50 RRangeType &value) const = 0;
52#define WrapperDiffusionModelMethods(POINT) \
53 virtual void source ( const POINT &x,\
54 const DRangeType &value,\
55 const DJacobianRangeType &gradient,\
56 RRangeType &flux ) const \
57 { impl().source(x, value, gradient, flux); } \
58 virtual void linSource ( const DRangeType& uBar,\
59 const DJacobianRangeType &gradientBar,\
61 const DRangeType &value,\
62 const DJacobianRangeType &gradient,\
63 RRangeType &flux ) const \
64 { impl().linSource(uBar, gradientBar, x, value, gradient, flux); } \
65 virtual void diffusiveFlux ( const POINT &x,\
66 const DRangeType &value,\
67 const DJacobianRangeType &gradient,\
68 RJacobianRangeType &flux ) const \
69 { impl().diffusiveFlux(x, value, gradient, flux); } \
70 virtual void linDiffusiveFlux ( const DRangeType& uBar,\
71 const DJacobianRangeType& gradientBar,\
73 const DRangeType &value,\
74 const DJacobianRangeType &gradient,\
75 RJacobianRangeType &flux ) const \
76 { impl().linDiffusiveFlux(uBar, gradientBar, x, value, gradient, flux); } \
77 virtual void fluxDivergence( const POINT &x,\
78 const DRangeType &value,\
79 const DJacobianRangeType &jacobian,\
80 const DHessianRangeType &hessian,\
81 RRangeType &flux) const \
82 { impl().fluxDivergence(x, value, jacobian, hessian, flux); } \
83 virtual void alpha(const POINT &x,\
84 const DRangeType &value,\
85 RRangeType &val) const \
86 { impl().alpha(x, value, val); } \
87 virtual void linAlpha(const DRangeType &uBar,\
89 const DRangeType &value,\
90 RRangeType &val) const \
91 { impl().linAlpha(uBar, x, value, val); } \
92 virtual void dirichlet( int bndId, const POINT &x,\
93 RRangeType &value) const \
94 { impl().dirichlet(bndId, x, value); }
96template<
class Gr
idPart,
int dimDomain,
int dimRange=dimDomain,
class RangeField =
double >
100 static const int dimD = dimDomain;
101 static const int dimR = dimRange;
119 typedef typename GridPartType::template Codim<0>::EntityType
EntityType;
124 template <
class F,
int d>
125 using Traits = Dune::FemPy::FempyQuadratureTraits<F,d>;
155 virtual std::string
name()
const = 0;
161 virtual void setTime(
const double t )
const = 0;
165 virtual double time()
const = 0;
181 typedef std::array<int, dimRange> DirichletComponentType;
190 class CheckTimeMethod
193 static std::true_type testSignature(
double& (T::*)());
196 static decltype(testSignature(&T::time)) test(std::nullptr_t);
199 static std::false_type test(...);
201 using type =
decltype(test<M>(
nullptr));
203 static const bool value = type::value;
207 template <
class Model,
bool>
210 static void setTime( Model&,
const double ) {}
211 static double time(
const Model& ) {
return -1.0; }
214 template <
class Model>
215 struct CallSetTime< Model, true >
217 static void setTime( Model& model,
const double t ) { model.time() = t; }
218 static double time (
const Model& model ) {
return model.time(); }
223template <
class ModelImpl >
225 typename ModelImpl::RRangeFieldType>
227 typedef ModelImpl Impl;
229 static const int dimD = ModelImpl::dimD;
230 static const int dimR = ModelImpl::dimR;
252 template<
class... Args, std::enable_if_t< std::is_constructible< ModelImpl, Args &&... >::value,
int > = 0 >
254 : impl_(
std::forward< Args >( args )... )
276 virtual std::string
name()
const
278 return impl().name();
285 detail::CallSetTime< ModelImpl, detail::CheckTimeMethod< ModelImpl >::value >
293 return detail::CallSetTime< ModelImpl, detail::CheckTimeMethod< ModelImpl >::value >
::time(
impl() );
299 return impl().hasDirichletBoundary();
303 return impl().hasNeumanBoundary();
307 return impl().isDirichletIntersection(inter, dirichletComponent);
311 return impl().init(entity);
325#define VirtualPenaltyMethods(POINT) \
326 virtual RRangeType penalty ( const POINT &x,\
327 const DRangeType &value ) const = 0;\
328 virtual RRangeType linPenalty ( const POINT &x,\
329 const DRangeType &value ) const = 0;
330#define WrapperPenaltyMethods(POINT) \
331 virtual RRangeType penalty ( const POINT &x,\
332 const DRangeType &value ) const \
333 { return impl().penalty(x, value); } \
334 virtual RRangeType linPenalty ( const POINT &x,\
335 const DRangeType &value ) const \
336 { return impl().linPenalty( x, value); }
338template<
class Gr
idPart,
int dimDomain,
int dimRange=dimDomain,
class RangeField =
double >
342 static const int dimD = dimDomain;
343 static const int dimR = dimRange;
361 typedef typename GridPartType::template Codim<0>::EntityType
EntityType;
366 template <
class F,
int d>
367 using Traits = Dune::FemPy::FempyQuadratureTraits<F,d>;
397 virtual std::string
name()
const = 0;
403 virtual void setTime(
const double t )
const = 0;
407 virtual double time()
const = 0;
423 typedef std::array<int, dimRange> DirichletComponentType;
443template <
class ModelImpl >
445 typename ModelImpl::RRangeFieldType>
447 typedef ModelImpl Impl;
449 static const int dimD = ModelImpl::dimD;
450 static const int dimR = ModelImpl::dimR;
472 template<
class... Args, std::enable_if_t< std::is_constructible< ModelImpl, Args &&... >::value,
int > = 0 >
474 : impl_(
std::forward< Args >( args )... )
508 virtual std::string
name()
const
510 return impl().name();
517 detail::CallSetTime< ModelImpl, detail::CheckTimeMethod< ModelImpl >::value >
525 return detail::CallSetTime< ModelImpl, detail::CheckTimeMethod< ModelImpl >::value >
::time(
impl() );
531 return impl().hasDirichletBoundary();
535 return impl().hasNeumanBoundary();
539 return impl().isDirichletIntersection(inter, dirichletComponent);
543 return impl().init(entity);
558template <
class ModelTraits >
559struct DiffusionModelEngine :
public DiffusionModel<typename ModelTraits::GridPartType, ModelTraits::dimD, ModelTraits::dimR, typename ModelTraits::RRangeFieldType>
561 typedef typename ModelTraits::GridPartType
GridPartType;
562 static const int dimD = ModelTraits::dimD;
563 static const int dimR = ModelTraits::dimR;
566 typedef typename Base::Point
Point;
581 typedef typename ModelTraits::ModelType ModelImpl;
583 DiffusionModelEngine(ModelImpl &model) : impl_(model) {}
584 ~DiffusionModelEngine()
595 virtual std::string
name()
const
597 return impl().name();
601 return impl().hasDirichletBoundary();
605 return impl().hasNeumanBoundary();
607 virtual bool isDirichletIntersection(
const IntersectionType& inter, std::array<int, dimR> &dirichletComponent )
const
609 return impl().isDirichletIntersection(inter, dirichletComponent);
611 virtual bool init(
const EntityType &entity)
const
613 return impl().init(entity);
615 const ModelImpl &impl()
const
#define WrapperDiffusionModelMethods(POINT)
Definition: diffusionmodel.hh:52
#define VirtualDiffusionModelMethods(POINT)
Definition: diffusionmodel.hh:16
#define VirtualPenaltyMethods(POINT)
Definition: diffusionmodel.hh:325
Definition: explicitfieldvector.hh:75
quadrature class supporting base function caching
Definition: cachingquadrature.hh:41
quadrature on the codim-0 reference element
Definition: elementquadrature.hh:18
Definition: diffusionmodel.hh:98
DiffusionModel()
Definition: diffusionmodel.hh:151
Dune::Fem::CachingQuadrature< GridPartType, 0, Traits >::QuadraturePointWrapperType Point
Definition: diffusionmodel.hh:127
Dune::Fem::CachingQuadrature< GridPartType, 1 >::QuadraturePointWrapperType OriginalIntersectionPoint
Definition: diffusionmodel.hh:139
Dune::Fem::ElementQuadrature< GridPartType, 1 >::QuadraturePointWrapperType OriginalElementIntersectionPoint
Definition: diffusionmodel.hh:143
DiffusionModel< GridPartType, dimD, dimR, RangeField > ModelType
Definition: diffusionmodel.hh:102
Dune::Fem::CachingQuadrature< GridPartType, 1, Traits >::QuadraturePointWrapperType IntersectionPoint
Definition: diffusionmodel.hh:129
DFunctionSpaceType::DomainFieldType DDomainFieldType
Definition: diffusionmodel.hh:113
RFunctionSpaceType::DomainFieldType rDomainFieldType
Definition: diffusionmodel.hh:117
Dune::Fem::ElementQuadrature< GridPartType, 0, Traits >::QuadraturePointWrapperType ElementPoint
Definition: diffusionmodel.hh:131
virtual std::string name() const =0
virtual bool hasNeumanBoundary() const =0
virtual bool init(const EntityType &entity) const =0
Dune::Fem::ElementQuadrature< GridPartType, 1, Traits >::QuadraturePointWrapperType ElementIntersectionPoint
Definition: diffusionmodel.hh:133
RFunctionSpaceType::JacobianRangeType RJacobianRangeType
Definition: diffusionmodel.hh:115
Dune::Fem::ElementQuadrature< GridPartType, 0 >::QuadraturePointWrapperType OriginalElementPoint
Definition: diffusionmodel.hh:141
virtual double time() const =0
virtual ~DiffusionModel()
Definition: diffusionmodel.hh:153
Dune::Fem::FunctionSpace< double, RangeFieldType, GridPart::dimensionworld, dimR > RFunctionSpaceType
Definition: diffusionmodel.hh:108
GridPartType::IntersectionType IntersectionType
Definition: diffusionmodel.hh:120
DFunctionSpaceType::JacobianRangeType DJacobianRangeType
Definition: diffusionmodel.hh:111
virtual bool isDirichletIntersection(const IntersectionType &inter, DirichletComponentType &dirichletComponent) const =0
GridPartType::template Codim< 0 >::EntityType EntityType
Definition: diffusionmodel.hh:119
RFunctionSpaceType::HessianRangeType RHessianRangeType
Definition: diffusionmodel.hh:116
static const int dimD
Definition: diffusionmodel.hh:100
Dune::Fem::CachingQuadrature< GridPartType, 0 >::QuadraturePointWrapperType OriginalPoint
Definition: diffusionmodel.hh:137
RangeField RangeFieldType
Definition: diffusionmodel.hh:103
EntityType::Geometry::LocalCoordinate LocalDomainType
Definition: diffusionmodel.hh:121
Dune::FemPy::FempyQuadratureTraits< F, d > Traits
Definition: diffusionmodel.hh:125
Dune::Fem::FunctionSpace< double, RangeFieldType, GridPart::dimensionworld, dimD > DFunctionSpaceType
Definition: diffusionmodel.hh:106
virtual VirtualDiffusionModelMethods(Point) VirtualDiffusionModelMethods(ElementPoint) VirtualDiffusionModelMethods(IntersectionPoint) VirtualDiffusionModelMethods(ElementIntersectionPoint) VirtualDiffusionModelMethods(OriginalPoint) VirtualDiffusionModelMethods(OriginalElementPoint) VirtualDiffusionModelMethods(OriginalIntersectionPoint) VirtualDiffusionModelMethods(OriginalElementIntersectionPoint) VirtualDiffusionModelMethods(LocalDomainType) typedef std bool hasDirichletBoundary() const =0
GridPart GridPartType
Definition: diffusionmodel.hh:99
DFunctionSpaceType::HessianRangeType DHessianRangeType
Definition: diffusionmodel.hh:112
DFunctionSpaceType::DomainType DomainType
Definition: diffusionmodel.hh:109
static const int dimR
Definition: diffusionmodel.hh:101
RFunctionSpaceType::RangeType RRangeType
Definition: diffusionmodel.hh:114
virtual void setTime(const double t) const =0
DFunctionSpaceType::RangeType DRangeType
Definition: diffusionmodel.hh:110
Definition: diffusionmodel.hh:226
Base::LocalDomainType LocalDomainType
Definition: diffusionmodel.hh:241
Base::DRangeType DRangeType
Definition: diffusionmodel.hh:243
DiffusionModelWrapper(Args &&... args)
Definition: diffusionmodel.hh:253
virtual bool hasNeumanBoundary() const
Definition: diffusionmodel.hh:301
WrapperDiffusionModelMethods(OriginalPoint)
virtual double time() const
Definition: diffusionmodel.hh:291
static const int dimR
Definition: diffusionmodel.hh:230
Base::RRangeType RRangeType
Definition: diffusionmodel.hh:246
virtual std::string name() const
Definition: diffusionmodel.hh:276
Base::DJacobianRangeType DJacobianRangeType
Definition: diffusionmodel.hh:244
~DiffusionModelWrapper()
Definition: diffusionmodel.hh:257
std::array< int, dimR > DirichletComponentType
Definition: diffusionmodel.hh:296
WrapperDiffusionModelMethods(OriginalIntersectionPoint)
Base::ElementPoint ElementPoint
Definition: diffusionmodel.hh:235
Base::RJacobianRangeType RJacobianRangeType
Definition: diffusionmodel.hh:247
Base::ElementIntersectionPoint ElementIntersectionPoint
Definition: diffusionmodel.hh:236
WrapperDiffusionModelMethods(IntersectionPoint)
Base::OriginalElementIntersectionPoint OriginalElementIntersectionPoint
Definition: diffusionmodel.hh:240
virtual bool hasDirichletBoundary() const
Definition: diffusionmodel.hh:297
const ModelImpl & impl() const
Definition: diffusionmodel.hh:313
WrapperDiffusionModelMethods(OriginalElementIntersectionPoint)
virtual bool init(const EntityType &entity) const
Definition: diffusionmodel.hh:309
Base::DomainType DomainType
Definition: diffusionmodel.hh:242
Base::EntityType EntityType
Definition: diffusionmodel.hh:249
virtual void setTime(const double t) const
Definition: diffusionmodel.hh:283
Base::OriginalElementPoint OriginalElementPoint
Definition: diffusionmodel.hh:239
Base::IntersectionType IntersectionType
Definition: diffusionmodel.hh:250
Base::Point Point
Definition: diffusionmodel.hh:233
WrapperDiffusionModelMethods(ElementPoint)
Base::IntersectionPoint IntersectionPoint
Definition: diffusionmodel.hh:234
WrapperDiffusionModelMethods(OriginalElementPoint)
WrapperDiffusionModelMethods(ElementIntersectionPoint)
ModelImpl & impl()
Definition: diffusionmodel.hh:317
Base::OriginalIntersectionPoint OriginalIntersectionPoint
Definition: diffusionmodel.hh:238
DiffusionModel< GridPartType, dimD, dimR, typename ModelImpl::RRangeFieldType > Base
Definition: diffusionmodel.hh:231
static const int dimD
Definition: diffusionmodel.hh:229
Base::OriginalPoint OriginalPoint
Definition: diffusionmodel.hh:237
WrapperDiffusionModelMethods(LocalDomainType)
virtual bool isDirichletIntersection(const IntersectionType &inter, DirichletComponentType &dirichletComponent) const
Definition: diffusionmodel.hh:305
Base::RHessianRangeType RHessianRangeType
Definition: diffusionmodel.hh:248
ModelImpl::GridPartType GridPartType
Definition: diffusionmodel.hh:228
Base::DHessianRangeType DHessianRangeType
Definition: diffusionmodel.hh:245
WrapperDiffusionModelMethods(Point)
Definition: diffusionmodel.hh:340
GridPartType::IntersectionType IntersectionType
Definition: diffusionmodel.hh:362
DFunctionSpaceType::RangeType DRangeType
Definition: diffusionmodel.hh:352
Dune::Fem::CachingQuadrature< GridPartType, 0, Traits >::QuadraturePointWrapperType Point
Definition: diffusionmodel.hh:369
Dune::Fem::ElementQuadrature< GridPartType, 1, Traits >::QuadraturePointWrapperType ElementIntersectionPoint
Definition: diffusionmodel.hh:375
DFunctionSpaceType::DomainFieldType DDomainFieldType
Definition: diffusionmodel.hh:355
virtual double time() const =0
Dune::Fem::FunctionSpace< double, RangeFieldType, GridPart::dimensionworld, dimD > DFunctionSpaceType
Definition: diffusionmodel.hh:348
GridPart GridPartType
Definition: diffusionmodel.hh:341
GridPartType::template Codim< 0 >::EntityType EntityType
Definition: diffusionmodel.hh:361
static const int dimR
Definition: diffusionmodel.hh:343
Dune::Fem::CachingQuadrature< GridPartType, 1, Traits >::QuadraturePointWrapperType IntersectionPoint
Definition: diffusionmodel.hh:371
virtual std::string name() const =0
RFunctionSpaceType::HessianRangeType RHessianRangeType
Definition: diffusionmodel.hh:358
RangeField RangeFieldType
Definition: diffusionmodel.hh:345
EntityType::Geometry::LocalCoordinate LocalDomainType
Definition: diffusionmodel.hh:363
virtual bool init(const EntityType &entity) const =0
DFunctionSpaceType::DomainType DomainType
Definition: diffusionmodel.hh:351
Dune::Fem::CachingQuadrature< GridPartType, 0 >::QuadraturePointWrapperType OriginalPoint
Definition: diffusionmodel.hh:379
Dune::Fem::ElementQuadrature< GridPartType, 0 >::QuadraturePointWrapperType OriginalElementPoint
Definition: diffusionmodel.hh:383
Dune::Fem::CachingQuadrature< GridPartType, 1 >::QuadraturePointWrapperType OriginalIntersectionPoint
Definition: diffusionmodel.hh:381
virtual ~DGDiffusionModel()
Definition: diffusionmodel.hh:395
Dune::FemPy::FempyQuadratureTraits< F, d > Traits
Definition: diffusionmodel.hh:367
Dune::Fem::ElementQuadrature< GridPartType, 0, Traits >::QuadraturePointWrapperType ElementPoint
Definition: diffusionmodel.hh:373
DFunctionSpaceType::JacobianRangeType DJacobianRangeType
Definition: diffusionmodel.hh:353
Dune::Fem::FunctionSpace< double, RangeFieldType, GridPart::dimensionworld, dimR > RFunctionSpaceType
Definition: diffusionmodel.hh:350
virtual bool hasNeumanBoundary() const =0
DGDiffusionModel< GridPartType, dimD, dimR, RangeField > ModelType
Definition: diffusionmodel.hh:344
RFunctionSpaceType::DomainFieldType rDomainFieldType
Definition: diffusionmodel.hh:359
static const int dimD
Definition: diffusionmodel.hh:342
virtual VirtualDiffusionModelMethods(Point) VirtualDiffusionModelMethods(ElementPoint) VirtualDiffusionModelMethods(IntersectionPoint) VirtualDiffusionModelMethods(ElementIntersectionPoint) VirtualDiffusionModelMethods(OriginalPoint) VirtualDiffusionModelMethods(OriginalElementPoint) VirtualDiffusionModelMethods(OriginalIntersectionPoint) VirtualDiffusionModelMethods(OriginalElementIntersectionPoint) VirtualDiffusionModelMethods(LocalDomainType) typedef std bool hasDirichletBoundary() const =0
Dune::Fem::ElementQuadrature< GridPartType, 1 >::QuadraturePointWrapperType OriginalElementIntersectionPoint
Definition: diffusionmodel.hh:385
DFunctionSpaceType::HessianRangeType DHessianRangeType
Definition: diffusionmodel.hh:354
virtual void setTime(const double t) const =0
RFunctionSpaceType::JacobianRangeType RJacobianRangeType
Definition: diffusionmodel.hh:357
virtual bool isDirichletIntersection(const IntersectionType &inter, DirichletComponentType &dirichletComponent) const =0
DGDiffusionModel()
Definition: diffusionmodel.hh:393
RFunctionSpaceType::RangeType RRangeType
Definition: diffusionmodel.hh:356
Definition: diffusionmodel.hh:446
static const int dimD
Definition: diffusionmodel.hh:449
static const int dimR
Definition: diffusionmodel.hh:450
virtual bool hasNeumanBoundary() const
Definition: diffusionmodel.hh:533
WrapperDiffusionModelMethods(OriginalElementIntersectionPoint)
WrapperDiffusionModelMethods(OriginalIntersectionPoint)
Base::LocalDomainType LocalDomainType
Definition: diffusionmodel.hh:461
Base::DRangeType DRangeType
Definition: diffusionmodel.hh:463
Base::EntityType EntityType
Definition: diffusionmodel.hh:469
Base::Point Point
Definition: diffusionmodel.hh:453
Base::ElementIntersectionPoint ElementIntersectionPoint
Definition: diffusionmodel.hh:456
WrapperPenaltyMethods(LocalDomainType)
Base::RRangeType RRangeType
Definition: diffusionmodel.hh:466
WrapperPenaltyMethods(OriginalIntersectionPoint)
WrapperDiffusionModelMethods(LocalDomainType)
Base::OriginalPoint OriginalPoint
Definition: diffusionmodel.hh:457
WrapperDiffusionModelMethods(IntersectionPoint)
DGDiffusionModelWrapper(Args &&... args)
Definition: diffusionmodel.hh:473
Base::OriginalElementIntersectionPoint OriginalElementIntersectionPoint
Definition: diffusionmodel.hh:460
WrapperDiffusionModelMethods(Point)
Base::OriginalElementPoint OriginalElementPoint
Definition: diffusionmodel.hh:459
Base::DomainType DomainType
Definition: diffusionmodel.hh:462
Base::IntersectionPoint IntersectionPoint
Definition: diffusionmodel.hh:454
~DGDiffusionModelWrapper()
Definition: diffusionmodel.hh:477
Base::IntersectionType IntersectionType
Definition: diffusionmodel.hh:470
std::array< int, dimR > DirichletComponentType
Definition: diffusionmodel.hh:528
WrapperDiffusionModelMethods(OriginalElementPoint)
Base::OriginalIntersectionPoint OriginalIntersectionPoint
Definition: diffusionmodel.hh:458
virtual bool hasDirichletBoundary() const
Definition: diffusionmodel.hh:529
virtual void setTime(const double t) const
Definition: diffusionmodel.hh:515
Base::RJacobianRangeType RJacobianRangeType
Definition: diffusionmodel.hh:467
DGDiffusionModel< GridPartType, dimD, dimR, typename ModelImpl::RRangeFieldType > Base
Definition: diffusionmodel.hh:451
WrapperPenaltyMethods(OriginalElementPoint)
ModelImpl::GridPartType GridPartType
Definition: diffusionmodel.hh:448
Base::DHessianRangeType DHessianRangeType
Definition: diffusionmodel.hh:465
virtual double time() const
Definition: diffusionmodel.hh:523
WrapperDiffusionModelMethods(ElementIntersectionPoint)
WrapperDiffusionModelMethods(OriginalPoint)
Base::ElementPoint ElementPoint
Definition: diffusionmodel.hh:455
virtual bool init(const EntityType &entity) const
Definition: diffusionmodel.hh:541
WrapperPenaltyMethods(OriginalElementIntersectionPoint)
const ModelImpl & impl() const
Definition: diffusionmodel.hh:545
Base::RHessianRangeType RHessianRangeType
Definition: diffusionmodel.hh:468
WrapperDiffusionModelMethods(ElementPoint)
virtual std::string name() const
Definition: diffusionmodel.hh:508
Base::DJacobianRangeType DJacobianRangeType
Definition: diffusionmodel.hh:464
ModelImpl & impl()
Definition: diffusionmodel.hh:549
WrapperPenaltyMethods(Point) WrapperPenaltyMethods(ElementPoint) WrapperPenaltyMethods(IntersectionPoint) WrapperPenaltyMethods(ElementIntersectionPoint) WrapperPenaltyMethods(OriginalPoint)
virtual bool isDirichletIntersection(const IntersectionType &inter, DirichletComponentType &dirichletComponent) const
Definition: diffusionmodel.hh:537
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