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