AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
ZeroOrderTestTrial.hpp
1#pragma once
2
3#include <type_traits>
4
5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/common/StaticSize.hpp>
7#include <amdis/typetree/FiniteElementType.hpp>
8
9#include <dune/common/typetree/nodeconcepts.hh>
10
11namespace AMDiS
12{
18 namespace tag
19 {
20 struct test_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> == 1,
35 "Expression must be of scalar type." );
36
37 const bool sameFE = std::is_same_v<FiniteElementType_t<RN>, FiniteElementType_t<CN>>;
38 const bool sameNode = rowNode.treeIndex() == colNode.treeIndex();
39
40 if (sameFE && sameNode)
41 getElementMatrixOptimized(contextGeo, quad, rowNode, colNode, localFct, elementMatrix);
42 else
43 getElementMatrixStandard(contextGeo, quad, rowNode, colNode, localFct, elementMatrix);
44 }
45
46 template <class CG, class Node, class Quad, class LocalFct, class Vec>
47 void assemble(CG const& contextGeo, Node const& node,
48 Quad const& quad, LocalFct const& localFct, Vec& elementVector) const
49 {
50 static_assert(static_size_v<typename LocalFct::Range> == 1,
51 "Expression must be of scalar type." );
52 static_assert(Dune::TypeTree::Concept::LeafTreeNode<Node>,
53 "Operator can be applied to Leaf-Nodes only");
54
55 std::size_t size = node.size();
56
57 for (auto const& qp : quad) {
58 // Position of the current quadrature point in the reference element
59 auto&& local = contextGeo.coordinateInElement(qp.position());
60
61 // The multiplicative factor in the integral transformation formula
62 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
63
64 auto const& shapeValues = node.localBasisValuesAt(local);
65 for (std::size_t i = 0; i < size; ++i) {
66 const auto local_i = node.localIndex(i);
67 elementVector[local_i] += factor * shapeValues[i];
68 }
69 }
70 }
71
72
73 protected:
74 template <class CG, class QR, class RN, class CN, class LocalFct, class Mat>
75 void getElementMatrixStandard(CG const& contextGeo, QR const& quad,
76 RN const& rowNode, CN const& colNode,
77 LocalFct const& localFct, Mat& elementMatrix) const
78 {
79 static_assert(Dune::TypeTree::Concept::LeafTreeNode<RN> && Dune::TypeTree::Concept::LeafTreeNode<CN>,
80 "Operator can be applied to Leaf-Nodes only.");
81
82 std::size_t rowSize = rowNode.size();
83 std::size_t colSize = colNode.size();
84
85 for (auto const& qp : quad) {
86 // Position of the current quadrature point in the reference element
87 auto&& local = contextGeo.coordinateInElement(qp.position());
88
89 // The multiplicative factor in the integral transformation formula
90 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
91
92 auto const& rowShapeValues = rowNode.localBasisValuesAt(local);
93 auto const& colShapeValues = colNode.localBasisValuesAt(local);
94
95 for (std::size_t i = 0; i < rowSize; ++i) {
96 const auto local_i = rowNode.localIndex(i);
97 const auto value = factor * rowShapeValues[i];
98
99 for (std::size_t j = 0; j < colSize; ++j) {
100 const auto local_j = colNode.localIndex(j);
101 elementMatrix[local_i][local_j] += value * colShapeValues[j];
102 }
103 }
104 }
105 }
106
107
108 template <class CG, class QR, class RN, class CN, class LocalFct, class Mat>
109 void getElementMatrixOptimized(CG const& contextGeo, QR const& quad,
110 RN const& node, CN const& /*colNode*/,
111 LocalFct const& localFct, Mat& elementMatrix) const
112 {
113 static_assert(Dune::TypeTree::Concept::LeafTreeNode<RN> && Dune::TypeTree::Concept::LeafTreeNode<CN>,
114 "Operator can be applied to Leaf-Nodes only.");
115
116 std::size_t size = node.size();
117
118 for (auto const& qp : quad) {
119 // Position of the current quadrature point in the reference element
120 auto&& local = contextGeo.coordinateInElement(qp.position());
121
122 // The multiplicative factor in the integral transformation formula
123 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
124
125 auto const& shapeValues = node.localBasisValuesAt(local);
126
127 for (std::size_t i = 0; i < size; ++i) {
128 const auto local_i = node.localIndex(i);
129
130 const auto value = factor * shapeValues[i];
131 elementMatrix[local_i][local_i] += value * shapeValues[i];
132
133 for (std::size_t j = i+1; j < size; ++j) {
134 const auto local_j = node.localIndex(j);
135
136 elementMatrix[local_i][local_j] += value * shapeValues[j];
137 elementMatrix[local_j][local_i] += value * shapeValues[j];
138 }
139 }
140 }
141 }
142 };
143
144 template <class LC>
145 struct GridFunctionOperatorRegistry<tag::test_trial, LC>
146 {
147 static constexpr int degree = 0;
148 using type = ZeroOrderTestTrial;
149 };
150
151
153 template <class Expr>
154 auto zot(Expr&& expr, int quadOrder = -1)
155 {
156 return makeOperator(tag::test_trial{}, FWD(expr), quadOrder);
157 }
158
161} // end namespace AMDiS
zero-order operator
Definition ZeroOrderTestTrial.hpp:26
auto makeOperator(Tag const &tag, Expr &&expr, int gridFctDeg=-1)
Definition GridFunctionOperator.hpp:235
auto zot(Expr &&expr, int quadOrder=-1)
Create a zero-order term.
Definition ZeroOrderTestTrial.hpp:154
Registry to specify a tag for each implementation type.
Definition GridFunctionOperator.hpp:216
Definition ZeroOrderTestTrial.hpp:20