AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
FirstOrderTestGradTrial.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 test_gradtrial {};
20 struct grad_trial {};
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> == CG::dow, "Expression must be of vector type.");
35 static_assert(Dune::TypeTree::Concept::LeafTreeNode<RN> && Dune::TypeTree::Concept::LeafTreeNode<CN>,
36 "Operator can be applied to Leaf-Nodes only.");
37
38 std::size_t rowSize = rowNode.size();
39 std::size_t colSize = colNode.size();
40
41 using RangeFieldType = typename CN::LocalBasis::Traits::RangeFieldType;
42 using WorldVector = FieldVector<RangeFieldType,CG::dow>;
43 std::vector<WorldVector> colGradients;
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 = contextGeo.integrationElement(qp.position()) * qp.weight();
54 const auto b = localFct(local);
55
56 // the values of the shape functions on the reference element at the quadrature point
57 auto const& shapeValues = rowNode.localBasisValuesAt(local);
58
59 // The gradients of the shape functions on the reference element
60 auto const& shapeGradients = colNode.localBasisJacobiansAt(local);
61
62 // Compute the shape function gradients on the real element
63 colGradients.resize(shapeGradients.size());
64
65 for (std::size_t i = 0; i < colGradients.size(); ++i)
66 jacobian.mv(shapeGradients[i][0], colGradients[i]);
67
68 for (std::size_t j = 0; j < colSize; ++j) {
69 const auto local_j = colNode.localIndex(j);
70 const auto value = factor * (b * colGradients[j]);
71 for (std::size_t i = 0; i < rowSize; ++i) {
72 const auto local_i = rowNode.localIndex(i);
73 elementMatrix[local_i][local_j] += value * shapeValues[i];
74 }
75 }
76 }
77
78 }
79 };
80
81 template <class LC>
82 struct GridFunctionOperatorRegistry<tag::test_gradtrial, LC>
83 {
84 static constexpr int degree = 1;
86 };
87
88
90 template <class Expr>
91 auto fot(Expr&& expr, tag::grad_trial, int quadOrder = -1)
92 {
93 return makeOperator(tag::test_gradtrial{}, FWD(expr), quadOrder);
94 }
95
98} // end namespace AMDiS
first-order operator
Definition FirstOrderTestGradTrial.hpp:26
auto fot(Expr &&expr, tag::grad_test, int quadOrder=-1)
Create a first-order term with derivative on trial-function.
Definition FirstOrderGradTestTrial.hpp:33
auto makeOperator(Tag const &tag, Expr &&expr, int gridFctDeg=-1)
Definition GridFunctionOperator.hpp:235
Registry to specify a tag for each implementation type.
Definition GridFunctionOperator.hpp:216
Definition FirstOrderTestGradTrial.hpp:20
Definition FirstOrderTestGradTrial.hpp:19