AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
ConvectionDiffusionOperator.hpp
1 #pragma once
2 
3 #include <type_traits>
4 
5 #include <dune/common/hybridutilities.hh>
6 
7 #include <amdis/LocalOperator.hpp>
8 #include <amdis/Output.hpp>
9 #include <amdis/common/TypeTraits.hpp>
10 #include <amdis/common/StaticSize.hpp>
11 
12 namespace AMDiS
13 {
19  template <class LC, class GridFctA, class GridFctB, class GridFctC, class GridFctF, bool conserving = true>
24  : public LocalOperator<ConvectionDiffusionOperator<LC, GridFctA, GridFctB, GridFctC, GridFctF, conserving>, LC>
25  {
26  static const int dow = ContextGeometry<LC>::dow;
27 
28  using A_range = typename GridFctA::Range;
29  static_assert( static_size_v<A_range> == 1 || (static_num_rows_v<A_range> == dow && static_num_cols_v<A_range> == dow),
30  "Expression A(x) must be of scalar or matrix type." );
31  using b_range = typename GridFctB::Range;
32  static_assert( static_size_v<b_range> == 1 || static_size_v<b_range> == dow,
33  "Expression b(x) must be of scalar or vector type." );
34  using c_range = typename GridFctC::Range;
35  static_assert( static_size_v<c_range> == 1,
36  "Expression c(x) must be of scalar type." );
37  using f_range = typename GridFctF::Range;
38  static_assert( static_size_v<f_range> == 1,
39  "Expression f(x) must be of scalar type." );
40 
41  public:
42  ConvectionDiffusionOperator(GridFctA const& gridFctA, GridFctB const& gridFctB,
43  GridFctC const& gridFctC, GridFctF const& gridFctF)
44  : gridFctA_(gridFctA)
45  , gridFctB_(gridFctB)
46  , gridFctC_(gridFctC)
47  , gridFctF_(gridFctF)
48  {}
49 
50  template <class CG, class RN, class CN, class Mat>
51  void getElementMatrix(CG const& contextGeo, RN const& rowNode, CN const& colNode, Mat& elementMatrix)
52  {
53  static_assert(RN::isLeaf && CN::isLeaf,
54  "Operator can be applied to Leaf-Nodes only." );
55 
56  static_assert(std::is_same_v<FiniteElementType_t<RN>, FiniteElementType_t<CN>>,
57  "Galerkin operator requires equal finite elements for test and trial space." );
58 
59  using RangeFieldType = typename RN::LocalBasis::Traits::RangeFieldType;
60  using WorldVector = FieldVector<RangeFieldType,CG::dow>;
61  std::vector<WorldVector> gradients;
62 
63  auto localFctA = localFunction(gridFctA_); localFctA.bind(contextGeo.element());
64  auto localFctB = localFunction(gridFctB_); localFctB.bind(contextGeo.element());
65  auto localFctC = localFunction(gridFctC_); localFctC.bind(contextGeo.element());
66 
67  std::size_t numLocalFe = rowNode.size();
68 
69  int quadDeg = std::max({this->getDegree(2,coeffOrder(localFctA),rowNode,colNode),
70  this->getDegree(1,coeffOrder(localFctB),rowNode,colNode),
71  this->getDegree(0,coeffOrder(localFctC),rowNode,colNode)});
72 
73  using QuadratureRules = Dune::QuadratureRules<typename CG::Geometry::ctype, CG::LocalContext::mydimension>;
74  auto const& quad = QuadratureRules::rule(contextGeo.type(), quadDeg);
75 
76  for (std::size_t iq = 0; iq < quad.size(); ++iq) {
77  // Position of the current quadrature point in the reference element
78  auto&& local = contextGeo.local(quad[iq].position());
79 
80  // The transposed inverse Jacobian of the map from the reference element to the element
81  const auto jacobian = contextGeo.geometry().jacobianInverseTransposed(local);
82 
83  // The multiplicative factor in the integral transformation formula
84  const auto factor = contextGeo.integrationElement(quad[iq].position()) * quad[iq].weight();
85 
86  // the values of the shape functions on the reference element at the quadrature point
87  auto const& shapeValues = rowNode.localBasisValuesAt(local);
88 
89  // The gradients of the shape functions on the reference element
90  auto const& shapeGradients = rowNode.localBasisJacobiansAt(local);
91 
92  // Compute the shape function gradients on the real element
93  gradients.resize(shapeGradients.size());
94 
95  for (std::size_t i = 0; i < gradients.size(); ++i)
96  jacobian.mv(shapeGradients[i][0], gradients[i]);
97 
98  const auto A = localFctA(local);
99  WorldVector b = makeB(localFctB(local));
100  const auto c = localFctC(local);
101 
102  if (conserving) {
103  WorldVector gradAi;
104  for (std::size_t i = 0; i < numLocalFe; ++i) {
105  const auto local_i = rowNode.localIndex(i);
106  gradAi = A * gradients[i];
107  auto gradBi = b * gradients[i];
108  for (std::size_t j = 0; j < numLocalFe; ++j) {
109  const auto local_j = colNode.localIndex(j);
110  elementMatrix[local_i][local_j] += (dot(gradAi, gradients[j]) + (c * shapeValues[i] - gradBi) * shapeValues[j]) * factor;
111  }
112  }
113  } else {
114  WorldVector grad_i;
115  for (std::size_t i = 0; i < numLocalFe; ++i) {
116  const auto local_i = rowNode.localIndex(i);
117  grad_i = A * gradients[i];
118  grad_i+= b * shapeValues[i];
119  for (std::size_t j = 0; j < numLocalFe; ++j) {
120  const auto local_j = colNode.localIndex(j);
121  elementMatrix[local_i][local_j] += (dot(grad_i, gradients[j]) + c * shapeValues[i] * shapeValues[j]) * factor;
122  }
123  }
124  }
125  }
126 
127  localFctA.unbind();
128  localFctB.unbind();
129  localFctC.unbind();
130  }
131 
132 
133  template <class CG, class Node, class Vec>
134  void getElementVector(CG const& contextGeo, Node const& node, Vec& elementVector)
135  {
136  static_assert(Node::isLeaf,
137  "Operator can be applied to Leaf-Nodes only." );
138 
139  auto localFctF = localFunction(gridFctF_); localFctF.bind(contextGeo.element());
140 
141  std::size_t numLocalFe = node.size();
142 
143  int quad_order = this->getDegree(0,coeffOrder(localFctF),node);
144 
145  using QuadratureRules = Dune::QuadratureRules<typename CG::Geometry::ctype, CG::LocalContext::dimension>;
146  auto const& quad = QuadratureRules::rule(contextGeo.type(), quad_order);
147 
148  for (std::size_t iq = 0; iq < quad.size(); ++iq) {
149  // Position of the current quadrature point in the reference element
150  auto&& local = contextGeo.local(quad[iq].position());
151 
152  // the values of the shape functions on the reference element at the quadrature point
153  auto const& shapeValues = node.localBasisValuesAt(local);
154 
155  // The multiplicative factor in the integral transformation formula
156  const auto factor = localFctF(local) * contextGeo.integrationElement(quad[iq].position()) * quad[iq].weight();
157 
158  for (std::size_t i = 0; i < numLocalFe; ++i) {
159  const auto local_i = node.localIndex(i);
160  elementVector[local_i] += shapeValues[i] * factor;
161  }
162  }
163 
164  localFctF.unbind();
165  }
166 
167  private:
168 
169  template <class LocalFct>
170  int coeffOrder(LocalFct const& localFct)
171  {
172  if constexpr (Concepts::Polynomial<LocalFct>)
173  return order(localFct);
174  else
175  return 0;
176  }
177 
178  template <class T, int N>
179  static FieldVector<T,dow> makeB(FieldVector<T,N> const& b) { return b; }
180 
181  template <class T, int N>
182  static FieldVector<T,dow> makeB(FieldVector<T,N>&& b) { return std::move(b); }
183 
184  template <class T>
185  static FieldVector<T,dow> makeB(FieldVector<T,1> const& b) { return FieldVector<T,dow>(T(b)); }
186 
187  template <class T>
188  static FieldVector<T,dow> makeB(FieldVector<T,1>&& b) { return FieldVector<T,dow>(T(b)); }
189 
190  template <class T, std::enable_if_t<std::is_arithmetic_v<T>, int> = 0>
191  static FieldVector<T,dow> makeB(T b) { return FieldVector<T,dow>(b); }
192 
193  private:
194  GridFctA gridFctA_;
195  GridFctB gridFctB_;
196  GridFctC gridFctC_;
197  GridFctF gridFctF_;
198  };
199 
200 #ifndef DOXYGEN
201  template <class PreGridFctA, class PreGridFctB, class PreGridFctC, class PreGridFctF, class c>
203  {
204  static constexpr bool conserving = c::value;
205  PreGridFctA gridFctA;
206  PreGridFctB gridFctB;
207  PreGridFctC gridFctC;
208  PreGridFctF gridFctF;
209  };
210 #endif
211 
212 
213  template <class PreGridFctA, class PreGridFctB, class PreGridFctC, class PreGridFctF, bool conserving = true>
214  auto convectionDiffusion(PreGridFctA const& gridFctA, PreGridFctB const& gridFctB,
215  PreGridFctC const& gridFctC, PreGridFctF const& gridFctF,
216  bool_t<conserving> = {})
217  {
219  return Pre{gridFctA, gridFctB, gridFctC, gridFctF};
220  }
221 
222 
223 #ifndef DOXYGEN
224  template <class Context, class... T, class GridView>
225  auto makeLocalOperator(PreConvectionDiffusionOperator<T...> pre, GridView const& gridView)
226  {
227  using Pre = PreConvectionDiffusionOperator<T...>;
228  auto gridFctA = makeGridFunction(std::move(pre.gridFctA), gridView);
229  auto gridFctB = makeGridFunction(std::move(pre.gridFctB), gridView);
230  auto gridFctC = makeGridFunction(std::move(pre.gridFctC), gridView);
231  auto gridFctF = makeGridFunction(std::move(pre.gridFctF), gridView);
232 
233  using GridFctOp = ConvectionDiffusionOperator<Context,
234  decltype(gridFctA), decltype(gridFctB), decltype(gridFctC), decltype(gridFctF), Pre::conserving>;
235 
236  GridFctOp localOperator{std::move(gridFctA), std::move(gridFctB), std::move(gridFctC), std::move(gridFctF)};
237  return localOperator;
238  }
239 #endif
240 
243 } // end namespace AMDiS
Definition: ConvectionDiffusionOperator.hpp:202
decltype(auto) makeGridFunction(PreGridFct const &preGridFct, GridView const &gridView)
Generator for Gridfunctions from Expressions (PreGridfunctions)
Definition: GridFunction.hpp:154
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
int getDegree(int derivOrder, int coeffDegree, RN const &rowNode, CN const &colNode) const
Return the quadrature degree for a matrix operator.
Definition: LocalOperator.hpp:130
auto order(F const &f) -> decltype(&F::operator(), f.order())
polynomial order of functions
Definition: Order.hpp:11
auto makeLocalOperator(PreGridFunctionOperator< Args... > op, GridView const &gridView)
Definition: GridFunctionOperator.hpp:258
The main implementation of an operator to be used in a Assembler.
Definition: LocalOperator.hpp:28
Wrapper class for element and geometry.
Definition: ContextGeometry.hpp:43
std::integral_constant< bool, B > bool_t
A wrapper for bool types.
Definition: Logical.hpp:12
Definition: ConvectionDiffusionOperator.hpp:23