CoordsGridFunction.hpp 3.14 KB
Newer Older
Praetorius, Simon's avatar
Praetorius, Simon committed
1
2
3
4
#pragma once

#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
5
#include <dune/common/diagonalmatrix.hh>
Praetorius, Simon's avatar
Praetorius, Simon committed
6
7
#include <dune/common/typeutilities.hh>

8
9
#include <amdis/gridfunctions/AnalyticGridFunction.hpp>

Praetorius, Simon's avatar
Praetorius, Simon committed
10
11
namespace AMDiS
{
12
  namespace Operation
Praetorius, Simon's avatar
Praetorius, Simon committed
13
  {
14
15
16
17
18
19
    /** \addtogroup operations
     *  @{
     **/

    /// A functor that evaluates to the global coordinates
    struct CoordsFunction
Praetorius, Simon's avatar
Praetorius, Simon committed
20
    {
21
22
23
24
25
26
27
28
29
30
31
      template <class T, int N>
      Dune::FieldVector<T, N> const& operator()(Dune::FieldVector<T, N> const& x) const
      {
        return x;
      }

      friend int order(CoordsFunction const& /*f*/, int /*d*/)
      {
        return 1;
      }

32
33
34
35
36
37
38
39
40
41
42
43
      struct Derivative
      {
        template <class T, int N>
        Dune::DiagonalMatrix<T, N> const& operator()(Dune::FieldVector<T, N> const& x) const
        {
          return Dune::DiagonalMatrix<T,N>{T(1)};
        }
      };
      friend Derivative partial(CoordsFunction const& /*f*/, index_t<0>)
      {
        return Derivative{};
      }
Praetorius, Simon's avatar
Praetorius, Simon committed
44
45
46
    };


47
48
    /// A functor that evaluates to a component of the global coordinates
    struct CoordsCompFunction
Praetorius, Simon's avatar
Praetorius, Simon committed
49
    {
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
      /// Constructor. Stores the component `comp` of the coordinates.
      explicit CoordsCompFunction(int comp)
        : comp_(comp)
      {}

      template <class T, int N>
      T const& operator()(Dune::FieldVector<T, N> const& x) const
      {
        return x[comp_];
      }

      friend int order(CoordsCompFunction const& /*f*/, int /*d*/)
      {
        return 1;
      }

66
67
68
69
70
      struct Derivative
      {
        explicit Derivative(int comp)
          : comp_(comp)
        {}
Praetorius, Simon's avatar
Praetorius, Simon committed
71

72
73
74
75
76
        template <class T, int N>
        Dune::FieldVector<T, N> operator()(Dune::FieldVector<T, N> const& x) const
        {
          Dune::FieldVector<T, N> result(0);
          result[comp_] = T(1);
Praetorius, Simon's avatar
Praetorius, Simon committed
77

78
79
          return result;
        }
Praetorius, Simon's avatar
Praetorius, Simon committed
80

81
82
83
84
85
      private:
        int comp_;
      };

      friend Derivative partial(CoordsCompFunction const& f, index_t<0>)
Praetorius, Simon's avatar
Praetorius, Simon committed
86
      {
87
88
        return Derivative{f.comp_};
      }
Praetorius, Simon's avatar
Praetorius, Simon committed
89
90
91
92
93

    private:
      int comp_;
    };

94
    /** @} **/
Praetorius, Simon's avatar
Praetorius, Simon committed
95

96
  } // end namespace Operation
Praetorius, Simon's avatar
Praetorius, Simon committed
97
98
99
100

  /// Generator for \ref CoordsFunction
  inline auto X()
  {
101
    return Operation::CoordsFunction{};
Praetorius, Simon's avatar
Praetorius, Simon committed
102
103
104
105
106
  }

  /// Generator for \ref CoordsCompFunction
  inline auto X(int comp)
  {
107
    return Operation::CoordsCompFunction{comp};
Praetorius, Simon's avatar
Praetorius, Simon committed
108
109
  }

110
111
112
113
114
115
116
117
118
119
120
  namespace Traits
  {
    template <>
    struct IsPreGridFunction<Operation::CoordsFunction>
      : std::true_type {};

    template <>
    struct IsPreGridFunction<Operation::CoordsCompFunction>
      : std::true_type {};
  }

121
122
  template <>
  struct GridFunctionCreator<Operation::CoordsFunction>
Praetorius, Simon's avatar
Praetorius, Simon committed
123
  {
124
125
    using PreGridFct = Operation::CoordsFunction;

126
    template <class GridView>
127
    static auto create(PreGridFct const& f, GridView const& gridView)
128
129
130
    {
      return makeAnalyticGridFunction(f, gridView);
    }
131
132
133
134
135
136
  };

  template <>
  struct GridFunctionCreator<Operation::CoordsCompFunction>
  {
    using PreGridFct = Operation::CoordsCompFunction;
Praetorius, Simon's avatar
Praetorius, Simon committed
137

138
    template <class GridView>
139
    static auto create(PreGridFct const& f, GridView const& gridView)
140
141
142
    {
      return makeAnalyticGridFunction(f, gridView);
    }
143
  };
Praetorius, Simon's avatar
Praetorius, Simon committed
144
145

} // end namespace AMDiS