FirstOrderGradTest.hpp 2.97 KB
Newer Older
1
2
3
4
5
6
#pragma once

#include <type_traits>

#include <amdis/GridFunctionOperator.hpp>
#include <amdis/Output.hpp>
7
8
#include <amdis/common/StaticSize.hpp>
#include <amdis/utility/LocalBasisCache.hpp>
9
10
11
12
13
14
15
16
17
18
19
20
21
22

namespace AMDiS
{
  /**
   * \addtogroup operators
   * @{
   **/

  namespace tag
  {
    struct gradtest {};
  }


23
  /// first-order operator \f$ \langle\nabla\psi, b\rangle \f$
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
  template <class LocalContext, class GridFct>
  class GridFunctionOperator<tag::gradtest, LocalContext, GridFct>
      : public GridFunctionOperatorBase<GridFunctionOperator<tag::gradtest, LocalContext, GridFct>,
                                        LocalContext, GridFct>
  {
    static const int dow = ContextGeometry<LocalContext>::dow;
    using Self = GridFunctionOperator;
    using Super = GridFunctionOperatorBase<Self, LocalContext, GridFct>;

    static_assert( Size_v<typename GridFct::Range> == dow, "Expression must be of vector type." );

  public:
    GridFunctionOperator(tag::gradtest, GridFct const& expr)
      : Super(expr, 1)
    {}

    template <class Context, class Node, class ElementVector>
    void getElementVector(Context const& context,
                          Node const& node,
                          ElementVector& elementVector)
    {
      static_assert(Node::isLeaf, "Node must be Leaf-Node.");

      auto const& quad = this->getQuadratureRule(context.type(), node);
      auto const& localFE = node.finiteElement();
      std::size_t feSize = localFE.size();

51
52
      using RangeFieldType = typename NodeQuadCache<Node>::Traits::RangeFieldType;
      NodeQuadCache<Node> cache(localFE.localBasis());
53

54
      auto const& shapeGradientsCache = cache.evaluateJacobianAtQP(context, quad);
55
56
57
58
59
60
61
62
63
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
      for (std::size_t iq = 0; iq < quad.size(); ++iq) {
        // Position of the current quadrature point in the reference element
        decltype(auto) local = context.local(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
        const auto factor = Super::coefficient(local);
        const auto dx = context.integrationElement(quad[iq].position()) * quad[iq].weight();

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

        // Compute the shape function gradients on the real element
        using WorldVector = FieldVector<RangeFieldType,Context::dow>;
        thread_local std::vector<WorldVector> gradients;
        gradients.resize(shapeGradients.size());

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

        for (std::size_t i = 0; i < feSize; ++i) {
          const auto local_i = node.localIndex(i);
          elementVector[local_i] += dx * (factor * gradients[i]);
        }
      }

    }
  };

  /** @} **/

} // end namespace AMDiS