AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
FirstOrderTestvecGradTrial.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 {
20 }
21
22
25 {
26 public:
28
29 template <class CG, class RN, class CN, class Quad, class LocalFct, class Mat>
30 void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode,
31 Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
32 {
33 static_assert(static_size_v<typename LocalFct::Range> == 1,
34 "Expression must be of scalar type.");
35 static_assert(Dune::TypeTree::Concept::UniformInnerTreeNode<RN> && Dune::TypeTree::Concept::LeafTreeNode<CN>,
36 "col-node must be Leaf-Node and row-node must be a Power-Node.");
37
38 assert( rowNode.degree() == CG::dow );
39
40 std::size_t rowSize = rowNode.child(0).size();
41 std::size_t colSize = colNode.size();
42
43 using RangeFieldType = typename CN::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 decltype(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.child(0).localBasisValuesAt(local);
59
60 // The gradients of the shape functions on the reference element
61 auto const& shapeGradients = colNode.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 value = factor * shapeValues[i];
71 for (std::size_t j = 0; j < colSize; ++j) {
72 const auto local_j = colNode.localIndex(j);
73 for (std::size_t k = 0; k < rowNode.degree(); ++k) {
74 const auto local_ki = rowNode.child(k).localIndex(i);
75 elementMatrix[local_ki][local_j] += value * colGradients[j][k];
76 }
77 }
78 }
79 }
80
81 }
82 };
83
84 template <class LC>
85 struct GridFunctionOperatorRegistry<tag::testvec_gradtrial, LC>
86 {
87 static constexpr int degree = 1;
89 };
90
93} // end namespace AMDiS
first-order operator
Definition FirstOrderTestvecGradTrial.hpp:25
Registry to specify a tag for each implementation type.
Definition GridFunctionOperator.hpp:216
Definition FirstOrderTestvecGradTrial.hpp:19