deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
flow_function.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2007 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_flow_function_h
16#define dealii_flow_function_h
17
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/base/mutex.h>
23#include <deal.II/base/point.h>
24
26
27namespace Functions
28{
46 template <int dim>
47 class FlowFunction : public Function<dim>
48 {
49 public:
54
58 virtual ~FlowFunction() override = default;
59
64 void
65 pressure_adjustment(double p);
66
72 virtual void
73 vector_values(const std::vector<Point<dim>> &points,
74 std::vector<std::vector<double>> &values) const override = 0;
80 virtual void
82 const std::vector<Point<dim>> &points,
83 std::vector<std::vector<Tensor<1, dim>>> &gradients) const override = 0;
92 virtual void
93 vector_laplacians(const std::vector<Point<dim>> &points,
94 std::vector<std::vector<double>> &values) const = 0;
95
96 virtual void
97 vector_value(const Point<dim> &points,
98 Vector<double> &value) const override;
99 virtual double
100 value(const Point<dim> &points,
101 const unsigned int component) const override;
102 virtual void
103 vector_value_list(const std::vector<Point<dim>> &points,
104 std::vector<Vector<double>> &values) const override;
105 virtual void
107 const std::vector<Point<dim>> &points,
108 std::vector<std::vector<Tensor<1, dim>>> &gradients) const override;
112 virtual void
113 vector_laplacian_list(const std::vector<Point<dim>> &points,
114 std::vector<Vector<double>> &values) const override;
115
119 virtual std::size_t
120 memory_consumption() const override;
121
122 protected:
127
128 private:
133
137 mutable std::vector<std::vector<double>> aux_values;
138
142 mutable std::vector<std::vector<Tensor<1, dim>>> aux_gradients;
143 };
144
152 template <int dim>
153 class PoisseuilleFlow : public FlowFunction<dim>
154 {
155 public:
160 PoisseuilleFlow(const double r, const double Re);
161
162 virtual ~PoisseuilleFlow() override = default;
163
164 virtual void
165 vector_values(const std::vector<Point<dim>> &points,
166 std::vector<std::vector<double>> &values) const override;
167 virtual void
169 const std::vector<Point<dim>> &points,
170 std::vector<std::vector<Tensor<1, dim>>> &gradients) const override;
171 virtual void
172 vector_laplacians(const std::vector<Point<dim>> &points,
173 std::vector<std::vector<double>> &values) const override;
174
175 private:
176 const double inv_sqr_radius;
177 const double Reynolds;
178 };
179
180
193 template <int dim>
194 class StokesCosine : public FlowFunction<dim>
195 {
196 public:
201 StokesCosine(const double viscosity = 1., const double reaction = 0.);
205 void
206 set_parameters(const double viscosity, const double reaction);
207
208 virtual ~StokesCosine() override = default;
209
210 virtual void
211 vector_values(const std::vector<Point<dim>> &points,
212 std::vector<std::vector<double>> &values) const override;
213 virtual void
215 const std::vector<Point<dim>> &points,
216 std::vector<std::vector<Tensor<1, dim>>> &gradients) const override;
217 virtual void
218 vector_laplacians(const std::vector<Point<dim>> &points,
219 std::vector<std::vector<double>> &values) const override;
220
221 private:
223 double viscosity;
225 double reaction;
226 };
227
228
246 {
247 public:
250
251 virtual void
252 vector_values(const std::vector<Point<2>> &points,
253 std::vector<std::vector<double>> &values) const override;
254 virtual void
256 const std::vector<Point<2>> &points,
257 std::vector<std::vector<Tensor<1, 2>>> &gradients) const override;
258 virtual void
259 vector_laplacians(const std::vector<Point<2>> &points,
260 std::vector<std::vector<double>> &values) const override;
261
262 private:
264 double
265 Psi(double phi) const;
267 double
268 Psi_1(double phi) const;
270 double
271 Psi_2(double phi) const;
273 double
274 Psi_3(double phi) const;
276 double
277 Psi_4(double phi) const;
279 const double omega;
282 static const double lambda;
284 const double coslo;
286 const double lp;
288 const double lm;
289 };
290
298 class Kovasznay : public FlowFunction<2>
299 {
300 public:
308 Kovasznay(const double Re, bool Stokes = false);
309
310 virtual ~Kovasznay() override = default;
311
312 virtual void
313 vector_values(const std::vector<Point<2>> &points,
314 std::vector<std::vector<double>> &values) const override;
315 virtual void
317 const std::vector<Point<2>> &points,
318 std::vector<std::vector<Tensor<1, 2>>> &gradients) const override;
319 virtual void
320 vector_laplacians(const std::vector<Point<2>> &points,
321 std::vector<std::vector<double>> &values) const override;
322
324 double
325 lambda() const;
326
327 private:
328 const double Reynolds;
329 double lbda;
330 double p_average;
331 const bool stokes;
332 };
333
334} // namespace Functions
335
337
338#endif
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
std::vector< std::vector< Tensor< 1, dim > > > aux_gradients
void pressure_adjustment(double p)
virtual void vector_gradient_list(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const override
virtual void vector_laplacian_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const override
virtual void vector_laplacians(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const =0
virtual void vector_gradients(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const override=0
virtual void vector_values(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const override=0
virtual std::size_t memory_consumption() const override
virtual ~FlowFunction() override=default
virtual void vector_value(const Point< dim > &points, Vector< double > &value) const override
std::vector< std::vector< double > > aux_values
virtual double value(const Point< dim > &points, const unsigned int component) const override
virtual ~Kovasznay() override=default
virtual void vector_values(const std::vector< Point< 2 > > &points, std::vector< std::vector< double > > &values) const override
virtual void vector_gradients(const std::vector< Point< 2 > > &points, std::vector< std::vector< Tensor< 1, 2 > > > &gradients) const override
virtual void vector_laplacians(const std::vector< Point< 2 > > &points, std::vector< std::vector< double > > &values) const override
double lambda() const
The value of lambda.
virtual void vector_gradients(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const override
virtual void vector_values(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const override
virtual void vector_laplacians(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const override
virtual ~PoisseuilleFlow() override=default
double reaction
The reaction parameter.
virtual void vector_values(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const override
virtual ~StokesCosine() override=default
double viscosity
The viscosity.
void set_parameters(const double viscosity, const double reaction)
virtual void vector_gradients(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const override
virtual void vector_laplacians(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const override
virtual void vector_gradients(const std::vector< Point< 2 > > &points, std::vector< std::vector< Tensor< 1, 2 > > > &gradients) const override
virtual void vector_values(const std::vector< Point< 2 > > &points, std::vector< std::vector< double > > &values) const override
double Psi_1(double phi) const
The derivative of Psi()
double Psi(double phi) const
The auxiliary function Psi.
const double lp
Auxiliary variable 1+lambda.
const double lm
Auxiliary variable 1-lambda.
const double omega
The angle of the reentrant corner, set to 3*pi/2.
virtual void vector_laplacians(const std::vector< Point< 2 > > &points, std::vector< std::vector< double > > &values) const override
double Psi_3(double phi) const
The 3rd derivative of Psi()
StokesLSingularity()
Constructor setting up some data.
double Psi_4(double phi) const
The 4th derivative of Psi()
double Psi_2(double phi) const
The 2nd derivative of Psi()
const double coslo
Cosine of lambda times omega.
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499