AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
SecondOrderGradTestGradTrial.hpp
1#pragma once
2
3#include <type_traits>
4
5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/Output.hpp>
7#include <amdis/common/StaticSize.hpp>
8#include <amdis/common/ValueCategory.hpp>
9#include <amdis/typetree/FiniteElementType.hpp>
10
11#include <dune/common/typetree/nodeconcepts.hh>
12
13namespace AMDiS
14{
20 namespace tag
21 {
23 }
24
25
28 {
29 public:
31
32 template <class CG, class RN, class CN, class Quad, class LocalFct, class Mat>
33 void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode,
34 Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
35 {
36 using expr_value_type = typename LocalFct::Range;
37 static_assert(static_size_v<expr_value_type> == 1 ||
38 (static_num_rows_v<expr_value_type> == CG::dow &&
39 static_num_cols_v<expr_value_type> == CG::dow),
40 "Expression must be of scalar or matrix type.");
41 static_assert(Dune::TypeTree::Concept::LeafTreeNode<RN> && Dune::TypeTree::Concept::LeafTreeNode<CN>,
42 "Operator can be applied to Leaf-Nodes only.");
43
44 const bool sameFEType = std::is_same_v<FiniteElementType_t<RN>, FiniteElementType_t<CN>>;
45 using Category = ValueCategory_t<typename LocalFct::Range>;
46
47 if constexpr (sameFEType) {
48 const bool sameNode = rowNode.treeIndex() == colNode.treeIndex();
49 if (sameNode)
50 return getElementMatrixOptimized(contextGeo, quad, rowNode, colNode, localFct, elementMatrix, Category{});
51 else {
52 const bool sameFEsize = rowNode.finiteElement().size() == colNode.finiteElement().size();
53 if (sameFEsize)
54 return getElementMatrixStandard(contextGeo, quad, rowNode, colNode, localFct, elementMatrix);
55 }
56 }
57
58 error_exit("Not implemented: currently only the implementation for equal fespaces available");
59 }
60
61 protected:
62
63 template <class CG, class QR, class RN, class CN, class LocalFct, class Mat>
64 void getElementMatrixStandard(CG const& contextGeo, QR const& quad,
65 RN const& rowNode, CN const& colNode,
66 LocalFct const& localFct, Mat& elementMatrix) const
67 {
68 std::size_t size = rowNode.size();
69
70 using RangeFieldType = typename RN::LocalBasis::Traits::RangeFieldType;
71 using WorldVector = Dune::FieldVector<RangeFieldType,CG::dow>;
72 std::vector<WorldVector> gradients;
73
74 for (auto const& qp : quad) {
75 // Position of the current quadrature point in the reference element
76 auto&& local = contextGeo.coordinateInElement(qp.position());
77
78 // The transposed inverse Jacobian of the map from the reference element to the element
79 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
80
81 // The multiplicative factor in the integral transformation formula
82 const auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
83 const auto exprValue = localFct(local);
84
85 // The gradients of the shape functions on the reference element
86 auto const& shapeGradients = rowNode.localBasisJacobiansAt(local);
87
88 // Compute the shape function gradients on the real element
89 gradients.resize(shapeGradients.size());
90
91 for (std::size_t i = 0; i < gradients.size(); ++i)
92 jacobian.mv(shapeGradients[i][0], gradients[i]);
93
94 for (std::size_t i = 0; i < size; ++i) {
95 const auto local_i = rowNode.localIndex(i);
96 for (std::size_t j = 0; j < size; ++j) {
97 const auto local_j = colNode.localIndex(j);
98 elementMatrix[local_i][local_j] += eval(exprValue, factor, gradients[i], gradients[j]);
99 }
100 }
101 }
102 }
103
104 template <class CG, class QR, class RN, class CN, class LocalFct, class Mat>
105 void getElementMatrixOptimized(CG const& contextGeo, QR const& quad,
106 RN const& node, CN const& /*colNode*/,
107 LocalFct const& localFct, Mat& elementMatrix, tag::scalar) const
108 {
109 std::size_t size = node.size();
110
111 using RangeFieldType = typename RN::LocalBasis::Traits::RangeFieldType;
112 using WorldVector = Dune::FieldVector<RangeFieldType,CG::dow>;
113 std::vector<WorldVector> gradients;
114
115 for (auto const& qp : quad) {
116 // Position of the current quadrature point in the reference element
117 auto&& local = contextGeo.coordinateInElement(qp.position());
118
119 // The transposed inverse Jacobian of the map from the reference element to the element
120 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
121
122 // The multiplicative factor in the integral transformation formula
123 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
124
125 // The gradients of the shape functions on the reference element
126 auto const& shapeGradients = node.localBasisJacobiansAt(local);
127
128 // Compute the shape function gradients on the real element
129 gradients.resize(shapeGradients.size());
130 for (std::size_t i = 0; i < gradients.size(); ++i)
131 jacobian.mv(shapeGradients[i][0], gradients[i]);
132
133 for (std::size_t i = 0; i < size; ++i) {
134 const auto local_i = node.localIndex(i);
135
136 elementMatrix[local_i][local_i] += factor * (gradients[i] * gradients[i]);
137
138 for (std::size_t j = i+1; j < size; ++j) {
139 const auto local_j = node.localIndex(j);
140 const auto value = factor * (gradients[i] * gradients[j]);
141
142 elementMatrix[local_i][local_j] += value;
143 elementMatrix[local_j][local_i] += value;
144 }
145 }
146 }
147 }
148
149 template <class CG, class QR, class RN, class CN, class LocalFct, class Mat>
150 void getElementMatrixOptimized(CG const& contextGeo, QR const& quad,
151 RN const& node, CN const& /*colNode*/,
152 LocalFct const& localFct, Mat& elementMatrix, tag::matrix) const
153 {
154 std::size_t size = node.size();
155
156 using RangeFieldType = typename RN::LocalBasis::Traits::RangeFieldType;
157 using WorldVector = Dune::FieldVector<RangeFieldType,CG::dow>;
158 std::vector<WorldVector> gradients;
159
160 for (auto const& qp : quad) {
161 // Position of the current quadrature point in the reference element
162 auto&& local = contextGeo.coordinateInElement(qp.position());
163
164 // The transposed inverse Jacobian of the map from the reference element to the element
165 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
166
167 // The multiplicative factor in the integral transformation formula
168 const auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
169 const auto exprValue = localFct(local);
170
171 // The gradients of the shape functions on the reference element
172 auto const& shapeGradients = node.localBasisJacobiansAt(local);
173
174 // Compute the shape function gradients on the real element
175 gradients.resize(shapeGradients.size());
176 for (std::size_t i = 0; i < gradients.size(); ++i)
177 jacobian.mv(shapeGradients[i][0], gradients[i]);
178
179 for (std::size_t i = 0; i < size; ++i) {
180 const auto local_i = node.localIndex(i);
181 for (std::size_t j = 0; j < size; ++j) {
182 const auto local_j = node.localIndex(j);
183 elementMatrix[local_i][local_j] += eval(exprValue, factor, gradients[i], gradients[j]);
184 }
185 }
186 }
187 }
188
189 protected:
190
191 template <class S, class F, class T, int dow,
192 std::enable_if_t<Category::Scalar<S>,int> = 0>
193 T eval(S const& scalar, F factor,
194 Dune::FieldVector<T,dow> const& grad_test,
195 Dune::FieldVector<T,dow> const& grad_trial) const
196 {
197 return (scalar * factor) * (grad_test * grad_trial);
198 }
199
200 template <class M, class F, class T, int dow,
201 std::enable_if_t<Category::Matrix<M>,int> = 0>
202 T eval(M const& mat, F factor,
203 Dune::FieldVector<T,dow> const& grad_test,
204 Dune::FieldVector<T,dow> const& grad_trial) const
205 {
206 return factor * (grad_test * (mat * grad_trial));
207 }
208 };
209
210 template <class LC>
211 struct GridFunctionOperatorRegistry<tag::gradtest_gradtrial, LC>
212 {
213 static constexpr int degree = 2;
215 };
216
218 template <class Expr>
219 auto sot(Expr&& expr, int quadOrder = -1)
220 {
221 return makeOperator(tag::gradtest_gradtrial{}, FWD(expr), quadOrder);
222 }
223
226} // end namespace AMDiS
second-order operator , or
Definition SecondOrderGradTestGradTrial.hpp:28
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 SecondOrderGradTestGradTrial.hpp:22
Definition Tags.hpp:11
Definition Tags.hpp:9