FirstOrderTestPartialTrial.hpp 4.11 KB
Newer Older
1
2
3
4
#pragma once

#include <type_traits>

5
#include <amdis/GridFunctionOperator.hpp>
6
#include <amdis/utility/LocalBasisCache.hpp>
7
#include <amdis/Output.hpp>
8
#include <amdis/common/StaticSize.hpp>
9
10
11

namespace AMDiS
{
12
13
14
15
16
  /**
   * \addtogroup operators
   * @{
   **/

17
18
19
20
  namespace tag
  {
    struct test_partialtrial
    {
21
      std::size_t comp;
22
    };
23
24
25
26
27

    struct partial_trial
    {
      std::size_t comp;
    };
28
29
30
  }


31
  /// first-order operator \f$ \langle\psi, c\,\partial_i\phi\rangle \f$
32
33
34
  template <class LC, class GridFct>
  class GridFunctionOperator<tag::test_partialtrial, LC, GridFct>
      : public GridFunctionOperatorBase<GridFunctionOperator<tag::test_partialtrial, LC, GridFct>, LC, GridFct>
35
  {
36
    using Self = GridFunctionOperator;
37
    using Super = GridFunctionOperatorBase<Self, LC, GridFct>;
38

39
    static_assert( Size_v<typename GridFct::Range> == 1, "Expression must be of scalar type." );
40
41

  public:
42
43
    GridFunctionOperator(tag::test_partialtrial tag, GridFct const& expr)
      : Super(expr, 1)
44
45
46
      , comp_(tag.comp)
    {}

47
48
    template <class CG, class RN, class CN, class Mat>
    void getElementMatrix(CG const& contextGeo, RN const& rowNode, CN const& colNode, Mat& elementMatrix)
49
    {
50
      static_assert(RN::isLeaf && CN::isLeaf,
51
        "Operator can be applied to Leaf-Nodes only.");
52

53
      auto const& quad = this->getQuadratureRule(contextGeo.type(), rowNode, colNode);
54
55
      auto const& rowLocalFE = rowNode.finiteElement();
      auto const& colLocalFE = colNode.finiteElement();
56
57
      std::size_t rowSize = rowLocalFE.size();
      std::size_t colSize = colLocalFE.size();
58

59
      using RangeFieldType = typename NodeQuadCache<CN>::Traits::RangeFieldType;
60

61
62
      NodeQuadCache<RN> rowCache(rowLocalFE.localBasis());
      NodeQuadCache<CN> colCache(colLocalFE.localBasis());
63

64
65
      auto const& shapeValuesCache = rowCache.evaluateFunctionAtQP(contextGeo, quad);
      auto const& shapeGradientsCache = colCache.evaluateJacobianAtQP(contextGeo, quad);
66
      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
67
        // Position of the current quadrature point in the reference element
68
        decltype(auto) local = contextGeo.local(quad[iq].position());
69
70

        // The transposed inverse Jacobian of the map from the reference element to the element
71
        const auto jacobian = contextGeo.geometry().jacobianInverseTransposed(local);
72
        auto&& jacobian_mat = Dune::MatVec::as_matrix(jacobian);
73
        assert(jacobian.M() == CG::dim);
74
75

        // The multiplicative factor in the integral transformation formula
76
        const auto factor = contextGeo.integrationElement(quad[iq].position()) * quad[iq].weight();
77
        const auto c = Super::coefficient(local);
78

79
80
        // the values of the shape functions on the reference element at the quadrature point
        auto const& shapeValues = shapeValuesCache[iq];
81
82

        // The gradients of the shape functions on the reference element
83
        auto const& shapeGradients = shapeGradientsCache[iq];
84
85

        // Compute the shape function gradients on the real element
86
87
        thread_local std::vector<RangeFieldType> colPartial;
        colPartial.resize(shapeGradients.size());
88
        for (std::size_t i = 0; i < colPartial.size(); ++i) {
89
90
91
          colPartial[i] = jacobian_mat[comp_][0] * shapeGradients[i][0][0];
          for (std::size_t j = 1; j < jacobian_mat.M(); ++j)
            colPartial[i] += jacobian_mat[comp_][j] * shapeGradients[i][0][j];
92
93
        }

94
        for (std::size_t j = 0; j < colSize; ++j) {
95
96
          const auto local_j = colNode.localIndex(j);
          const auto value = factor * (c * colPartial[j]);
97
          for (std::size_t i = 0; i < rowSize; ++i) {
98
99
            const auto local_i = rowNode.localIndex(i);
            elementMatrix[local_i][local_j] += value * shapeValues[i];
100
101
102
103
104
105
          }
        }
      }
    }

  private:
106
    std::size_t comp_;
107
108
  };

109
110
111
112
113
114
115
116

  /// Create a first-order term with derivative on trial-function
  template <class Expr, class... QuadratureArgs>
  auto fot(Expr&& expr, tag::partial_trial t, QuadratureArgs&&... args)
  {
    return makeOperator(tag::test_partialtrial{t.comp}, FWD(expr), FWD(args)...);
  }

117
118
  /** @} **/

119
} // end namespace AMDiS