FirstOrderTestPartialTrial.hpp 4.25 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
35
  template <class LocalContext, class GridFct>
  class GridFunctionOperator<tag::test_partialtrial, LocalContext, GridFct>
      : public GridFunctionOperatorBase<GridFunctionOperator<tag::test_partialtrial, LocalContext, GridFct>,
                                        LocalContext, GridFct>
36
  {
37
    using Self = GridFunctionOperator;
38
    using Super = GridFunctionOperatorBase<Self, LocalContext, GridFct>;
39

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

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

48
49
50
51
52
    template <class Context, class RowNode, class ColNode, class ElementMatrix>
    void getElementMatrix(Context const& context,
                          RowNode const& rowNode,
                          ColNode const& colNode,
                          ElementMatrix& elementMatrix)
53
    {
54
55
      static_assert(RowNode::isLeaf && ColNode::isLeaf,
        "Operator can be applied to Leaf-Nodes only.");
56

57
      auto const& quad = this->getQuadratureRule(context.type(), rowNode, colNode);
58
59
      auto const& rowLocalFE = rowNode.finiteElement();
      auto const& colLocalFE = colNode.finiteElement();
60
61
      std::size_t rowSize = rowLocalFE.size();
      std::size_t colSize = colLocalFE.size();
62

63
      using RangeFieldType = typename NodeQuadCache<ColNode>::Traits::RangeFieldType;
64

65
66
      NodeQuadCache<RowNode> rowCache(rowLocalFE.localBasis());
      NodeQuadCache<ColNode> colCache(colLocalFE.localBasis());
67

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

        // The transposed inverse Jacobian of the map from the reference element to the element
75
        const auto jacobian = context.geometry().jacobianInverseTransposed(local);
76
        assert(jacobian.M() == Context::dim);
77
78

        // The multiplicative factor in the integral transformation formula
79
        const auto factor = context.integrationElement(quad[iq].position()) * quad[iq].weight();
80
        const auto c = Super::coefficient(local);
81

82
83
        // the values of the shape functions on the reference element at the quadrature point
        auto const& shapeValues = shapeValuesCache[iq];
84
85

        // The gradients of the shape functions on the reference element
86
        auto const& shapeGradients = shapeGradientsCache[iq];
87
88

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

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

  private:
109
    std::size_t comp_;
110
111
  };

112
113
114
115
116
117
118
119

  /// 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)...);
  }

120
121
  /** @} **/

122
} // end namespace AMDiS