13#ifndef FUNCTIONVIEWS_HH
14#define FUNCTIONVIEWS_HH
58 namespace FunctionViews
68 template <
typename Function>
72 typedef typename Function::Space
Space;
81 int order(SFS
const& sfs)
const {
return 2*f.
order(sfs); }
83 template<
class EPtr,
class V>
86 return f.
value(ci,v).two_norm2();
91 return f.
value(evaluator).two_norm2();
113 template <
typename Function,
bool asVector=Function::components==1>
119 typedef typename Function::Space
Space;
120 static int const dim = Space::dim;
125 using ValueType = std::conditional_t<asVector, Dune::FieldVector<Scalar,dim*m>,
148 return convertToValueType(f.derivative(evaluator));
155 template<
class Cell,
class Position>
158 return convertToValueType(f.derivative(cell,xi));
182 template <
class Function>
193 template <
class Function>
207 template <
typename Function>
211 typedef typename Function::Space
Space;
231 return f.derivative(evaluator).frobenius_norm2();
234 template <
class Cell,
class Position>
237 f.derivative(cell,xi).frobenius_norm2();
251 template<
template<
class Fct1,
class Fct2>
class View,
class Fct1,
class Fct2>
252 View<Fct1, Fct2>
makeView(Fct1
const& f1, Fct2
const& f2)
254 return View<Fct1,Fct2>(f1,f2);
261 template<
template<
typename Function>
class View,
typename Function>
264 return View<Function>(f1);
285 namespace WeakFunctionViews
291 int order(SFS
const& sfs)
const {
return 0; }
301 int order(SFS
const& sfs)
const {
return 0; }
303 template<
class EPtr,
class V>
316 template<
class EPtr,
class V>
317 double value(EPtr
const& ci, V
const&)
const
319 return ci.mightBeCoarsened()/ci.geometry().volume();
326 template<
class EPtr,
class V>
327 double value(EPtr
const& ci, V
const&)
const
329 return ci.wasRefined()/ci.geometry().volume();
336 template<
class EPtr,
class V>
337 double value(EPtr
const& ci, V
const&)
const
339 return ci.level()/ci.geometry().volume();
346 template<
class EPtr,
class V>
347 double value(EPtr
const& ci, V
const&)
const
349 return ci.isRegular()/ci.geometry().volume();
355 template<
class Gr
id,
class T>
class CellData;
356 template<
class Space>
class LocalIntegral;
359 template<
class View,
class Space>
364 errorIndicator(localIntegral(View(),space));
365 return errorIndicator;
376 template <
class Function>
389 template <
class Cell>
398 template <
class Cell>
402 return f1.
value(cell,localCoordinate) * f2.
value(cell,localCoordinate);
Function is the interface for evaluatable functions .
int order(Cell const &cell) const
static int const components
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.
AbsSquare(Function const &f_)
Dune::FieldVector< Scalar, 1 > ValueType
Space const & space() const
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
Space const & space() const
GradientAbsSquare(Function const &f_)
int order(SFS const &sfs) 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
Space const & space() const
ValueType value(typename Function::Space::Evaluator const &evaluator) const
evaluates the derivative
Gradient(Function const &f_)
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.
Dune::FieldVector< T, n > max(Dune::FieldVector< T, n > x, Dune::FieldVector< T, n > const &y)
Componentwise maximum.
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.
CellData< typenameSpace::Grid, Dune::FieldVector< double, 1 > >::CellDataVector evalCellProperties(Space const &space)
Evaluate WeakFunctionViews and construct CellData.
int order(SFS const &sfs) const
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()
int order(SFS const &sfs) const