SecondOrderPartialTestPartialTrial.hpp 4.13 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
27
28
29
  /// second-order operator \f$ \langle\partial_i\psi, c\,\partial_j\phi\rangle \f$
  template <class GridFct, class QuadCreator>
  class GridFunctionOperator<tag::partialtest_partialtrial, GridFct, QuadCreator>
      : public GridFunctionOperatorBase<GridFct, QuadCreator>
30
  {
31
    using Super = GridFunctionOperatorBase<GridFct, QuadCreator>;
32

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

  public:
36
37
    GridFunctionOperator(tag::partialtest_partialtrial tag, GridFct const& expr, QuadCreator const& quadCreator)
      : Super(expr, quadCreator, 2)
38
39
40
41
42
43
44
45
46
47
48
49
50
51
      , compTest_(tag.comp_test)
      , compTrial_(tag.comp_trial)
    {}

    template <class Context, class QuadratureRule,
              class ElementMatrix, class RowNode, class ColNode, bool sameFE, bool sameNode>
    void calculateElementMatrix(Context const& context,
                                QuadratureRule const& quad,
                                ElementMatrix& elementMatrix,
                                RowNode const& rowNode,
                                ColNode const& colNode,
                                std::integral_constant<bool, sameFE>,
                                std::integral_constant<bool, sameNode>)
    {
52
53
      static_assert(RowNode::isLeaf && ColNode::isLeaf,
        "Operator can be applied to Leaf-Nodes only.");
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68

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

      std::vector<Dune::FieldMatrix<double,1,Context::dim> > rowReferenceGradients;
      std::vector<Dune::FieldMatrix<double,1,Context::dim> > colReferenceGradients;

      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
        // Position of the current quadrature point in the reference element
        decltype(auto) local = context.position(quad[iq].position());

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

        // The multiplicative factor in the integral transformation formula
69
70
        double factor = context.integrationElement(quad[iq].position()) * quad[iq].weight();
        double c = Super::coefficient(local);
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105

        // The gradients of the shape functions on the reference element
        rowLocalFE.localBasis().evaluateJacobian(local, rowReferenceGradients);
        colLocalFE.localBasis().evaluateJacobian(local, colReferenceGradients);

        // Compute the shape function gradients on the real element
        std::vector<double> rowPartial(rowReferenceGradients.size());
        std::vector<double> colPartial(colReferenceGradients.size());

        // transform gradients to global element
        for (std::size_t i = 0; i < rowPartial.size(); ++i) {
          rowPartial[i] = jacobian[compTest_][0] * rowReferenceGradients[i][0][0];
          for (std::size_t j = 1; j < jacobian.cols(); ++j)
            rowPartial[i] += jacobian[compTest_][j] * rowReferenceGradients[i][0][j];
        }

        for (std::size_t i = 0; i < colPartial.size(); ++i) {
          colPartial[i] = jacobian[compTrial_][0] * colReferenceGradients[i][0][0];
          for (std::size_t j = 1; j < jacobian.cols(); ++j)
            colPartial[i] += jacobian[compTrial_][j] * colReferenceGradients[i][0][j];
        }

        for (std::size_t j = 0; j < colLocalFE.size(); ++j) {
          const int local_j = colNode.localIndex(j);
          double value = factor * (c * colPartial[j]);
          for (std::size_t i = 0; i < rowLocalFE.size(); ++i) {
            const int local_i = rowNode.localIndex(i);
            elementMatrix[local_i][local_j] += value * rowPartial[i];
          }
        }
      }

    }

  private:
106
107
    std::size_t compTest_;
    std::size_t compTrial_;
108
109
  };

110
111
  /** @} **/

112
} // end namespace AMDiS