AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
ZeroOrderTestvecTrialvec.hpp
1#pragma once
2
3#include <type_traits>
4
5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/common/StaticSize.hpp>
7#include <amdis/common/ValueCategory.hpp>
8#include <amdis/typetree/FiniteElementType.hpp>
9
10#include <dune/common/typetree/nodeconcepts.hh>
11
12namespace AMDiS
13{
19 namespace tag
20 {
22 }
23
24
27 {
28 public:
30
31 template <class CG, class RN, class CN, class Quad, class LocalFct, class Mat>
32 void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode,
33 Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
34 {
35 using expr_value_type = typename LocalFct::Range;
36 static_assert(static_size_v<expr_value_type> == 1 ||
37 (static_num_rows_v<expr_value_type> == CG::dow &&
38 static_num_cols_v<expr_value_type> == CG::dow),
39 "Expression must be of scalar or matrix type." );
40 static_assert(Dune::TypeTree::Concept::UniformInnerTreeNode<RN> && Dune::TypeTree::Concept::UniformInnerTreeNode<CN>,
41 "Operator can be applied to Power-Nodes only.");
42 assert(rowNode.degree() == colNode.degree());
43
44 const bool sameFE = std::is_same_v<FiniteElementType_t<RN>, FiniteElementType_t<CN>>;
45 const bool sameNode = rowNode.treeIndex() == colNode.treeIndex();
46
47 using Category = ValueCategory_t<expr_value_type>;
48
49 if (sameFE && sameNode && std::is_same_v<Category,tag::scalar>)
50 getElementMatrixOptimized(contextGeo, quad, rowNode, colNode, localFct, elementMatrix, Category{});
51 else
52 getElementMatrixStandard(contextGeo, quad, rowNode, colNode, localFct, elementMatrix, Category{});
53 }
54
55 protected: // specialization for different value-categories
56
57 template <class CG, class QR, class RN, class CN, class LocalFct, class Mat>
58 void getElementMatrixStandard(CG const& contextGeo, QR const& quad,
59 RN const& rowNode, CN const& colNode,
60 LocalFct const& localFct, Mat& elementMatrix, tag::scalar) const
61 {
62 std::size_t rowSize = rowNode.child(0).size();
63 std::size_t colSize = colNode.child(0).size();
64
65 for (auto const& qp : quad) {
66 // Position of the current quadrature point in the reference element
67 auto&& local = contextGeo.coordinateInElement(qp.position());
68
69 // The multiplicative factor in the integral transformation formula
70 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
71
72 auto const& rowShapeValues = rowNode.child(0).localBasisValuesAt(local);
73 auto const& colShapeValues = colNode.child(0).localBasisValuesAt(local);
74
75 for (std::size_t i = 0; i < rowSize; ++i) {
76 const auto value = factor * rowShapeValues[i];
77
78 for (std::size_t j = 0; j < colSize; ++j) {
79 const auto value0 = value * colShapeValues[j];
80
81 for (std::size_t k = 0; k < rowNode.degree(); ++k) {
82 const auto local_ki = rowNode.child(k).localIndex(i);
83 const auto local_kj = colNode.child(k).localIndex(j);
84
85 elementMatrix[local_ki][local_kj] += value0;
86 }
87 }
88 }
89 }
90 }
91
92
93 template <class CG, class QR, class RN, class CN, class LocalFct, class Mat>
94 void getElementMatrixOptimized(CG const& contextGeo, QR const& quad,
95 RN const& node, CN const& /*colNode*/,
96 LocalFct const& localFct, Mat& elementMatrix, tag::scalar) const
97 {
98 std::size_t size = node.child(0).size();
99
100 for (auto const& qp : quad) {
101 // Position of the current quadrature point in the reference element
102 auto&& local = contextGeo.coordinateInElement(qp.position());
103
104 // The multiplicative factor in the integral transformation formula
105 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
106
107 auto const& shapeValues = node.child(0).localBasisValuesAt(local);
108
109 for (std::size_t i = 0; i < size; ++i) {
110 const auto value = factor * shapeValues[i];
111
112 for (std::size_t k = 0; k < node.degree(); ++k) {
113 const auto local_ki = node.child(k).localIndex(i);
114 elementMatrix[local_ki][local_ki] += value * shapeValues[i];
115 }
116
117 for (std::size_t j = i+1; j < size; ++j) {
118 const auto value0 = value * shapeValues[j];
119
120 for (std::size_t k = 0; k < node.degree(); ++k) {
121 const auto local_ki = node.child(k).localIndex(i);
122 const auto local_kj = node.child(k).localIndex(j);
123
124 elementMatrix[local_ki][local_kj] += value0;
125 elementMatrix[local_kj][local_ki] += value0;
126 }
127 }
128 }
129 }
130 }
131
132 template <class CG, class QR, class RN, class CN, class LocalFct, class Mat>
133 void getElementMatrixStandard(CG const& contextGeo, QR const& quad,
134 RN const& rowNode, CN const& colNode,
135 LocalFct const& localFct, Mat& elementMatrix, tag::matrix) const
136 {
137 using expr_value_type [[maybe_unused]] = typename LocalFct::Range;
138 assert(static_num_rows_v<expr_value_type> == rowNode.degree() &&
139 static_num_cols_v<expr_value_type> == rowNode.degree());
140
141 std::size_t rowSize = rowNode.child(0).size();
142 std::size_t colSize = colNode.child(0).size();
143
144 for (auto const& qp : quad) {
145 // Position of the current quadrature point in the reference element
146 auto&& local = contextGeo.coordinateInElement(qp.position());
147
148 // The multiplicative factor in the integral transformation formula
149 const auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
150 const auto exprValue = localFct(local);
151
152 auto const& rowShapeValues = rowNode.child(0).localBasisValuesAt(local);
153 auto const& colShapeValues = colNode.child(0).localBasisValuesAt(local);
154
155 for (std::size_t i = 0; i < rowSize; ++i) {
156 const auto value0 = factor * rowShapeValues[i];
157
158 for (std::size_t j = 0; j < colSize; ++j) {
159 const auto value = value0 * colShapeValues[j];
160 const auto mat = exprValue * value;
161
162 for (std::size_t k0 = 0; k0 < rowNode.degree(); ++k0) {
163 const auto local_ki = rowNode.child(k0).localIndex(i);
164 for (std::size_t k1 = 0; k1 < rowNode.degree(); ++k1) {
165 const auto local_kj = colNode.child(k1).localIndex(j);
166
167 elementMatrix[local_ki][local_kj] += mat[k0][k1];
168 }
169 }
170 }
171 }
172 }
173 }
174
175 template <class CG, class QR, class RN, class CN, class LocalFct, class Mat>
176 void getElementMatrixOptimized(CG const&, QR const&, RN const&, CN const&,
177 LocalFct const&, Mat&, tag::matrix) const
178 {
179 DUNE_THROW(Dune::NotImplemented,
180 "Optimized LocalOperator with matrix coefficients not implemented.");
181 }
182 };
183
184 template <class LC>
185 struct GridFunctionOperatorRegistry<tag::testvec_trialvec, LC>
186 {
187 static constexpr int degree = 0;
189 };
190
193} // end namespace AMDiS
zero-order operator , or
Definition ZeroOrderTestvecTrialvec.hpp:27
Registry to specify a tag for each implementation type.
Definition GridFunctionOperator.hpp:216
Definition Tags.hpp:11
Definition Tags.hpp:9
Definition ZeroOrderTestvecTrialvec.hpp:21