AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
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
9#include <dune/common/typetree/nodeconcepts.hh>
10
11namespace AMDiS
12{
18 namespace tag
19 {
21 {
22 std::size_t comp;
23 };
24
26 {
27 std::size_t comp;
28 };
29 }
30
31
34 {
35 int comp_;
36
37 public:
39 : comp_(tag.comp)
40 {}
41
42 template <class CG, class RN, class CN, class Quad, class LocalFct, class Mat>
43 void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode,
44 Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
45 {
46 static_assert(static_size_v<typename LocalFct::Range> == 1,
47 "Expression must be of scalar type.");
48 static_assert(Dune::TypeTree::Concept::LeafTreeNode<RN> && Dune::TypeTree::Concept::LeafTreeNode<CN>,
49 "Operator can be applied to Leaf-Nodes only.");
50
51 std::size_t rowSize = rowNode.size();
52 std::size_t colSize = colNode.size();
53
54 LocalToGlobalBasisAdapter colLocalBasis(colNode, contextGeo.elementGeometry());
55 for (auto const& qp : quad) {
56 // Position of the current quadrature point in the reference element
57 decltype(auto) local = contextGeo.coordinateInElement(qp.position());
58
59 // The multiplicative factor in the integral transformation formula
60 auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
61 factor *= localFct(local);
62
63 // the values of the shape functions on the reference element at the quadrature point
64 auto const& shapeValues = rowNode.localBasisValuesAt(local);
65
66 // Compute the shape function gradients on the real element
67 auto const& colPartial = colLocalBasis.partialsAt(local, comp_);
68
69 for (std::size_t j = 0; j < colSize; ++j) {
70 const auto local_j = colNode.localIndex(j);
71 const auto value = factor * colPartial[j];
72 for (std::size_t i = 0; i < rowSize; ++i) {
73 const auto local_i = rowNode.localIndex(i);
74 elementMatrix[local_i][local_j] += value * shapeValues[i];
75 }
76 }
77 }
78 }
79 };
80
81 template <class LC>
82 struct GridFunctionOperatorRegistry<tag::test_partialtrial, LC>
83 {
84 static constexpr int degree = 1;
86 };
87
88
90 template <class Expr>
91 auto fot(Expr&& expr, tag::partial_trial t, int quadOrder = -1)
92 {
93 return makeOperator(tag::test_partialtrial{t.comp}, FWD(expr), quadOrder);
94 }
95
98} // end namespace AMDiS
first-order operator
Definition FirstOrderTestPartialTrial.hpp:34
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::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 FirstOrderTestPartialTrial.hpp:26
Definition FirstOrderTestPartialTrial.hpp:21