AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
ConvectionDiffusionOperator.hpp
1#pragma once
2
3#include <type_traits>
4
5#include <dune/common/hybridutilities.hh>
6#include <dune/common/typetree/nodeconcepts.hh>
7
8#include <amdis/LocalOperator.hpp>
9#include <amdis/VariadicGridFunctionOperator.hpp>
10#include <amdis/common/Order.hpp>
11#include <amdis/common/Logical.hpp>
12#include <amdis/common/StaticSize.hpp>
13#include <amdis/common/TypeTraits.hpp>
14#include <amdis/typetree/FiniteElementType.hpp>
15
16namespace AMDiS
17{
18 namespace tag{
19 template<bool conserving>
21 }
29 template <bool conserving>
31 {
32 public:
34 {}
35
36 template <class CG, class RN, class CN, class Quad, class LocalFctA, class LocalFctB, class LocalFctC, class LocalFctF, class Mat>
37 void assemble(CG const& contextGeo, RN const& rowNode, CN const& colNode, Quad const& quad,
38 LocalFctA const& localFctA_, LocalFctB const& localFctB_, LocalFctC const& localFctC_, LocalFctF const& localFctF_,
39 Mat& elementMatrix) const
40 {
41 using A_range = typename LocalFctA::Range;
42 static_assert(static_size_v<A_range> == 1 || (static_num_rows_v<A_range> == CG::dow && static_num_cols_v<A_range> == CG::dow),
43 "Expression A(x) must be of scalar or matrix type.");
44 using b_range = typename LocalFctB::Range;
45 static_assert(static_size_v<b_range> == 1 || static_size_v<b_range> == CG::dow,
46 "Expression b(x) must be of scalar or vector type.");
47 using c_range = typename LocalFctC::Range;
48 static_assert(static_size_v<c_range> == 1,
49 "Expression c(x) must be of scalar type.");
50 static_assert(Dune::TypeTree::Concept::LeafTreeNode<RN> && Dune::TypeTree::Concept::LeafTreeNode<CN>,
51 "Operator can be applied to Leaf-Nodes only." );
52 static_assert(std::is_same_v<FiniteElementType_t<RN>, FiniteElementType_t<CN>>,
53 "Galerkin operator requires equal finite elements for test and trial space." );
54
55 using RangeFieldType = typename RN::LocalBasis::Traits::RangeFieldType;
56 using WorldVector = FieldVector<RangeFieldType,CG::dow>;
57 std::vector<WorldVector> gradients;
58
59 std::size_t numLocalFe = rowNode.size();
60
61 auto contextGeometry = contextGeo.geometry();
62
63 for (std::size_t iq = 0; iq < quad.size(); ++iq) {
64 // Position of the current quadrature point in the reference element
65 auto&& local = contextGeo.coordinateInElement(quad[iq].position());
66
67 // The transposed inverse Jacobian of the map from the reference element to the element
68 const auto jacobian = contextGeo.elementGeometry().jacobianInverseTransposed(local);
69
70 // The multiplicative factor in the integral transformation formula
71 const auto factor = contextGeo.integrationElement(quad[iq].position()) * quad[iq].weight();
72
73 // the values of the shape functions on the reference element at the quadrature point
74 auto const& shapeValues = rowNode.localBasisValuesAt(local);
75
76 // The gradients of the shape functions on the reference element
77 auto const& shapeGradients = rowNode.localBasisJacobiansAt(local);
78
79 // Compute the shape function gradients on the real element
80 gradients.resize(shapeGradients.size());
81
82 for (std::size_t i = 0; i < gradients.size(); ++i)
83 jacobian.mv(shapeGradients[i][0], gradients[i]);
84
85 const auto A = localFctA_(local);
86 WorldVector b = makeB<CG::dow>(localFctB_(local));
87 const auto c = localFctC_(local);
88
89 if (conserving) {
90 WorldVector gradAi;
91 for (std::size_t i = 0; i < numLocalFe; ++i) {
92 const auto local_i = rowNode.localIndex(i);
93 multiply(Dune::MatVec::as_scalar(A), gradients[i], gradAi);
94 auto gradBi = b * gradients[i];
95 for (std::size_t j = 0; j < numLocalFe; ++j) {
96 const auto local_j = colNode.localIndex(j);
97 elementMatrix[local_i][local_j] += (dot(gradAi, gradients[j]) + (c * shapeValues[i] - gradBi) * shapeValues[j]) * factor;
98 }
99 }
100 } else {
101 WorldVector grad_i;
102 for (std::size_t i = 0; i < numLocalFe; ++i) {
103 const auto local_i = rowNode.localIndex(i);
104 multiply(Dune::MatVec::as_scalar(A), gradients[i], grad_i);
105 grad_i.axpy(shapeValues[i], b);
106 for (std::size_t j = 0; j < numLocalFe; ++j) {
107 const auto local_j = colNode.localIndex(j);
108 elementMatrix[local_i][local_j] += (dot(grad_i, gradients[j]) + c * shapeValues[i] * shapeValues[j]) * factor;
109 }
110 }
111 }
112 }
113 }
114
115
116 template <class CG, class Node, class Quad, class LocalFctA, class LocalFctB, class LocalFctC, class LocalFctF, class Vec>
117 void assemble(CG const& contextGeo, Node const& node, Quad const& quad, LocalFctA const& localFctA_, LocalFctB const& localFctB_, LocalFctC const& localFctC_, LocalFctF const& localFctF_, Vec& elementVector) const
118 {
119 using f_range = typename LocalFctF::Range;
120 static_assert( static_size_v<f_range> == 1,
121 "Expression f(x) must be of scalar type." );
122 static_assert(Dune::TypeTree::Concept::LeafTreeNode<Node>,
123 "Operator can be applied to Leaf-Nodes only." );
124
125 std::size_t numLocalFe = node.size();
126
127 for (std::size_t iq = 0; iq < quad.size(); ++iq) {
128 // Position of the current quadrature point in the reference element
129 auto&& local = contextGeo.coordinateInElement(quad[iq].position());
130
131 // the values of the shape functions on the reference element at the quadrature point
132 auto const& shapeValues = node.localBasisValuesAt(local);
133
134 // The multiplicative factor in the integral transformation formula
135 const auto factor = localFctF_(local) * contextGeo.integrationElement(quad[iq].position()) * quad[iq].weight();
136
137 for (std::size_t i = 0; i < numLocalFe; ++i) {
138 const auto local_i = node.localIndex(i);
139 elementVector[local_i] += shapeValues[i] * factor;
140 }
141 }
142 }
143
144 private:
145
146 template <class T1, int M, int N, class T2, class T3>
147 static void multiply(FieldMatrix<T1,M,N> const& A, FieldVector<T2,N> const& x, FieldVector<T3,M>& y)
148 {
149 A.mv(x, y);
150 }
151
152 template <class T1, int N, class T2, class T3,
153 std::enable_if_t<Dune::IsNumber<T1>::value, int> = 0>
154 static void multiply(T1 const& alpha, FieldVector<T2,N> const& x, FieldVector<T3,N>& y)
155 {
156 y = x;
157 y *= alpha;
158 }
159
160 template <int dow, class T, int N>
161 static FieldVector<T,dow> makeB(FieldVector<T,N> const& b) { return b; }
162
163 template <int dow, class T>
164 static FieldVector<T,dow> makeB(FieldVector<T,1> const& b) { return FieldVector<T,dow>(T(b)); }
165
166 template <int dow, class T,
167 std::enable_if_t<Dune::IsNumber<T>::value, int> = 0>
168 static FieldVector<T,dow> makeB(T b) { return FieldVector<T,dow>(b); }
169
170 private:
171 int quadOrder_;
172 };
173
174 template <class LC, bool conserving>
176 {
177 static constexpr int degree = 2;
179 };
180
182 template <bool conserving = true,
183 class PreGridFctA, class PreGridFctB, class PreGridFctC, class PreGridFctF>
184 auto convectionDiffusion(PreGridFctA const& gridFctA, PreGridFctB const& gridFctB,
185 PreGridFctC const& gridFctC, PreGridFctF const& gridFctF,
186 int gridFctDegree = -1, bool_t<conserving> = {})
187 {
188 return makeOperator(tag::convectionDiffusion<conserving>{}, gridFctDegree, gridFctA,gridFctB,gridFctC, gridFctF);
189 }
190
193} // end namespace AMDiS
Definition ConvectionDiffusionOperator.hpp:31
auto convectionDiffusion(PreGridFctA const &gridFctA, PreGridFctB const &gridFctB, PreGridFctC const &gridFctC, PreGridFctF const &gridFctF, int gridFctDegree=-1, bool_t< conserving >={})
Generator function for constructing a Convection-Diffusion Operator.
Definition ConvectionDiffusionOperator.hpp:184
auto makeOperator(Tag const &tag, Expr &&expr, int gridFctDeg=-1)
Definition GridFunctionOperator.hpp:235
Registry to specify a tag for each implementation type.
Definition VariadicGridFunctionOperator.hpp:246
Definition ConvectionDiffusionOperator.hpp:20