1 #ifndef DUNE_PDELAB_LOCALOPERATOR_STOKESPARAMETER_HH 2 #define DUNE_PDELAB_LOCALOPERATOR_STOKESPARAMETER_HH 4 #include <dune/common/parametertree.hh> 44 template<
typename GV,
typename RF>
53 dimDomain = GV::dimension
60 typedef Dune::FieldVector<DomainField,dimDomain>
Domain;
78 typedef typename GV::Traits::template Codim<0>::Entity
Element;
89 template<
typename GF,
typename Entity,
typename Domain>
90 typename GF::Traits::RangeType
91 evaluateVelocityGridFunction(
const GF& gf,
95 static_assert(
int(GF::Traits::dimRange) ==
int(Domain::dimension),
"dimension of function range does not match grid dimension");
96 typename GF::Traits::RangeType y;
105 template<
typename GF,
typename Entity,
typename Domain>
106 FieldVector<typename GF::template Child<0>::Type::Traits::RangeFieldType,GF::CHILDREN>
107 evaluateVelocityGridFunction(
const GF& gf,
111 static_assert(Domain::dimension == GF::CHILDREN,
"dimension of function range does not match grid dimension");
112 FieldVector<typename GF::template Child<0>::Type::Traits::RangeFieldType,GF::CHILDREN> y;
113 typename GF::template Child<0>::Type::Traits::RangeType cy;
114 for (
int d = 0; d < Domain::dimension; ++d)
116 gf.child(d).evaluate(e,x,cy);
142 template <
typename GV,
typename RF,
typename F,
typename B,
typename V,
typename J,
bool navier = false,
bool tensor = false>
147 static const bool assemble_navier = navier;
148 static const bool assemble_full_tensor = tensor;
159 : _rho(config.get<double>(
"rho"))
160 , _mu(config.get<double>(
"mu"))
183 template<
typename EG>
187 typename F::Traits::RangeType fvalue;
188 return evaluateVelocityGridFunction(_f,e.entity(),x);
192 template<
typename IG>
193 typename Traits::BoundaryCondition::Type
197 typename B::Traits::RangeType y;
203 template<
typename EG>
211 template<
typename IG>
219 template<
typename EG>
227 template<
typename IG>
235 template<
typename EG>
239 typename V::Traits::RangeType y;
245 template<
typename IG>
247 DUNE_DEPRECATED_MSG(
"Deprecated in DUNE-PDELab 2.4, use entity-based version instead!")
248 g(const IG&
ig, const typename Traits::IntersectionDomain& x)
const 250 auto e = ig.inside();
251 typename V::Traits::RangeType y;
252 _v.evaluate(
e,ig.geometryInInside().global(x),y);
257 template<
typename EG>
267 template<
typename IG>
276 template<
typename IG>
278 J::Traits::dimRange == 1 &&
279 (GV::dimension > 1) &&
287 typename J::Traits::RangeType r;
288 auto e = ig.inside();
289 _j.evaluate(
e,ig.geometryInInside().global(x),r);
295 template<
typename IG>
297 J::Traits::dimRange == GV::dimension &&
305 auto e = ig.inside();
306 typename J::Traits::RangeType y;
307 _j.evaluate(
e,ig.geometryInInside().global(x),y);
334 template<
typename PRM>
350 const Dune::FieldVector<typename I::ctype, I::dimension-1> & coord)
const 361 template<
typename PRM>
377 const Dune::FieldVector<typename I::ctype, I::dimension-1> & coord)
const 388 template<
typename PRM,
int rangeDim>
389 class NavierStokesFunctionAdapterBase
391 Dune::PDELab::GridFunctionTraits<
392 typename PRM::Traits::GridView,
393 typename PRM::Traits::RangeField,
395 Dune::FieldVector<typename PRM::Traits::RangeField,rangeDim>
397 NavierStokesFunctionAdapterBase<PRM,rangeDim>
403 typename PRM::Traits::GridView,
404 typename PRM::Traits::RangeField,
406 Dune::FieldVector<typename PRM::Traits::RangeField,rangeDim>
410 NavierStokesFunctionAdapterBase(PRM& prm)
414 void setTime(
const double time)
419 const PRM& parameters()
const 425 const typename Traits::GridViewType& getGridView ()
const 427 return _prm.gridView();
443 template<
typename PRM>
445 :
public NavierStokesFunctionAdapterBase<PRM,PRM::Traits::dimDomain>
449 typedef NavierStokesFunctionAdapterBase<PRM,PRM::Traits::dimDomain> Base;
451 using Base::parameters;
464 const typename Traits::DomainType& x,
465 typename Traits::RangeType& y)
const 467 y = parameters().g(e,x);
476 template <
typename PRM >
477 class NavierStokesDirichletFunctionAdapterFactory
481 NavierStokesDirichletFunctionAdapter<PRM>,
482 NavierStokesPressureDirichletFunctionAdapter<PRM> >
483 BoundaryDirichletFunction;
485 NavierStokesDirichletFunctionAdapterFactory(PRM & prm)
486 : v(prm),
p(prm), df(v,
p)
489 BoundaryDirichletFunction & dirichletFunction()
495 NavierStokesVelocityDirichletFunctionAdapter<PRM> v;
496 NavierStokesPressureDirichletFunctionAdapter<PRM>
p;
497 BoundaryDirichletFunction df;
Traits::BoundaryCondition::Type bctype(const IG &is, const typename Traits::IntersectionDomain &x) const
boundary condition type from local intersection coordinate
Definition: stokesparameter.hh:194
void evaluate(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeType &y) const
Evaluate dirichlet function.
Definition: stokesparameter.hh:463
Dune::FieldVector< DomainField, dimDomain-1 > IntersectionDomain
domain type
Definition: stokesparameter.hh:63
static const unsigned int value
Definition: gridfunctionspace/tags.hh:177
composite functions
Definition: function.hh:799
Definition: stokesparameter.hh:444
StokesBoundaryCondition BoundaryCondition
boundary type value
Definition: stokesparameter.hh:75
GV::Intersection Intersection
grid intersection type
Definition: stokesparameter.hh:80
RF RangeField
Export type for range field.
Definition: stokesparameter.hh:66
const IG & ig
Definition: constraints.hh:148
StokesVelocityDirichletConstraints(const PRM &_prm)
Constructor.
Definition: stokesparameter.hh:344
GV GridView
the grid view
Definition: stokesparameter.hh:48
bool isDirichlet(const I &intersection, const Dune::FieldVector< typename I::ctype, I::dimension-1 > &coord) const
Definition: stokesparameter.hh:376
Definition: adaptivity.hh:27
Definition: stokesparameter.hh:45
Definition: stokesparameter.hh:32
NavierStokesDefaultParameters(const RF &mu, const RF &rho, F &f, B &b, V &v, J &j)
Definition: stokesparameter.hh:167
bool isDirichlet(const I &intersection, const Dune::FieldVector< typename I::ctype, I::dimension-1 > &coord) const
Definition: stokesparameter.hh:349
const Entity & e
Definition: localfunctionspace.hh:111
Base::Traits Traits
Definition: stokesparameter.hh:455
Traits::VelocityRange f(const EG &e, const typename Traits::Domain &x) const
source term
Definition: stokesparameter.hh:185
GV::Grid::ctype DomainField
Export type for domain field.
Definition: stokesparameter.hh:57
Definition: stokesparameter.hh:36
Definition: stokesparameter.hh:335
StokesPressureDirichletConstraints(const PRM &_prm)
Constructor.
Definition: stokesparameter.hh:371
Definition: stokesparameter.hh:37
NavierStokesParameterTraits< GV, RF > Traits
Type traits.
Definition: stokesparameter.hh:151
Traits::RangeField rho(const IG &ig, const typename Traits::IntersectionDomain &x) const
Density value from local intersection coordinate.
Definition: stokesparameter.hh:229
Definition: stokesparameter.hh:35
Definition: stokesparameter.hh:34
Definition: stokesparameter.hh:143
Definition: constraintsparameters.hh:24
Traits::RangeField rho(const EG &eg, const typename Traits::Domain &x) const
Density value from local cell coordinate.
Definition: stokesparameter.hh:221
Dune::FieldVector< RF, 1 > PressureRange
pressure range type
Definition: stokesparameter.hh:72
const P & p
Definition: constraints.hh:147
Traits::VelocityRange g(const EG &e, const typename Traits::Domain &x) const
Dirichlet boundary condition value from local cell coordinate.
Definition: stokesparameter.hh:237
Type
Definition: stokesparameter.hh:33
GV::Traits::template Codim< 0 >::Entity Element
grid element type
Definition: stokesparameter.hh:78
Traits::RangeField mu(const IG &ig, const typename Traits::IntersectionDomain &x) const
Dynamic viscosity value from local intersection coordinate.
Definition: stokesparameter.hh:213
void setTime(RF time)
Definition: stokesparameter.hh:313
NavierStokesDefaultParameters(const Dune::ParameterTree &config, F &f, B &b, V &v, J &j)
Constructor.
Definition: stokesparameter.hh:154
Traits::RangeField mu(const EG &e, const typename Traits::Domain &x) const
Dynamic viscosity value from local cell coordinate.
Definition: stokesparameter.hh:205
Dune::FieldVector< DomainField, dimDomain > Domain
domain type
Definition: stokesparameter.hh:60
Dune::FieldVector< RF, GV::dimensionworld > VelocityRange
deformation range type
Definition: stokesparameter.hh:69
Definition: stokesparameter.hh:362
leaf of a function tree
Definition: function.hh:576
Traits::RangeField g2(const EG &e, const typename Traits::Domain &x) const
pressure source term
Definition: stokesparameter.hh:259
NavierStokesVelocityFunctionAdapter(PRM &prm)
Constructor.
Definition: stokesparameter.hh:458
traits class holding the function signature, same as in local function
Definition: function.hh:175