AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
SecondOrderDivTestvecDivTrialvec.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/childaccess.hh>
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 static_assert(static_size_v<typename LocalFct::Range> == 1,
36 "Expression must be of scalar type.");
37 static_assert(Dune::TypeTree::Concept::UniformInnerTreeNode<RN> && Dune::TypeTree::Concept::UniformInnerTreeNode<CN>,
38 "Operator can be applied to Power-Nodes only.");
39
40 const bool sameFE = std::is_same_v<FiniteElementType_t<RN>, FiniteElementType_t<CN>>;
41 const bool sameNode = rowNode.treeIndex() == colNode.treeIndex();
42
43 if (sameFE && sameNode)
44 getElementMatrixOptimized(contextGeo, quad, rowNode, colNode, localFct, elementMatrix);
45 else
46 getElementMatrixStandard(contextGeo, quad, rowNode, colNode, localFct, elementMatrix);
47 }
48
49 protected:
50 template <class CG, class QR, class RN, class CN, class LocalFct, class Mat>
51 void getElementMatrixStandard(CG const& contextGeo, QR const& quad,
52 RN const& rowNode, CN const& colNode,
53 LocalFct const& localFct, Mat& elementMatrix) const
54 {
55 assert( rowNode.degree() == colNode.degree() );
56
57 std::size_t rowSize = rowNode.child(0).size();
58 std::size_t colSize = colNode.child(0).size();
59
60 using RowFieldType = typename Dune::TypeTree::Child<RN,0>::LocalBasis::Traits::RangeFieldType;
61 using RowWorldVector = FieldVector<RowFieldType,CG::dow>;
62 std::vector<RowWorldVector> rowGradients;
63
64 using ColFieldType = typename Dune::TypeTree::Child<CN,0>::LocalBasis::Traits::RangeFieldType;
65 using ColWorldVector = FieldVector<ColFieldType,CG::dow>;
66 std::vector<ColWorldVector> colGradients;
67
68 for (auto const& qp : quad) {
69 // Position of the current quadrature point in the reference element
70 decltype(auto) local = contextGeo.coordinateInElement(qp.position());
71
72 // The transposed inverse Jacobian of the map from the reference element to the element
73 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
74
75 // The multiplicative factor in the integral transformation formula
76 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
77
78 // The gradients of the shape functions on the reference element
79 auto const& rowShapeGradients = rowNode.child(0).localBasisJacobiansAt(local);
80 auto const& colShapeGradients = colNode.child(0).localBasisJacobiansAt(local);
81
82 // Compute the shape function gradients on the real element
83 rowGradients.resize(rowShapeGradients.size());
84
85 for (std::size_t i = 0; i < rowGradients.size(); ++i)
86 jacobian.mv(rowShapeGradients[i][0], rowGradients[i]);
87
88 colGradients.resize(colShapeGradients.size());
89
90 for (std::size_t i = 0; i < colGradients.size(); ++i)
91 jacobian.mv(colShapeGradients[i][0], colGradients[i]);
92
93 for (std::size_t i = 0; i < rowSize; ++i) {
94 for (std::size_t j = 0; j < colSize; ++j) {
95 for (std::size_t k = 0; k < CG::dow; ++k) {
96 for (std::size_t l = 0; l < CG::dow; ++l) {
97 const auto local_ki = rowNode.child(k).localIndex(i);
98 const auto local_kj = colNode.child(l).localIndex(j);
99 elementMatrix[local_ki][local_kj] += factor * rowGradients[i][k] * colGradients[j][l];
100 }
101 }
102 }
103 }
104 }
105 }
106
107
108 template <class CG, class QR, class Node, class LocalFct, class Mat>
109 void getElementMatrixOptimized(CG const& contextGeo, QR const& quad,
110 Node const& node, Node const& /*colNode*/,
111 LocalFct const& localFct, Mat& elementMatrix) const
112 {
113 std::size_t size = node.child(0).size();
114
115 using RangeFieldType = typename Dune::TypeTree::Child<Node,0>::LocalBasis::Traits::RangeFieldType;
116 using WorldVector = FieldVector<RangeFieldType,CG::dow>;
117 std::vector<WorldVector> gradients;
118
119 for (auto const& qp : quad) {
120 // Position of the current quadrature point in the reference element
121 auto&& local = contextGeo.coordinateInElement(qp.position());
122
123 // The transposed inverse Jacobian of the map from the reference element to the element
124 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
125
126 // The multiplicative factor in the integral transformation formula
127 const auto factor = localFct(local) * contextGeo.integrationElement(qp.position()) * qp.weight();
128
129 // The gradients of the shape functions on the reference element
130 auto const& shapeGradients = node.child(0).localBasisJacobiansAt(local);
131
132 // Compute the shape function gradients on the real element
133 gradients.resize(shapeGradients.size());
134
135 for (std::size_t i = 0; i < gradients.size(); ++i)
136 jacobian.mv(shapeGradients[i][0], gradients[i]);
137
138 for (std::size_t k = 0; k < CG::dow; ++k) {
139 for (std::size_t l = 0; l < CG::dow; ++l) {
140 for (std::size_t i = 0; i < size; ++i) {
141 const auto local_ki = node.child(k).localIndex(i);
142 const auto value = factor * gradients[i][k];
143 elementMatrix[local_ki][local_ki] += value * gradients[i][k];
144
145 for (std::size_t j = i + 1; j < size; ++j) {
146 const auto local_kj = node.child(l).localIndex(j);
147 elementMatrix[local_ki][local_kj] += value * gradients[j][l];
148 elementMatrix[local_kj][local_ki] += value * gradients[j][l];
149 }
150 }
151 }
152 }
153 }
154 }
155 };
156
157 template <class LC>
158 struct GridFunctionOperatorRegistry<tag::divtestvec_divtrialvec, LC>
159 {
160 static constexpr int degree = 2;
162 };
163
166} // end namespace AMDiS
second-order operator
Definition SecondOrderDivTestvecDivTrialvec.hpp:27
Registry to specify a tag for each implementation type.
Definition GridFunctionOperator.hpp:216
Definition SecondOrderDivTestvecDivTrialvec.hpp:21