AMDiS 2.10
The Adaptive Multi-Dimensional Simulation Toolbox
FirstOrderTestPartialTrial.hpp
1#pragma once
2
3#include <type_traits>
4
5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/common/StaticSize.hpp>
7#include <amdis/utility/LocalToGlobalAdapter.hpp>
8
9namespace AMDiS
10{
16 namespace tag
17 {
19 {
20 std::size_t comp;
21 };
22
24 {
25 std::size_t comp;
26 };
27 }
28
29
32 {
33 int comp_;
34
35 public:
37 : comp_(tag.comp)
38 {}
39
40 template <class CG, class RN, class CN, class Quad, class LocalFct, class Mat>
41 void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode,
42 Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
43 {
44 static_assert(static_size_v<typename LocalFct::Range> == 1,
45 "Expression must be of scalar type.");
46 static_assert(RN::isLeaf && CN::isLeaf,
47 "Operator can be applied to Leaf-Nodes only.");
48
49 std::size_t rowSize = rowNode.size();
50 std::size_t colSize = colNode.size();
51
52 LocalToGlobalBasisAdapter colLocalBasis(colNode, contextGeo.elementGeometry());
53 for (auto const& qp : quad) {
54 // Position of the current quadrature point in the reference element
55 decltype(auto) local = contextGeo.coordinateInElement(qp.position());
56
57 // The multiplicative factor in the integral transformation formula
58 auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
59 factor *= localFct(local);
60
61 // the values of the shape functions on the reference element at the quadrature point
62 auto const& shapeValues = rowNode.localBasisValuesAt(local);
63
64 // Compute the shape function gradients on the real element
65 auto const& colPartial = colLocalBasis.partialsAt(local, comp_);
66
67 for (std::size_t j = 0; j < colSize; ++j) {
68 const auto local_j = colNode.localIndex(j);
69 const auto value = factor * colPartial[j];
70 for (std::size_t i = 0; i < rowSize; ++i) {
71 const auto local_i = rowNode.localIndex(i);
72 elementMatrix[local_i][local_j] += value * shapeValues[i];
73 }
74 }
75 }
76 }
77 };
78
79 template <class LC>
80 struct GridFunctionOperatorRegistry<tag::test_partialtrial, LC>
81 {
82 static constexpr int degree = 1;
84 };
85
86
88 template <class Expr>
89 auto fot(Expr&& expr, tag::partial_trial t, int quadOrder = -1)
90 {
91 return makeOperator(tag::test_partialtrial{t.comp}, FWD(expr), quadOrder);
92 }
93
96} // end namespace AMDiS
first-order operator
Definition: FirstOrderTestPartialTrial.hpp:32
Convert a simple (scalar) local basis into a global basis.
Definition: LocalToGlobalAdapter.hpp:65
auto const & partialsAt(typename Traits::DomainLocal const &x, std::size_t comp) const
Definition: LocalToGlobalAdapter.hpp:175
auto fot(Expr &&expr, tag::partial_trial t, int quadOrder=-1)
Create a first-order term with derivative on trial-function.
Definition: FirstOrderTestPartialTrial.hpp:89
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: FirstOrderTestPartialTrial.hpp:24
Definition: FirstOrderTestPartialTrial.hpp:19