KASKADE 7 development version
functionviews.hh
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the library KASKADE 7 */
4/* https://www.zib.de/research/projects/kaskade7-finite-element-toolbox */
5/* */
6/* Copyright (C) 2002-2021 Zuse Institute Berlin */
7/* */
8/* KASKADE 7 is distributed under the terms of the ZIB Academic License. */
9/* see $KASKADE/academic.txt */
10/* */
11/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12
13#ifndef FUNCTIONVIEWS_HH
14#define FUNCTIONVIEWS_HH
15
16#include <type_traits>
17
18namespace Kaskade
19{
58 namespace FunctionViews
59 {
68 template <typename Function>
70 {
71 public:
72 typedef typename Function::Space Space;
73 typedef typename Function::Scalar Scalar;
75
76 AbsSquare(Function const& f_) : f(f_) {}
77
78 Space const& space() const {return f.space();}
79
80 template<class SFS>
81 int order(SFS const& sfs) const {return 2*f.order(sfs); }
82
83 template<class EPtr, class V>
84 ValueType value(EPtr const& ci, V const& v) const
85 {
86 return f.value(ci,v).two_norm2();
87 }
88
89 ValueType value(typename Function::Space::Evaluator const& evaluator) const
90 {
91 return f.value(evaluator).two_norm2();
92 }
93
94 private:
95 Function const& f;
96 };
97
98 // --------------------------------------------------------------------------------------------
99
113 template <typename Function, bool asVector=Function::components==1>
115 {
116 static int const m = Function::components;
117
118 public:
119 typedef typename Function::Space Space;
120 static int const dim = Space::dim;
121
122 typedef typename Function::Scalar Scalar;
123
124
125 using ValueType = std::conditional_t<asVector, Dune::FieldVector<Scalar,dim*m>,
126 /*otherwise*/ Dune::FieldMatrix<Scalar,m,dim>>;
127
128 Gradient(Function const& f_) : f(f_)
129 {}
130
131 Space const& space() const
132 {
133 return f.space();
134 }
135
136 template<class SFS>
137 int order(SFS const& sfs) const
138 {
139 return std::max(f.order(sfs)-1,0);
140 }
141
146 ValueType value(typename Function::Space::Evaluator const& evaluator) const
147 {
148 return convertToValueType(f.derivative(evaluator));
149 }
150
155 template<class Cell, class Position>
156 ValueType value(Cell const& cell, Position const& xi) const
157 {
158 return convertToValueType(f.derivative(cell,xi));
159 }
160
161
162 private:
163 Function const& f;
164
165 ValueType convertToValueType(Dune::FieldMatrix<Scalar,m,dim> const& fx) const
166 {
167 if constexpr (m==1)
168 return fx[0];
169 else if constexpr (asVector)
170 return Dune::asVector(fx);
171 else
172 return fx;
173 }
174 };
175
176
182 template <class Function>
184 {
185 return Gradient<Function,false>(f);
186 }
187
193 template <class Function>
195 {
196 return Gradient<Function,true>(f);
197 }
198
199 // --------------------------------------------------------------------------------------------
200
207 template <typename Function>
209 {
210 public:
211 typedef typename Function::Space Space;
212 typedef typename Function::Scalar Scalar;
214
215 GradientAbsSquare(Function const& f_) : f(f_)
216 { }
217
218 Space const& space() const
219 {
220 return f.space();
221 }
222
223 template <class SFS>
224 int order(SFS const& sfs) const
225 {
226 return 2*std::max(f.order(sfs)-1,0);
227 }
228
229 ValueType value(typename Function::Space::Evaluator const& evaluator) const
230 {
231 return f.derivative(evaluator).frobenius_norm2();
232 }
233
234 template <class Cell, class Position>
235 ValueType value(Cell const& cell, Position const& xi) const
236 {
237 f.derivative(cell,xi).frobenius_norm2();
238 }
239
240 private:
241 Function const& f;
242 };
243
244
245 }
246
251 template<template<class Fct1,class Fct2> class View, class Fct1, class Fct2>
252 View<Fct1, Fct2> makeView(Fct1 const& f1, Fct2 const& f2)
253 {
254 return View<Fct1,Fct2>(f1,f2);
255 }
256
261 template<template<typename Function> class View, typename Function>
262 View<Function> makeView(Function const& f1)
263 {
264 return View<Function>(f1);
265 }
266
285 namespace WeakFunctionViews
286 {
289 {
290 template<class SFS>
291 int order(SFS const& sfs) const { return 0; }
292
293 virtual ~ZeroOrderView() {};
294 };
295
297 {
299
300 template<class SFS>
301 int order(SFS const& sfs) const { return 0; }
302
303 template<class EPtr, class V>
304 ValueType value(EPtr const& ci, V const&) const
305 {
306 ValueType v;
307 v[0]=0;
308 v[1]=1;
309 return v;
310 }
311 };
312
315 {
316 template<class EPtr, class V>
317 double value(EPtr const& ci, V const&) const
318 {
319 return ci.mightBeCoarsened()/ci.geometry().volume();
320 }
321 };
322
324 struct WasRefined : public ZeroOrderView
325 {
326 template<class EPtr, class V>
327 double value(EPtr const& ci, V const&) const
328 {
329 return ci.wasRefined()/ci.geometry().volume();
330 }
331 };
332
334 struct GridLevel : public ZeroOrderView
335 {
336 template<class EPtr, class V>
337 double value(EPtr const& ci, V const&) const
338 {
339 return ci.level()/ci.geometry().volume();
340 }
341 };
342
344 struct IsRegular : public ZeroOrderView
345 {
346 template<class EPtr, class V>
347 double value(EPtr const& ci, V const&) const
348 {
349 return ci.isRegular()/ci.geometry().volume();
350 }
351 };
352
353 }
354
355 template<class Grid, class T> class CellData;
356 template<class Space> class LocalIntegral;
357
359 template<class View, class Space>
361 {
362 LocalIntegral<Space> localIntegral;
364 errorIndicator(localIntegral(View(),space));
365 return errorIndicator;
366 }
367
376 template <class Function>
378 {
379 public:
383 ScalarProductView(Function const& f1_, Function const& f2_): f1(f1_), f2(f2_) {}
384
385 typedef typename Function::Scalar Scalar;
388
389 template <class Cell>
390 int order(Cell const& cell) const
391 {
392 if (std::max(f1.order(cell),f2.order(cell)) < std::numeric_limits<int>::max())
393 return f1.order(cell)+f2.order(cell);
394 else
396 }
397
398 template <class Cell>
399 ValueType value(Cell const& cell,
401 {
402 return f1.value(cell,localCoordinate) * f2.value(cell,localCoordinate);
403 }
404
405
406 private:
407 Function const& f1;
408 Function const& f2;
409 };
410} /* end of namespace Kaskade */
411
412#endif
Function is the interface for evaluatable functions .
Definition: concepts.hh:324
int order(Cell const &cell) const
static int const components
Definition: concepts.hh:335
unspecified Scalar
Definition: concepts.hh:330
ValueType value(Cell const &cell, Dune::FieldVector< typename Cell::Geometry::ctype, Cell::dimension > const &localCoordinate) const
Class that stores information for each cell of a grid.
Definition: celldata.hh:37
Dune::FieldVector< Scalar, 1 > ValueType
ValueType value(EPtr const &ci, V const &v) const
int order(SFS const &sfs) const
ValueType value(typename Function::Space::Evaluator const &evaluator) const
A function view returning the square of the Frobenius norm of a function's derivative.
ValueType value(Cell const &cell, Position const &xi) const
ValueType value(typename Function::Space::Evaluator const &evaluator) const
Derivative of a given finite element function.
ValueType value(Cell const &cell, Position const &xi) const
evaluates the derivative
ValueType value(typename Function::Space::Evaluator const &evaluator) const
evaluates the derivative
std::conditional_t< asVector, Dune::FieldVector< Scalar, dim *m >, Dune::FieldMatrix< Scalar, m, dim > > ValueType
int order(SFS const &sfs) const
Create a CellData by computing local integrals over each Cell.
(Scalar) product of two functions.
Dune::FieldVector< Scalar, components > ValueType
static int const components
ValueType value(Cell const &cell, Dune::FieldVector< typename Cell::Geometry::ctype, Cell::dimension > const &localCoordinate) const
ScalarProductView(Function const &f1_, Function const &f2_)
int order(Cell const &cell) const
Gradient< Function, false > gradient(Function const &f)
creates a function view of the gradient
Gradient< Function, true > gradientAsVector(Function const &f)
creates a function view of the gradient
View< Fct1, Fct2 > makeView(Fct1 const &f1, Fct2 const &f2)
Construct a View on two Functions without having to state the type of these explicitely.
typename GridView::template Codim< 0 >::Entity Cell
The type of cells (entities of codimension 0) in the grid view.
Definition: gridBasics.hh:35
Dune::FieldVector< T, n > max(Dune::FieldVector< T, n > x, Dune::FieldVector< T, n > const &y)
Componentwise maximum.
Definition: fixdune.hh:110
Dune::FieldVector< T, n *m > asVector(Dune::FieldMatrix< T, n, m > const &x)
Converts a matrix of size n x m to a vector of size n*m by concatenating all the columns.
Definition: fixdune.hh:135
CellData< typenameSpace::Grid, Dune::FieldVector< double, 1 > >::CellDataVector evalCellProperties(Space const &space)
Evaluate WeakFunctionViews and construct CellData.
Dune::FieldVector< double, 2 > ValueType
ValueType value(EPtr const &ci, V const &) const
Get level() of entity in Grid.
double value(EPtr const &ci, V const &) const
Return 1, if entity is constructed from red refinement.
double value(EPtr const &ci, V const &) const
Return 1, if entity mightbecoarsened.
double value(EPtr const &ci, V const &) const
Return 1, if entity was refined.
double value(EPtr const &ci, V const &) const
Base class providing int order()