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