DUNE-ACFEM (2.5.1)

basicfunctions.hh
1 #ifndef __DUNE_ACFEM_BASICFUNCTIONS_HH__
2 #define __DUNE_ACFEM_BASICFUNCTIONS_HH__
3 
4 #include <dune/fem/function/localfunction/mutable.hh>
5 #include <dune/fem/function/localfunction/const.hh>
6 
7 #include "../functions/gridfunctionexpression.hh"
8 #include "../functions/boundaryfunctionexpression.hh"
9 
10 #include "../functions/constantfunction.hh"
11 #include "../functions/coordinatefunctions.hh"
12 #include "../functions/divergencefunction.hh"
13 #include "../functions/gradientfunction.hh"
14 #include "../functions/piecewisefunction.hh"
15 #include "../functions/righthandsidefunction.hh"
16 
17 namespace Dune {
18 
19  namespace ACFem {
20 
60  template<class GridPart>
61  static inline
62  auto identityGridFunction(const GridPart& gridPart)
64  -> decltype(wrapToGridFunction(std::declval<std::string>(),
65  std::declval<IdentityFunction<typename GridPart::ctype,
66  GridPart::dimensionworld> >(),
67  gridPart,
68  1))
70  {
72 
73  return wrapToGridFunction("x -> x", FunctionType(), gridPart, 1);
74  }
75 
81  template<unsigned N>
82  struct X
83  : public std::integral_constant<unsigned, N>
84  {};
85 
89  template<class GridPart, unsigned N>
90  static inline
91  auto coordinateGridFunction(const GridPart& gridPart, const std::integral_constant<unsigned, N>& dummy)
93  -> decltype(wrapToGridFunction(std::declval<std::string>(),
94  std::declval<CoordinateFunction<typename GridPart::ctype,
95  GridPart::dimensionworld, N> >(),
96  gridPart,
97  1))
99  {
101 
102  // could use stringstream here.
103  return wrapToGridFunction("(x -> x_" + std::to_string(N) + ")", FunctionType(), gridPart, 1);
104  }
105 
107  template<class GridPart, unsigned N>
108  static inline
109  ConstantGridFunction<Fem::FunctionSpace<typename GridPart::ctype,
110  typename GridPart::ctype,
111  GridPart::dimensionworld,
112  GridPart::dimensionworld>,
113  GridPart>
114  coordinateVectorFunction(const GridPart& gridPart, const std::integral_constant<unsigned, N>& dummy)
115  {
116  const int dimWorld = GridPart::dimensionworld;
117  typedef typename GridPart::ctype FieldType;
118  typedef Fem::FunctionSpace<FieldType, FieldType, dimWorld, dimWorld> FunctionSpaceType;
120  typename FunctionSpaceType::RangeType e_N;
121 
122  e_N = 0.;
123  e_N[N] = 1.0;
124 
125  return FunctionType(e_N, gridPart);
126  }
127 
129  template<class GridFunction>
130  LocalFunctionWrapper<LocalDivergenceAdapter<GridFunction>, typename GridFunction::GridPartType>
131  divergence(const Fem::Function<typename GridFunction::FunctionSpaceType, GridFunction>& f_, const std::string& name = "")
132  {
133  typedef LocalDivergenceAdapter<GridFunction> LocalFunctionType;
135 
136  const GridFunction& f(static_cast<const GridFunction&>(f_));
137  LocalFunctionType div(f, name);
138 
139  return GridFunctionType(div.name(), div, f.space().gridPart(), f.space().order() - 1);
140  }
141 
143  template<class GridFunction>
144  LocalFunctionWrapper<LocalGradientAdapter<GridFunction>, typename GridFunction::GridPartType>
145  gradient(const Fem::Function<typename GridFunction::FunctionSpaceType, GridFunction>& f_, const std::string& name = "")
146  {
147  typedef LocalGradientAdapter<GridFunction> LocalFunctionType;
149 
150  const GridFunction& f(static_cast<const GridFunction&>(f_));
151  LocalFunctionType div(f, name);
152 
153  return GridFunctionType(div.name(), div, f.space().gridPart(), f.space().order() - 1);
154  }
155 
175  template<class CharacteristicFunction, class GridFunction0, class GridFunction1>
176  LocalFunctionWrapper<LocalPiecewiseAdapter<CharacteristicFunction, GridFunction0, GridFunction1>, typename GridFunction0::GridPartType>
177  piecewise(const Fem::Function<typename CharacteristicFunction::FunctionSpaceType, CharacteristicFunction>& characteristic,
178  const Fem::Function<typename GridFunction0::FunctionSpaceType, GridFunction0>& function0,
179  const Fem::Function<typename GridFunction1::FunctionSpaceType, GridFunction1>& function1,
180  const std::string& name = "")
181  {
184 
185  const auto& realFunction0 = static_cast<const GridFunction0&>(function0);
186  const auto& realFunction1 = static_cast<const GridFunction1&>(function1);
187  const auto& gridPart = realFunction0.space().gridPart();
188 
189  const auto& c = static_cast<const CharacteristicFunction&>(characteristic);
190 
191  LocalFunctionType pw(c, function0, function1, name);
192 
193  return GridFunctionType(pw.name(), pw, gridPart, std::max(realFunction0.space().order(), realFunction1.space().order()));
194  }
195 
196 
198  template<class CharacteristicFunction, class GridFunction1>
199  LocalFunctionWrapper<LocalPiecewiseAdapter<CharacteristicFunction, ConstantGridFunction<typename GridFunction1::FunctionSpaceType, typename GridFunction1::GridPartType>, GridFunction1>, typename GridFunction1::GridPartType>
200  piecewise(const Fem::Function<typename CharacteristicFunction::FunctionSpaceType, CharacteristicFunction>& characteristic,
201  const typename GridFunction1::FunctionSpaceType::RangeType value0,
202  const Fem::Function<typename GridFunction1::FunctionSpaceType, GridFunction1>& function1,
203  const std::string& name = "")
204  {
208 
209  const auto& realFunction1 = static_cast<const GridFunction1&>(function1);
210  const auto& gridPart = realFunction1.space().gridPart();
211 
212  const GridFunction0 realFunction0(value0, gridPart);
213 
214  const auto& c = static_cast<const CharacteristicFunction&>(characteristic);
215 
216  LocalFunctionType pw(c, realFunction0, function1, name);
217 
218  return GridFunctionType(pw.name(), pw, gridPart, std::max(realFunction0.space().order(), realFunction1.space().order()));
219  }
220 
222  template<class CharacteristicFunction, class GridFunction0>
223  LocalFunctionWrapper<LocalPiecewiseAdapter<CharacteristicFunction, GridFunction0, ConstantGridFunction<typename GridFunction0::FunctionSpaceType, typename GridFunction0::GridPartType> >, typename GridFunction0::GridPartType>
224  piecewise(const Fem::Function<typename CharacteristicFunction::FunctionSpaceType, CharacteristicFunction>& characteristic,
225  const Fem::Function<typename GridFunction0::FunctionSpaceType, GridFunction0>& function0,
226  const typename GridFunction0::FunctionSpaceType::RangeType value1,
227  const std::string& name = "")
228  {
232 
233  const auto& realFunction0 = static_cast<const GridFunction0&>(function0);
234  const auto& gridPart = realFunction0.space().gridPart();
235  const GridFunction1 realFunction1(value1, gridPart);
236 
237  const auto& c = static_cast<const CharacteristicFunction&>(characteristic);
238 
239  LocalFunctionType pw(c, function0, realFunction1, name);
240 
241  return GridFunctionType(pw.name(), pw, gridPart, std::max(realFunction0.space().order(), realFunction1.space().order()));
242  }
244  template<class CharacteristicFunction, class FunctionSpaceType>
245  LocalFunctionWrapper<LocalPiecewiseAdapter<CharacteristicFunction, ConstantGridFunction<FunctionSpaceType, typename CharacteristicFunction::GridPartType>, ConstantGridFunction<FunctionSpaceType, typename CharacteristicFunction::GridPartType> >, typename CharacteristicFunction::GridPartType>
246  piecewise(const Fem::Function<typename CharacteristicFunction::FunctionSpaceType, CharacteristicFunction>& characteristic,
247  const typename FunctionSpaceType::RangeType value0,
248  const typename FunctionSpaceType::RangeType value1,
249  const std::string& name = "")
250  {
255 
256  const auto& c = static_cast<const CharacteristicFunction&>(characteristic);
257 
258  const auto& gridPart = c.space().gridPart();
259  const GridFunction0 realFunction0(value0, gridPart);
260  const GridFunction1 realFunction1(value1, gridPart);
261 
262  LocalFunctionType pw(c, realFunction0, realFunction1, name);
263 
264  return GridFunctionType(pw.name(), pw, gridPart, std::max(realFunction0.space().order(), realFunction1.space().order()));
265  }
266 
267  template<class Traits>
268  typename Fem::ConstLocalFunction<typename Traits::DiscreteFunctionType>
269  localFunction(const Fem::DiscreteFunctionInterface<Traits>& f_)
270  {
271  typedef typename Traits::DiscreteFunctionType DiscreteFunctionType;
272  const auto& f = static_cast<const DiscreteFunctionType&>(f_);
273 
274  return Fem::ConstLocalFunction<DiscreteFunctionType>(f);
275  }
276 
277  template<class Traits>
278  typename Fem::MutableLocalFunction<typename Traits::DiscreteFunctionType>
279  localFunction(Fem::DiscreteFunctionInterface<Traits>& f_)
280  {
281  typedef typename Traits::DiscreteFunctionType DiscreteFunctionType;
282  auto& f = static_cast<DiscreteFunctionType&>(f_);
283 
284  return Fem::MutableLocalFunction<DiscreteFunctionType>(f);
285  }
286 
287 #if 0
288  // Does not make sense without special "const" object.
289  template<class GridFunction>
290  const typename GridFunction::LocalFunctionType
291  localFunction(const Fem::Function<typename GridFunction::FunctionSpaceType, GridFunction>& f_)
292  {
293  typedef typename GridFunction::LocalFunctionType LocalFunctionType;
294  const auto& f = asImp(f_);
295 
296  return LocalFunctionType(f);
297  }
298 #endif
299 
300  template<class GridFunction>
301  typename GridFunction::LocalFunctionType
302  localFunction(const Fem::Function<typename GridFunction::FunctionSpaceType, GridFunction>& f_)
303  {
304  typedef typename GridFunction::LocalFunctionType LocalFunctionType;
305  auto& f = asImp(f_);
306 
307  return LocalFunctionType(f);
308  }
309 
311 
313 
315 
316  } // ACFem::
317 
318 } // Dune::
319 
320 #endif // __DUNE_ACFEM_BASICFUNCTIONS_HH__
ConstantGridFunction implements a constant function.
Definition: constantfunction.hh:108
Coordinate function returns a specific component of the world point, x_N, so to say.
Definition: coordinatefunctions.hh:63
A function which returns its argument.
Definition: coordinatefunctions.hh:19
An adapter class to compute the divergence of a GridFunction.
Definition: divergencefunction.hh:21
LocalFunctionWrapper wraps a class with a local evaluate method into a grid function.
Definition: localfunctionwrapper.hh:71
An adapter class to compute the gradient of a scalar GridFunction.
Definition: gradientfunction.hh:25
An adapter class to form a piecewise defined function where the domain decomposition is given by the ...
Definition: piecewisefunction.hh:34
const Implementation & asImp(const Fem::BartonNackmanInterface< Interface, Implementation > &arg)
Up-cast to the implementation for any Fem::BartonNackmanInterface.
Definition: expressionoperations.hh:71
static auto wrapToGridFunction(const std::string &name, const FunctionImp &f, const GridPart &g, unsigned order=111)
Possibly wrap a function into a GridFunctionWrapper.
Definition: gridfunctionwrapper.hh:342
static auto coordinateGridFunction(const GridPart &gridPart, const std::integral_constant< unsigned, N > &dummy)
A grid-function which returns the N-th component of the world coordinates of its argument.
Definition: basicfunctions.hh:91
LocalFunctionWrapper< LocalPiecewiseAdapter< CharacteristicFunction, GridFunction0, GridFunction1 >, typename GridFunction0::GridPartType > piecewise(const Fem::Function< typename CharacteristicFunction::FunctionSpaceType, CharacteristicFunction > &characteristic, const Fem::Function< typename GridFunction0::FunctionSpaceType, GridFunction0 > &function0, const Fem::Function< typename GridFunction1::FunctionSpaceType, GridFunction1 > &function1, const std::string &name="")
Form a piecewise defined function where the domain decomposition is defined by the indicator function...
Definition: basicfunctions.hh:177
static ConstantGridFunction< Fem::FunctionSpace< typename GridPart::ctype, typename GridPart::ctype, GridPart::dimensionworld, GridPart::dimensionworld >, GridPart > coordinateVectorFunction(const GridPart &gridPart, const std::integral_constant< unsigned, N > &dummy)
Construct a unit coordinate vector.
Definition: basicfunctions.hh:114
static auto identityGridFunction(const GridPart &gridPart)
A grid-function which returns the world coordinates of its argument.
Definition: basicfunctions.hh:62
LocalFunctionWrapper< LocalDivergenceAdapter< GridFunction >, typename GridFunction::GridPartType > divergence(const Fem::Function< typename GridFunction::FunctionSpaceType, GridFunction > &f_, const std::string &name="")
Take the divergence of a given function.
Definition: basicfunctions.hh:131
LocalFunctionWrapper< LocalGradientAdapter< GridFunction >, typename GridFunction::GridPartType > gradient(const Fem::Function< typename GridFunction::FunctionSpaceType, GridFunction > &f_, const std::string &name="")
Take the gradient of a given function.
Definition: basicfunctions.hh:145
LocalFunctionWrapper< LocalMaxAdapter< GridFunction1, GridFunction2 >, typename GridFunction1::GridPartType > max(const Fem::Function< typename GridFunction1::FunctionSpaceType, GridFunction1 > &f1, const Fem::Function< typename GridFunction2::FunctionSpaceType, GridFunction2 > &f2, const std::string &name="")
Pointwise maximum of two given functions.
Definition: maxfunction.hh:121
Used as "tag" for coordinateGridFunction() and coordinateVectorFunction().
Definition: basicfunctions.hh:84
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)