6 #include <dune/common/typetraits.hh>
7 #include <dune/common/version.hh>
16 template <
class T,
int N>
19 template <
class T,
int N,
int M>
25 template <
class Gr
idView>
28 using Element =
typename GridView::template Codim<0>::Entity;
29 using LocalDomain =
typename Element::Geometry::LocalCoordinate;
32 using IsGridFunction = decltype(
localFunction(std::declval<GF>()));
34 template <
class LocalFunction,
class LF = std::decay_t<LocalFunction>>
35 using IsLocalFunction = decltype((
36 std::declval<LF&>().bind(std::declval<Element>()),
37 std::declval<LF&>().unbind(),
38 std::declval<LF>()(std::declval<LocalDomain>()),
41 template <
class F,
class D>
42 using Range = std::decay_t<std::result_of_t<F(D)>>;
46 template <
class T,
int N>
47 static auto sizeOfImpl (
FieldVector<T,N>) -> std::integral_constant<int, N> {
return {}; }
49 template <
class T,
int N,
int M>
50 static auto sizeOfImpl (
FieldMatrix<T,N,M>) -> std::integral_constant<int, N*M> {
return {}; }
52 static auto sizeOfImpl (...) -> std::integral_constant<int, 1> {
return {}; }
55 static constexpr
int sizeOf () {
return decltype(sizeOfImpl(std::declval<T>()))::value; }
57 static std::vector<int> allComponents(
int n)
59 std::vector<int> components(n);
60 std::iota(components.begin(), components.end(), 0);
77 template <
class LF,
class... Args,
78 class = IsLocalFunction<LF>>
79 Function (LF&& localFct, std::string
name, std::vector<int> components, Args&&... args)
80 : localFct_(std::forward<LF>(localFct))
81 , name_(std::move(
name))
99 template <
class LF,
class... Args,
100 class = IsLocalFunction<LF>>
102 :
Function(std::forward<LF>(localFct), std::move(
name), allComponents(ncomps),
103 std::forward<Args>(args)...)
112 template <
class LF,
class... Args,
113 class = IsLocalFunction<LF>,
114 class R = Range<LF,LocalDomain> >
116 :
Function(std::forward<LF>(localFct), std::move(
name), sizeOf<R>(),
117 std::forward<Args>(args)...)
121 template <
class... Args>
124 getArg<std::string, char const*>(args..., fct.name_),
125 getArg<int,unsigned int,long,unsigned long,std::vector<int>>(args..., fct.components_),
139 template <
class GF,
class... Args,
140 disableCopyMove<Function, GF> = 0,
141 class = IsGridFunction<GF> >
156 explicit Function (std::shared_ptr<VTKFunction<GridView>
const>
const& fct, ...)
171 return self.localFct_;
175 std::string
const&
name ()
const
183 name_ = std::move(
name);
197 components_ = components;
198 localFct_.setComponents(components_);
245 std::vector<int> components_;
Definition: datacollectorinterface.hh:9
Vtk::DataTypes dataTypeOf(Dune::VTK::Precision p)
Definition: types.cc:98
RangeTypes
Type used to determine whether to limit output components to e.g. 3 (vector), or 9 (tensor)
Definition: types.hh:34
DataTypes
Definition: types.hh:50
@ FLOAT64
Definition: types.hh:57
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:17
Definition: function.hh:20
Wrapper class for functions allowing local evaluations.
Definition: function.hh:27
Function(LF &&localFct, std::string name, int ncomps, Args &&... args)
(2) Construct from a LocalFunction directly
Definition: function.hh:101
friend Vtk::LocalFunction< GridView > localFunction(Function const &self)
Create a LocalFunction.
Definition: function.hh:169
Vtk::DataTypes dataType() const
Return the VTK Datatype associated with the functions range type.
Definition: function.hh:208
void setName(std::string name)
Set the function name.
Definition: function.hh:181
Function(GF &&fct, std::string name, Args &&... args)
(5) Construct from a GridFunction
Definition: function.hh:142
Function()=default
(8) Default constructor. After construction, the function is an an invalid state.
std::string const & name() const
Return a name associated with the function.
Definition: function.hh:175
void setRangeType(Vtk::RangeTypes type, std::size_t ncomp=1)
Set the category of the range, SCALAR, VECTOR, TENSOR, or UNSPECIFIED.
Definition: function.hh:226
int numComponents() const
Return the number of components of the Range as it is written to the file.
Definition: function.hh:187
void setDataType(Vtk::DataTypes type)
Set the data-type for the components.
Definition: function.hh:214
Function(F &&fct, Vtk::FieldInfo info,...)
(6) Constructor that forwards the number of components and data type to the other constructor
Definition: function.hh:148
void setComponents(int ncomps)
Set the number of components of the Range and generate component range [0...ncomps)
Definition: function.hh:202
Function(std::shared_ptr< VTKFunction< GridView > const > const &fct,...)
(7) Construct from legacy VTKFunction
Definition: function.hh:156
Vtk::RangeTypes rangeType() const
The category of the range, SCALAR, VECTOR, TENSOR, or UNSPECIFIED.
Definition: function.hh:220
Function(Function< GridView > const &fct, Args &&... args)
(4) Construct from a Vtk::Function
Definition: function.hh:122
void setComponents(std::vector< int > components)
Set the components of the Range to visualize.
Definition: function.hh:195
Function(LF &&localFct, std::string name, Args &&... args)
(3) Construct from a LocalFunction directly.
Definition: function.hh:115
Function(LF &&localFct, std::string name, std::vector< int > components, Args &&... args)
(1) Construct from a LocalFunction directly
Definition: function.hh:79
void setFieldInfo(Vtk::FieldInfo info)
Set all the parameters from a FieldInfo object.
Definition: function.hh:234
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:245
std::string const & name() const
The name of the data field.
Definition: types.hh:239
Vtk::RangeTypes rangeType() const
Return the category of the stored range.
Definition: types.hh:251
Vtk::DataTypes dataType() const
Return the data tpe of the data field.
Definition: types.hh:257