AMDiS 2.10
The Adaptive Multi-Dimensional Simulation Toolbox
StokesOperator.hpp
1#pragma once
2
3#include <type_traits>
4#include <vector>
5
6#include <amdis/GridFunctionOperator.hpp>
7#include <amdis/common/StaticSize.hpp>
8
9namespace AMDiS
10{
16 namespace tag
17 {
18 struct stokes
19 {
20 bool constViscosity = false;
21 };
22 }
23
24
27 {
28 bool constViscosity_;
29
30 public:
32 : constViscosity_(t.constViscosity)
33 {}
34
35 template <class CG, class Node, class Quad, class LocalFct, class Mat>
36 void assemble(CG const& contextGeo, Node const& tree, Node const& /*colNode*/,
37 Quad const& quad, LocalFct const& localFct, Mat& elementMatrix) const
38 {
39 static_assert(static_size_v<typename LocalFct::Range> == 1,
40 "Viscosity must be of scalar type.");
41
42 using VelocityNode = Dune::TypeTree::Child<Node,0>;
43 using VelocityCompNode = Dune::TypeTree::Child<Node,0,0>;
44 using PressureNode = Dune::TypeTree::Child<Node,1>;
45
46 static_assert(VelocityNode::isPower && PressureNode::isLeaf, "");
47
48 using namespace Dune::Indices;
49
50 std::size_t numVelocityLocalFE = tree.child(_0,0).size();
51 std::size_t numPressureLocalFE = tree.child(_1).size();
52
53 using RangeFieldType = typename VelocityNode::ChildType::LocalBasis::Traits::RangeFieldType;
54 using WorldVector = Dune::FieldVector<RangeFieldType,CG::dow>;
55 std::vector<WorldVector> gradients;
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 transposed inverse Jacobian of the map from the reference element to the element
62 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
63
64 // The multiplicative factor in the integral transformation formula
65 const auto factor = contextGeo.integrationElement(qp.position()) * qp.weight();
66 const auto vel_factor = localFct(local) * factor;
67
68 auto const& shapeGradients = tree.child(_0,0).localBasisJacobiansAt(local);
69
70 // Compute the shape function gradients on the real element
71 gradients.resize(shapeGradients.size());
72 for (std::size_t i = 0; i < gradients.size(); ++i)
73 jacobian.mv(shapeGradients[i][0], gradients[i]);
74
75 auto const& pressureValues = tree.child(_1).localBasisValuesAt(local);
76
77 // <viscosity * grad(u), grad(v)>
78 for (std::size_t i = 0; i < numVelocityLocalFE; ++i) {
79 const auto value_ii = vel_factor * (gradients[i] * gradients[i]);
80 for (std::size_t k = 0; k < CG::dow; ++k) {
81 const auto local_ki = tree.child(_0,k).localIndex(i);
82 elementMatrix[local_ki][local_ki] += value_ii;
83 }
84
85 for (std::size_t j = i+1; j < numVelocityLocalFE; ++j) {
86 const auto value = vel_factor * (gradients[i] * gradients[j]);
87 for (std::size_t k = 0; k < CG::dow; ++k) {
88 const auto local_ki = tree.child(_0,k).localIndex(i);
89 const auto local_kj = tree.child(_0,k).localIndex(j);
90 elementMatrix[local_ki][local_kj] += value;
91 elementMatrix[local_kj][local_ki] += value;
92 }
93 }
94 }
95
96 if (!constViscosity_) {
97 // <viscosity * d_i u_j, d_j v_i>
98 for (std::size_t i = 0; i < numVelocityLocalFE; ++i) {
99 for (std::size_t kj = 0; kj < CG::dow; ++kj) {
100 const auto value_i_kj = vel_factor * gradients[i][kj];
101 for (std::size_t j = 0; j < numVelocityLocalFE; ++j) {
102 for (std::size_t ki = 0; ki < CG::dow; ++ki) {
103 const auto local_ki = tree.child(_0,ki).localIndex(i);
104 const auto local_kj = tree.child(_0,kj).localIndex(j);
105 elementMatrix[local_ki][local_kj] += value_i_kj * gradients[j][ki];
106 }
107 }
108 }
109 }
110 }
111
112 // <p, div(v)> + <div(u), q>
113 for (std::size_t i = 0; i < numVelocityLocalFE; ++i) {
114 for (std::size_t j = 0; j < numPressureLocalFE; ++j) {
115 const auto value = factor * pressureValues[j];
116 for (std::size_t k = 0; k < CG::dow; ++k) {
117 const auto local_v = tree.child(_0,k).localIndex(i);
118 const auto local_p = tree.child(_1).localIndex(j);
119
120 elementMatrix[local_v][local_p] += gradients[i][k] * value;
121 elementMatrix[local_p][local_v] += gradients[i][k] * value;
122 }
123 }
124 }
125 }
126
127 } // end assemble
128 };
129
130 template <class LC>
131 struct GridFunctionOperatorRegistry<tag::stokes, LC>
132 {
133 static constexpr int degree = 2;
134 using type = StokesOperator;
135 };
136
139} // end namespace AMDiS
Stokes operator .
Definition: StokesOperator.hpp:27
Registry to specify a tag for each implementation type.
Definition: GridFunctionOperator.hpp:216
Definition: StokesOperator.hpp:19