SecondOrderPartialTestPartialTrial.hpp 4.46 KB
Newer Older
1
2
3
4
#pragma once

#include <type_traits>

5
6
7
#include <amdis/GridFunctionOperator.hpp>
#include <amdis/Output.hpp>
#include <amdis/common/ValueCategory.hpp>
8
9
10

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

16
17
18
19
  namespace tag
  {
    struct partialtest_partialtrial
    {
20
21
      std::size_t comp_test; // i
      std::size_t comp_trial; // j
22
23
24
25
    };
  }


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

35
    static_assert( Category::Scalar<typename GridFct::Range>, "Expression must be of scalar type." );
36
37

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

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

      auto const& rowLocalFE = rowNode.finiteElement();
      auto const& colLocalFE = colNode.finiteElement();

56
57
      using RowLocalBasisType = typename std::decay_t<decltype(rowLocalFE)>::Traits::LocalBasisType;
      using ColLocalBasisType = typename std::decay_t<decltype(colLocalFE)>::Traits::LocalBasisType;
58

59
60
61
62
63
64
65
      using RowFieldType = typename RowLocalBasisType::Traits::RangeFieldType;
      using ColFieldType = typename ColLocalBasisType::Traits::RangeFieldType;
      using RowJacobianType = typename RowLocalBasisType::Traits::JacobianType;
      using ColJacobianType = typename ColLocalBasisType::Traits::JacobianType;

      auto const& quad = this->getQuadratureRule(context.type(), rowNode, colNode);
      for (auto const& qp : quad) {
66
        // Position of the current quadrature point in the reference element
67
        decltype(auto) local = context.local(qp.position());
68
69

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

        // The multiplicative factor in the integral transformation formula
73
74
        const auto factor = context.integrationElement(qp.position()) * qp.weight();
        const auto c = Super::coefficient(local);
75
76

        // The gradients of the shape functions on the reference element
77
78
79
80
        thread_local std::vector<RowJacobianType> rowShapeGradients;
        rowLocalFE.localBasis().evaluateJacobian(local, rowShapeGradients);
        thread_local std::vector<ColJacobianType> colShapeGradients;
        colLocalFE.localBasis().evaluateJacobian(local, colShapeGradients);
81
82

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

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

        for (std::size_t j = 0; j < colLocalFE.size(); ++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 < rowLocalFE.size(); ++i) {
101
            const auto local_i = rowNode.localIndex(i);
102
103
104
105
106
107
108
109
            elementMatrix[local_i][local_j] += value * rowPartial[i];
          }
        }
      }

    }

  private:
110
111
    std::size_t compTest_;
    std::size_t compTrial_;
112
113
  };

114
115
  /** @} **/

116
} // end namespace AMDiS