AMDiS 2.10
The Adaptive Multi-Dimensional Simulation Toolbox
SecondOrderPartialTestPartialTrial.hpp
1#pragma once
2
3#include <type_traits>
4
5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/common/StaticSize.hpp>
7
8namespace AMDiS
9{
15 namespace tag
16 {
18 {
19 std::size_t comp_test; // i
20 std::size_t comp_trial; // j
21 };
22 }
23
24
27 {
28 int compTest_;
29 int compTrial_;
30
31 public:
33 : compTest_(tag.comp_test)
34 , compTrial_(tag.comp_trial)
35 {}
36
37 template <class CG, class RN, class CN, class Quad, class LocalFct, class Mat>
38 void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode,
39 Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
40 {
41 static_assert(static_size_v<typename LocalFct::Range> == 1,
42 "Expression must be of scalar type." );
43 static_assert(RN::isLeaf && CN::isLeaf,
44 "Operator can be applied to Leaf-Nodes only.");
45
46 LocalToGlobalBasisAdapter rowLocalBasis(rowNode, contextGeo.elementGeometry());
47 LocalToGlobalBasisAdapter colLocalBasis(colNode, contextGeo.elementGeometry());
48
49 for (auto const& qp : quad) {
50 // Position of the current quadrature point in the reference element
51 decltype(auto) local = contextGeo.coordinateInElement(qp.position());
52
53 // The multiplicative factor in the integral transformation formula
54 auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
55 factor *= localFct(local);
56
57 auto const& rowPartial = rowLocalBasis.partialsAt(local, compTest_);
58 auto const& colPartial = colLocalBasis.partialsAt(local, compTrial_);
59
60 for (std::size_t j = 0; j < colLocalBasis.size(); ++j) {
61 const auto local_j = colNode.localIndex(j);
62 const auto value = factor * colPartial[j];
63 for (std::size_t i = 0; i < rowLocalBasis.size(); ++i) {
64 const auto local_i = rowNode.localIndex(i);
65 elementMatrix[local_i][local_j] += value * rowPartial[i];
66 }
67 }
68 }
69
70 }
71 };
72
73 template <class LC>
74 struct GridFunctionOperatorRegistry<tag::partialtest_partialtrial, LC>
75 {
76 static constexpr int degree = 2;
78 };
79
80
82 template <class Expr>
83 auto sot_ij(Expr&& expr, std::size_t comp_test, std::size_t comp_trial, int quadOrder = -1)
84 {
85 return makeOperator(tag::partialtest_partialtrial{comp_test, comp_trial}, FWD(expr), quadOrder);
86 }
87
90} // end namespace AMDiS
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
std::size_t size() const
Return the number of local basis functions.
Definition: LocalToGlobalAdapter.hpp:96
second-order operator
Definition: SecondOrderPartialTestPartialTrial.hpp:27
auto sot_ij(Expr &&expr, std::size_t comp_test, std::size_t comp_trial, int quadOrder=-1)
Create a second-order term of partial derivatives.
Definition: SecondOrderPartialTestPartialTrial.hpp:83
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: SecondOrderPartialTestPartialTrial.hpp:18