6#include <dune/common/typetraits.hh>
7#include <dune/common/version.hh>
17 template <
class T,
int N>
20 template <
class T,
int N,
int M>
26 template <
class Gr
idView>
29 using Element =
typename GridView::template Codim<0>::Entity;
30 using LocalDomain =
typename Element::Geometry::LocalCoordinate;
32 template <
class F,
class D>
33 using Range = std::decay_t<std::result_of_t<F(D)>>;
37 template <
class T,
int N>
38 static auto sizeOfImpl (
FieldVector<T,N>) -> std::integral_constant<int, N> {
return {}; }
40 template <
class T,
int N,
int M>
41 static auto sizeOfImpl (
FieldMatrix<T,N,M>) -> std::integral_constant<int, N*M> {
return {}; }
43 static auto sizeOfImpl (...) -> std::integral_constant<int, 1> {
return {}; }
46 static constexpr int sizeOf () {
return decltype(sizeOfImpl(std::declval<T>()))::value; }
48 static std::vector<int> allComponents(
int n)
50 std::vector<int> components(n);
51 std::iota(components.begin(), components.end(), 0);
68 template <
class LF,
class... Args,
70 Function (LF&& localFct, std::string
name, std::vector<int> components, Args&&... args)
71 : localFct_(std::forward<LF>(localFct))
72 , name_(std::move(
name))
90 template <
class LF,
class... Args,
92 Function (LF&& localFct, std::string
name,
int ncomps, Args&&... args)
93 :
Function(std::forward<LF>(localFct), std::move(
name), allComponents(ncomps),
94 std::forward<Args>(args)...)
103 template <
class LF,
class... Args,
105 class R = Range<LF,LocalDomain> >
107 :
Function(std::forward<LF>(localFct), std::move(
name), sizeOf<R>(),
108 std::forward<Args>(args)...)
112 template <
class... Args>
115 getArg<std::string, char const*>(args..., fct.name_),
116 getArg<int,unsigned int,long,unsigned long,std::vector<int>>(args..., fct.components_),
130 template <
class GF,
class... Args,
131 disableCopyMove<Function, GF> = 0,
144 template <
class F,
class... Args,
145 disableCopyMove<Function, F> = 0,
147 class =
decltype(std::declval<F>().name()),
148 class =
decltype(std::declval<F>().numComponents()),
149 class =
decltype(std::declval<F>().dataType()) >
159 explicit Function (std::shared_ptr<VTKFunction<GridView>
const>
const& fct, ...)
174 return self.localFct_;
178 std::string
const&
name ()
const
186 name_ = std::move(
name);
200 components_ = components;
201 localFct_.setComponents(components_);
248 std::vector<int> components_;
decltype((std::declval< LF & >().bind(std::declval< LocalContext >()), std::declval< LF & >().unbind(), std::declval< LF >()(std::declval< typename LocalContext::Geometry::LocalCoordinate >()), true)) IsLocalFunction
Definition: concepts.hh:39
Vtk::DataTypes dataTypeOf()
Definition: types.hh:79
decltype((localFunction(std::declval< GF const & >()), true)) IsGridFunction
Definition: concepts.hh:32
RangeTypes
Type used to determine whether to limit output components to e.g. 3 (vector), or 9 (tensor)
Definition: types.hh:35
DataTypes
Definition: types.hh:52
decltype(auto) getArg(Arg0 &&arg0, Args &&... args)
Definition: arguments.hh:29
Vtk::RangeTypes rangeTypeOf(Dune::VTK::FieldInfo::Type t)
Definition: types.cc:56
Definition: function.hh:18
Definition: function.hh:21
Wrapper class for functions allowing local evaluations.
Definition: function.hh:28
Function(LF &&localFct, std::string name, int ncomps, Args &&... args)
(2) Construct from a LocalFunction directly
Definition: function.hh:92
Vtk::DataTypes dataType() const
Return the VTK Datatype associated with the functions range type.
Definition: function.hh:211
void setName(std::string name)
Set the function name.
Definition: function.hh:184
friend Vtk::LocalFunction< GridView > localFunction(Function const &self)
Create a LocalFunction.
Definition: function.hh:172
std::string const & name() const
Return a name associated with the function.
Definition: function.hh:178
Function(GF &&fct, std::string name, Args &&... args)
(5) Construct from a GridFunction
Definition: function.hh:133
Function()=default
(9) Default constructor. After construction, the function is an an invalid state.
void setRangeType(Vtk::RangeTypes type, std::size_t ncomp=1)
Set the category of the range, SCALAR, VECTOR, TENSOR, or UNSPECIFIED.
Definition: function.hh:229
int numComponents() const
Return the number of components of the Range as it is written to the file.
Definition: function.hh:190
void setDataType(Vtk::DataTypes type)
Set the data-type for the components.
Definition: function.hh:217
Function(F &&fct, Vtk::FieldInfo info,...)
(6) Constructor that forwards the number of components and data type to the other constructor
Definition: function.hh:139
void setComponents(int ncomps)
Set the number of components of the Range and generate component range [0...ncomps)
Definition: function.hh:205
Function(std::shared_ptr< VTKFunction< GridView > const > const &fct,...)
(8) Construct from legacy VTKFunction
Definition: function.hh:159
Vtk::RangeTypes rangeType() const
The category of the range, SCALAR, VECTOR, TENSOR, or UNSPECIFIED.
Definition: function.hh:223
Function(Function< GridView > const &fct, Args &&... args)
(4) Construct from a Vtk::Function
Definition: function.hh:113
void setComponents(std::vector< int > components)
Set the components of the Range to visualize.
Definition: function.hh:198
Function(LF &&localFct, std::string name, Args &&... args)
(3) Construct from a LocalFunction directly.
Definition: function.hh:106
Function(LF &&localFct, std::string name, std::vector< int > components, Args &&... args)
(1) Construct from a LocalFunction directly
Definition: function.hh:70
void setFieldInfo(Vtk::FieldInfo info)
Set all the parameters from a FieldInfo object.
Definition: function.hh:237
A Vtk::LocalFunction is a function-like object that can be bound to a grid element an that provides a...
Definition: localfunction.hh:23
int size() const
The number of components in the data field.
Definition: types.hh:253
std::string const & name() const
The name of the data field.
Definition: types.hh:247
Vtk::RangeTypes rangeType() const
Return the category of the stored range.
Definition: types.hh:259
Vtk::DataTypes dataType() const
Return the data tpe of the data field.
Definition: types.hh:265