SecondOrderPartialTestPartialTrial.hpp 4.75 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 partialtest_partialtrial
    {
21
22
      std::size_t comp_test; // i
      std::size_t comp_trial; // j
23
24
25
26
    };
  }


27
  /// second-order operator \f$ \langle\partial_i\psi, c\,\partial_j\phi\rangle \f$
28
29
30
31
  template <class LocalContext, class GridFct>
  class GridFunctionOperator<tag::partialtest_partialtrial, LocalContext, GridFct>
      : public GridFunctionOperatorBase<GridFunctionOperator<tag::partialtest_partialtrial, LocalContext, GridFct>,
                                        LocalContext, GridFct>
32
  {
33
    using Self = GridFunctionOperator;
34
    using Super = GridFunctionOperatorBase<Self, LocalContext, GridFct>;
35

36
    static_assert( Size_v<typename GridFct::Range> == 1, "Expression must be of scalar type." );
37
38

  public:
39
40
    GridFunctionOperator(tag::partialtest_partialtrial tag, GridFct const& expr)
      : Super(expr, 2)
41
42
43
44
      , compTest_(tag.comp_test)
      , compTrial_(tag.comp_trial)
    {}

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

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

60
61
      using RowFieldType = typename NodeQuadCache<RowNode>::Traits::RangeFieldType;
      using ColFieldType = typename NodeQuadCache<ColNode>::Traits::RangeFieldType;
62

63
64
      NodeQuadCache<RowNode> rowCache(rowLocalFE.localBasis());
      NodeQuadCache<ColNode> colCache(colLocalFE.localBasis());
65

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

        // The transposed inverse Jacobian of the map from the reference element to the element
73
        const auto jacobian = context.geometry().jacobianInverseTransposed(local);
74
75

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

        // The gradients of the shape functions on the reference element
80
81
        auto const& rowShapeGradients = rowGradientsCache[iq];
        auto const& colShapeGradients = colGradientsCache[iq];
82
83

        // Compute the shape function gradients on the real element
84
85
        thread_local std::vector<RowFieldType> rowPartial;
        rowPartial.resize(rowShapeGradients.size());
86
        for (std::size_t i = 0; i < rowPartial.size(); ++i) {
87
          rowPartial[i] = jacobian[compTest_][0] * rowShapeGradients[i][0][0];
88
          for (std::size_t j = 1; j < jacobian.cols; ++j)
89
            rowPartial[i] += jacobian[compTest_][j] * rowShapeGradients[i][0][j];
90
91
        }

92
93
        thread_local std::vector<ColFieldType> colPartial;
        colPartial.resize(colShapeGradients.size());
94
        for (std::size_t i = 0; i < colPartial.size(); ++i) {
95
          colPartial[i] = jacobian[compTrial_][0] * colShapeGradients[i][0][0];
96
          for (std::size_t j = 1; j < jacobian.cols; ++j)
97
            colPartial[i] += jacobian[compTrial_][j] * colShapeGradients[i][0][j];
98
99
        }

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

    }

  private:
113
114
    std::size_t compTest_;
    std::size_t compTrial_;
115
116
  };

117
118
119
120
121
122
123
124

  /// Create a second-order term of partial derivatives
  template <class Expr, class... QuadratureArgs>
  auto sot_ij(Expr&& expr, std::size_t comp_test, std::size_t comp_trial, QuadratureArgs&&... args)
  {
    return makeOperator(tag::partialtest_partialtrial{comp_test, comp_trial}, FWD(expr), FWD(args)...);
  }

125
126
  /** @} **/

127
} // end namespace AMDiS