DUNE PDELab (2.8)

convectiondiffusionparameter.hh
1// -*- tab-width: 4; indent-tabs-mode: nil -*-
2#ifndef DUNE_PDELAB_LOCALOPERATOR_CONVECTIONDIFFUSIONPARAMETER_HH
3#define DUNE_PDELAB_LOCALOPERATOR_CONVECTIONDIFFUSIONPARAMETER_HH
4
5#include<vector>
6#include<type_traits>
7
12#include<dune/geometry/referenceelements.hh>
14#include<dune/pdelab/common/geometrywrapper.hh>
15#include<dune/pdelab/common/function.hh>
16#include<dune/pdelab/common/functionutilities.hh>
17#include<dune/pdelab/constraints/common/constraintsparameters.hh>
18
19namespace Dune {
20 namespace PDELab {
21
22 #ifndef DOXYGEN
23
24 // helper construct for backwards-compatible addition of hasPermeabilityIsConstantPerCell()
25
26 namespace Impl {
27
28 template<typename T, typename = void>
29 struct hasPermeabilityIsConstantPerCell
30 : public std::false_type
31 {};
32
33 template<typename T>
34 struct hasPermeabilityIsConstantPerCell<
35 T,
36 void_t<decltype(std::declval<T>().permeabilityIsConstantPerCell())>
37 >
38 : public std::true_type
39 {};
40
41 template<typename T>
42 [[deprecated("Starting from PDELab 2.6, parameter classes must have a method `bool permeabilityIsConstantPerCell()`. For now, we assume a default value of true.")]]
43 constexpr
44 std::enable_if_t<
45 not hasPermeabilityIsConstantPerCell<T>::value,
46 bool
47 >
48 permeabilityIsConstantPerCell(const T& param)
49 {
50 return true;
51 }
52
53 template<typename T>
54 constexpr
55 std::enable_if_t<
56 hasPermeabilityIsConstantPerCell<T>::value,
57 bool
58 >
59 permeabilityIsConstantPerCell(const T& param)
60 {
61 return param.permeabilityIsConstantPerCell();
62 }
63
64 } // namespace Impl
65
66 #endif // DOXYGEN
67
74 template<typename GV, typename RF>
76 {
78 typedef GV GridViewType;
79
81 enum {
83 dimDomain = GV::dimension
84 };
85
87 typedef typename GV::Grid::ctype DomainFieldType;
88
91
94
96 typedef RF RangeFieldType;
97
100
103
105 typedef typename GV::Traits::template Codim<0>::Entity ElementType;
106 typedef typename GV::Intersection IntersectionType;
107 };
108
112 {
113 enum Type { Dirichlet=1, Neumann=-1, Outflow=-2, None=-3 }; // BC requiring constraints must be >0 if
114 // constraints assembler coming with PDELab is used
115 };
116
134 template<typename GV, typename RF>
136 {
137 typedef ConvectionDiffusionBoundaryConditions::Type BCType;
138
139 public:
141
143 static constexpr bool permeabilityIsConstantPerCell()
144 {
145 return true;
146 }
147
150 A (const typename Traits::ElementType& e, const typename Traits::DomainType& x) const
151 {
152 typename Traits::PermTensorType I;
153 for (std::size_t i=0; i<Traits::dimDomain; i++)
154 for (std::size_t j=0; j<Traits::dimDomain; j++)
155 I[i][j] = (i==j) ? 1 : 0;
156 return I;
157 }
158
160 typename Traits::RangeType
161 b (const typename Traits::ElementType& e, const typename Traits::DomainType& x) const
162 {
163 typename Traits::RangeType v(0.0);
164 return v;
165 }
166
169 c (const typename Traits::ElementType& e, const typename Traits::DomainType& x) const
170 {
171 return 0.0;
172 }
173
176 f (const typename Traits::ElementType& e, const typename Traits::DomainType& x) const
177 {
178 return 0.0;
179 }
180
182 BCType
183 bctype (const typename Traits::IntersectionType& is, const typename Traits::IntersectionDomainType& x) const
184 {
185 return ConvectionDiffusionBoundaryConditions::Dirichlet;
186 }
187
190 g (const typename Traits::ElementType& e, const typename Traits::DomainType& x) const
191 {
192 typename Traits::DomainType xglobal = e.geometry().global(x);
193 return xglobal.two_norm();
194 }
195
198 j (const typename Traits::IntersectionType& is, const typename Traits::IntersectionDomainType& x) const
199 {
200 return 0.0;
201 }
202
205 o (const typename Traits::IntersectionType& is, const typename Traits::IntersectionDomainType& x) const
206 {
207 return 0.0;
208 }
209 };
210
211
216 template<typename T>
218 :
220 public Dune::PDELab::DirichletConstraintsParameters /*@\label{bcp:base}@*/
221 {
222 const T& t;
223
224 public:
225
226 ConvectionDiffusionBoundaryConditionAdapter(const typename T::Traits::GridViewType& gv_,
227 const T& t_ )
228 : t( t_ )
229 {}
230
232 : t( t_ )
233 {}
234
235 template<typename I>
236 bool isDirichlet(const I & ig /*@\label{bcp:name}@*/
238 ) const
239 {
240 return( t.bctype( ig.intersection(), coord )
241 == ConvectionDiffusionBoundaryConditions::Dirichlet );
242 }
243
244 template<typename I>
245 bool isNeumann(const I & ig, /*@\label{bcp:name}@*/
247 ) const
248 {
249 return !isDirichlet( ig, coord );
250 }
251
252 };
253
254
255
256
257
262 template<typename T>
264 : public Dune::PDELab::AnalyticGridFunctionBase<Dune::PDELab::AnalyticGridFunctionTraits
265 <typename T::Traits::GridViewType,
266 typename T::Traits::RangeFieldType,
267 T::Traits::GridViewType::dimension>,
268 ConvectionDiffusionVelocityExtensionAdapter<T> >
269 {
270 public:
271 typedef Dune::PDELab::AnalyticGridFunctionTraits<typename T::Traits::GridViewType,
272 typename T::Traits::RangeFieldType,
273 T::Traits::GridViewType::dimension> Traits;
275
276
279 : BaseT(gv_), gv(gv_), t(t_)
280 {}
281
282 inline void evaluateGlobal (const typename Traits::DomainType& x,
283 typename Traits::RangeType& y) const
284 {
285 y = t.b(x);
286 }
287
289 inline void evaluate (const typename Traits::ElementType& e,
290 const typename Traits::DomainType& x,
291 typename Traits::RangeType& y) const
292 {
293 y = t.b(e,x);
294 }
295
296 inline const typename Traits::GridViewType& getGridView () const
297 {
298 return gv;
299 }
300
301 inline void setTime(double time_)
302 {
303 t.setTime(time_);
304 }
305
306 private:
307 const typename Traits::GridViewType gv;
308 T& t;
309 };
310
311
312
313
314
319 template<typename T>
321 : public Dune::PDELab::GridFunctionBase<Dune::PDELab::GridFunctionTraits<typename T::Traits::GridViewType,
322 typename T::Traits::RangeFieldType,
323 1,Dune::FieldVector<typename T::Traits::RangeFieldType,1> >
324 ,ConvectionDiffusionDirichletExtensionAdapter<T> >
325 {
326 public:
327 typedef Dune::PDELab::GridFunctionTraits<typename T::Traits::GridViewType,
328 typename T::Traits::RangeFieldType,
330
333 : g(g_), t(t_)
334 {}
335
337 inline void evaluate (const typename Traits::ElementType& e,
338 const typename Traits::DomainType& x,
339 typename Traits::RangeType& y) const
340 {
341 y = t.g(e,x);
342 }
343
344 inline const typename Traits::GridViewType& getGridView () const
345 {
346 return g;
347 }
348
349 inline void setTime(double time_)
350 {
351 t.setTime(time_);
352 }
353
354 private:
355 const typename Traits::GridViewType g;
356 T& t;
357 };
358
359
364template<typename T>
366 : public Dune::PDELab::GridFunctionBase<Dune::PDELab::GridFunctionTraits<typename T::Traits::GridViewType,
367 typename T::Traits::RangeFieldType,
368 T::Traits::GridViewType::dimension,Dune::FieldVector<typename T::Traits::RangeFieldType,T::Traits::GridViewType::dimension> >
369 ,ConvectionDiffusionExactGradientAdapter<T> >
370{
371public:
372 typedef Dune::PDELab::GridFunctionTraits<typename T::Traits::GridViewType,
373 typename T::Traits::RangeFieldType,
375
377 ConvectionDiffusionExactGradientAdapter (const typename Traits::GridViewType& g_, const T& t_) : g(g_), t(t_) {}
378
380 inline void evaluate (const typename Traits::ElementType& e,
381 const typename Traits::DomainType& x,
382 typename Traits::RangeType& y) const
383 {
384 y = t.gradient(e,x);
385 }
386
387 inline const typename Traits::GridViewType& getGridView () const
388 {
389 return g;
390 }
391
392private:
393 const typename Traits::GridViewType g;
394 const T& t;
395};
396 }
397}
398
399
400#endif // DUNE_PDELAB_LOCALOPERATOR_CONVECTIONDIFFUSIONPARAMETER_HH
FieldTraits< value_type >::real_type two_norm() const
two norm sqrt(sum over squared values of entries)
Definition: densevector.hh:642
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:95
an analytic grid function
Definition: function.hh:660
Definition: convectiondiffusionparameter.hh:221
Definition: convectiondiffusionparameter.hh:325
ConvectionDiffusionDirichletExtensionAdapter(const typename Traits::GridViewType &g_, T &t_)
constructor
Definition: convectiondiffusionparameter.hh:332
void evaluate(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeType &y) const
Evaluate the GridFunction at given position.
Definition: convectiondiffusionparameter.hh:337
Definition: convectiondiffusionparameter.hh:370
void evaluate(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeType &y) const
Evaluate the GridFunction at given position.
Definition: convectiondiffusionparameter.hh:380
ConvectionDiffusionExactGradientAdapter(const typename Traits::GridViewType &g_, const T &t_)
constructor
Definition: convectiondiffusionparameter.hh:377
Parameter class for solving the linear convection-diffusion equation.
Definition: convectiondiffusionparameter.hh:136
Traits::PermTensorType A(const typename Traits::ElementType &e, const typename Traits::DomainType &x) const
tensor diffusion coefficient
Definition: convectiondiffusionparameter.hh:150
Traits::RangeFieldType o(const typename Traits::IntersectionType &is, const typename Traits::IntersectionDomainType &x) const
outflow boundary condition
Definition: convectiondiffusionparameter.hh:205
Traits::RangeFieldType f(const typename Traits::ElementType &e, const typename Traits::DomainType &x) const
source term
Definition: convectiondiffusionparameter.hh:176
BCType bctype(const typename Traits::IntersectionType &is, const typename Traits::IntersectionDomainType &x) const
boundary condition type function
Definition: convectiondiffusionparameter.hh:183
Traits::RangeType b(const typename Traits::ElementType &e, const typename Traits::DomainType &x) const
velocity field
Definition: convectiondiffusionparameter.hh:161
static constexpr bool permeabilityIsConstantPerCell()
tensor diffusion constant per cell? return false if you want more than one evaluation of A per cell.
Definition: convectiondiffusionparameter.hh:143
Traits::RangeFieldType j(const typename Traits::IntersectionType &is, const typename Traits::IntersectionDomainType &x) const
Neumann boundary condition.
Definition: convectiondiffusionparameter.hh:198
Traits::RangeFieldType g(const typename Traits::ElementType &e, const typename Traits::DomainType &x) const
Dirichlet boundary condition value.
Definition: convectiondiffusionparameter.hh:190
Traits::RangeFieldType c(const typename Traits::ElementType &e, const typename Traits::DomainType &x) const
sink term
Definition: convectiondiffusionparameter.hh:169
Definition: convectiondiffusionparameter.hh:269
ConvectionDiffusionVelocityExtensionAdapter(const typename Traits::GridViewType &gv_, T &t_)
constructor
Definition: convectiondiffusionparameter.hh:278
void evaluate(const typename Traits::ElementType &e, const typename Traits::DomainType &x, typename Traits::RangeType &y) const
Evaluate the GridFunction at given position.
Definition: convectiondiffusionparameter.hh:289
leaf of a function tree
Definition: function.hh:302
Definition of the DUNE_DEPRECATED macro for the case that config.h is not available.
A few common exception classes.
Implements a vector constructed from a given type representing a field and a compile-time given size.
typename Impl::voider< Types... >::type void_t
Is void for all valid input types. The workhorse for C++11 SFINAE-techniques.
Definition: typetraits.hh:38
Dune namespace.
Definition: alignedallocator.hh:11
STL namespace.
Static tag representing a codimension.
Definition: dimension.hh:22
function signature for analytic functions on a grid
Definition: function.hh:642
Class to define the boundary condition types.
Definition: convectiondiffusionparameter.hh:112
Traits class for convection diffusion parameters.
Definition: convectiondiffusionparameter.hh:76
Dune::FieldVector< DomainFieldType, dimDomain-1 > IntersectionDomainType
domain type
Definition: convectiondiffusionparameter.hh:93
GV::Grid::ctype DomainFieldType
Export type for domain field.
Definition: convectiondiffusionparameter.hh:87
@ dimDomain
dimension of the domain
Definition: convectiondiffusionparameter.hh:83
Dune::FieldMatrix< RangeFieldType, dimDomain, dimDomain > PermTensorType
permeability tensor type
Definition: convectiondiffusionparameter.hh:102
Dune::FieldVector< RF, GV::dimensionworld > RangeType
range type
Definition: convectiondiffusionparameter.hh:99
GV::Traits::template Codim< 0 >::Entity ElementType
grid types
Definition: convectiondiffusionparameter.hh:105
RF RangeFieldType
Export type for range field.
Definition: convectiondiffusionparameter.hh:96
GV GridViewType
the grid view
Definition: convectiondiffusionparameter.hh:78
Dune::FieldVector< DomainFieldType, dimDomain > DomainType
domain type
Definition: convectiondiffusionparameter.hh:90
Definition: constraintsparameters.hh:26
Definition: constraintsparameters.hh:122
Dune::FieldVector< GV::Grid::ctype, GV::dimension > DomainType
domain type in dim-size coordinates
Definition: function.hh:50
traits class holding the function signature, same as in local function
Definition: function.hh:183
GV::Traits::template Codim< 0 >::Entity ElementType
codim 0 entity
Definition: function.hh:119
GV GridViewType
The type of the grid view the function lives on.
Definition: function.hh:116
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)