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