AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
FirstOrderGradTest.hpp
1#pragma once
2
3#include <type_traits>
4
5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/common/StaticSize.hpp>
7
8#include <dune/common/typetree/nodeconcepts.hh>
9
10namespace AMDiS
11{
17 namespace tag
18 {
19 struct gradtest {};
20 }
21
22
25 {
26 public:
28
29 template <class CG, class Node, class Quad, class LocalFct, class Vec>
30 void assemble(CG const& contextGeo, Node const& node, Quad const& quad,
31 LocalFct const& localFct, Vec& elementVector) const
32 {
33 static_assert(static_size_v<typename LocalFct::Range> == CG::dow, "Expression must be of vector type.");
34 static_assert(Dune::TypeTree::Concept::LeafTreeNode<Node>, "Node must be Leaf-Node.");
35
36 std::size_t feSize = node.size();
37
38 using RangeFieldType = typename Node::LocalBasis::Traits::RangeFieldType;
39 using WorldVector = FieldVector<RangeFieldType,CG::dow>;
40 std::vector<WorldVector> gradients;
41
42 for (auto const& qp : quad) {
43 // Position of the current quadrature point in the reference element
44 auto&& local = contextGeo.coordinateInElement(qp.position());
45
46 // The transposed inverse Jacobian of the map from the reference element to the element
47 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
48
49 // The multiplicative factor in the integral transformation formula
50 const auto factor = localFct(local);
51 const auto dx = contextGeo.integrationElement(qp.position()) * qp.weight();
52
53 // The gradients of the shape functions on the reference element
54 auto const& shapeGradients = node.localBasisJacobiansAt(local);
55
56 // Compute the shape function gradients on the real element
57 gradients.resize(shapeGradients.size());
58
59 for (std::size_t i = 0; i < gradients.size(); ++i)
60 jacobian.mv(shapeGradients[i][0], gradients[i]);
61
62 for (std::size_t i = 0; i < feSize; ++i) {
63 const auto local_i = node.localIndex(i);
64 elementVector[local_i] += dx * (factor * gradients[i]);
65 }
66 }
67 }
68 };
69
70
71 template <class LC>
72 struct GridFunctionOperatorRegistry<tag::gradtest, LC>
73 {
74 static constexpr int degree = 1;
76 };
77
80} // end namespace AMDiS
first-order operator
Definition FirstOrderGradTest.hpp:25
Registry to specify a tag for each implementation type.
Definition GridFunctionOperator.hpp:216
Definition FirstOrderGradTest.hpp:19