FirstOrderTestGradTrial.hpp 3.32 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
20
21
  namespace tag
  {
    struct test_gradtrial {};
  }


22
23
24
25
  /// first-order operator \f$ \langle\psi, \mathbf{b}\cdot\nabla\phi\rangle \f$
  template <class GridFct, class QuadCreator>
  class GridFunctionOperator<tag::test_gradtrial, GridFct, QuadCreator>
      : public GridFunctionOperatorBase<GridFct, QuadCreator>
26
  {
27
    using Super = GridFunctionOperatorBase<GridFct, QuadCreator>;
28

29
    static_assert( Category::Vector<typename GridFct::Range>, "Expression must be of vector type." );
30
31

  public:
32
33
    GridFunctionOperator(tag::test_gradtrial, GridFct const& expr, QuadCreator const& quadCreator)
      : Super(expr, quadCreator, 1)
34
35
36
37
38
39
40
41
42
43
44
45
    {}

    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>)
    {
46
47
      static_assert(RowNode::isLeaf && ColNode::isLeaf,
        "Operator can be applied to Leaf-Nodes only.");
48
49
50
51
52
53
54
55
56
57
58
59
60
61

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

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

      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
62
63
        double factor = context.integrationElement(quad[iq].position()) * quad[iq].weight();
        const auto b = Super::coefficient(local);
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91

        rowLocalFE.localBasis().evaluateFunction(local, rowShapeValues_);

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

        // Compute the shape function gradients on the real element
        std::vector<Dune::FieldVector<double,Context::dow> > colGradients(referenceGradients.size());

        for (std::size_t i = 0; i < colGradients.size(); ++i)
          jacobian.mv(referenceGradients[i][0], colGradients[i]);

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

    }

  private:
    std::vector<Dune::FieldVector<double,1>> rowShapeValues_;
  };

92
93
  /** @} **/

94
} // end namespace AMDiS