FirstOrderPartialTest.hpp 3.16 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
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
51
52
53
54
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
89
90
91
92
93
94
95
#pragma once

#include <type_traits>

#include <amdis/assembler/FirstOrderTestPartialTrial.hpp>

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

  namespace tag
  {
    struct partialtest
    {
      std::size_t comp;
    };
  }


  /// first-order operator \f$ \langle\partial_i\psi, c\rangle \f$

  template <class LocalContext, class GridFct>
  class GridFunctionOperator<tag::partialtest, LocalContext, GridFct>
      : public GridFunctionOperatorBase<GridFunctionOperator<tag::partialtest, LocalContext, GridFct>,
                                        LocalContext, GridFct>
  {
    using Self = GridFunctionOperator;
    using Super = GridFunctionOperatorBase<Self, LocalContext, GridFct>;

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

  public:
    GridFunctionOperator(tag::partialtest tag, GridFct const& expr)
      : Super(expr, 1)
      , comp_(tag.comp)
    {}

    template <class Context, class Node, class ElementVector>
    void getElementVector(Context const& context,
                          Node const& node,
                          ElementVector& elementVector)
    {
      static_assert(Node::isLeaf, "Operator can be applied to Leaf-Nodes only.");

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

      using LocalBasisType = typename std::decay_t<decltype(localFE)>::Traits::LocalBasisType;
      using RangeFieldType = typename LocalBasisType::Traits::RangeFieldType;

      LocalBasisCache<LocalBasisType> cache(localFE.localBasis());

      auto const& shapeGradientsCache = cache.evaluateJacobianAtQp(context, quad);
      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);
        assert(jacobian.M() == Context::dim);

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

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

        // Compute the shape function gradients on the real element
        thread_local std::vector<RangeFieldType> partial;
        partial.resize(shapeGradients.size());
        for (std::size_t i = 0; i < partial.size(); ++i) {
          partial[i] = jacobian[comp_][0] * shapeGradients[i][0][0];
          for (std::size_t j = 1; j < jacobian.M(); ++j)
            partial[i] += jacobian[comp_][j] * shapeGradients[i][0][j];
        }

        for (std::size_t j = 0; j < feSize; ++j) {
          const auto local_j = node.localIndex(j);
          elementVector[local_j] += dx * (c * partial[j]);
        }
      }
    }

  private:
    std::size_t comp_;
  };

  /** @} **/

} // end namespace AMDiS