Reference documentation for deal.II version 8.5.1
flow_function.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2007 - 2015 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__flow_function_h
17 #define dealii__flow_function_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/function.h>
22 #include <deal.II/base/point.h>
23 #include <deal.II/base/thread_management.h>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 namespace Functions
28 {
47  template <int dim>
48  class FlowFunction : public Function<dim>
49  {
50  public:
54  FlowFunction();
55 
59  virtual ~FlowFunction();
60 
65  void pressure_adjustment(double p);
66 
72  virtual void vector_values (const std::vector<Point<dim> > &points,
73  std::vector<std::vector<double> > &values) const = 0;
79  virtual void vector_gradients (const std::vector<Point<dim> > &points,
80  std::vector<std::vector<Tensor<1,dim> > > &gradients) const = 0;
89  virtual void vector_laplacians (const std::vector<Point<dim> > &points,
90  std::vector<std::vector<double> > &values) const = 0;
91 
92  virtual void vector_value (const Point<dim> &points, Vector<double> &value) const;
93  virtual double value (const Point<dim> &points, const unsigned int component) const;
94  virtual void vector_value_list (const std::vector<Point<dim> > &points,
95  std::vector<Vector<double> > &values) const;
96  virtual void vector_gradient_list (const std::vector<Point<dim> > &points,
97  std::vector<std::vector<Tensor<1,dim> > > &gradients) const;
101  virtual void vector_laplacian_list (const std::vector<Point<dim> > &points,
102  std::vector<Vector<double> > &values) const;
103 
104  std::size_t memory_consumption () const;
105 
106  protected:
111 
112  private:
113 
118 
122  mutable std::vector<std::vector<double> > aux_values;
123 
127  mutable std::vector<std::vector<Tensor<1,dim> > > aux_gradients;
128  };
129 
138  template <int dim>
139  class PoisseuilleFlow : public FlowFunction<dim>
140  {
141  public:
146  PoisseuilleFlow<dim> (const double r,
147  const double Re);
148  virtual ~PoisseuilleFlow();
149 
150  virtual void vector_values (const std::vector<Point<dim> > &points,
151  std::vector<std::vector<double> > &values) const;
152  virtual void vector_gradients (const std::vector<Point<dim> > &points,
153  std::vector<std::vector<Tensor<1,dim> > > &gradients) const;
154  virtual void vector_laplacians (const std::vector<Point<dim> > &points,
155  std::vector<std::vector<double> > &values) const;
156 
157  private:
158  const double radius;
159  const double Reynolds;
160  };
161 
162 
176  template <int dim>
177  class StokesCosine :
178  public FlowFunction<dim>
179  {
180  public:
185  StokesCosine (const double viscosity = 1., const double reaction = 0.);
189  void set_parameters (const double viscosity, const double reaction);
190  virtual ~StokesCosine();
191 
192  virtual void vector_values (const std::vector<Point<dim> > &points,
193  std::vector<std::vector<double> > &values) const;
194  virtual void vector_gradients (const std::vector<Point<dim> > &points,
195  std::vector<std::vector<Tensor<1,dim> > > &gradients) const;
196  virtual void vector_laplacians (const std::vector<Point<dim> > &points,
197  std::vector<std::vector<double> > &values) const;
198 
199  private:
201  double viscosity;
203  double reaction;
204  };
205 
206 
216  {
217  public:
220 
221  virtual void vector_values (const std::vector<Point<2> > &points,
222  std::vector<std::vector<double> > &values) const;
223  virtual void vector_gradients (const std::vector<Point<2> > &points,
224  std::vector<std::vector<Tensor<1,2> > > &gradients) const;
225  virtual void vector_laplacians (const std::vector<Point<2> > &points,
226  std::vector<std::vector<double> > &values) const;
227  private:
229  double Psi(double phi) const;
231  double Psi_1(double phi) const;
233  double Psi_2(double phi) const;
235  double Psi_3(double phi) const;
237  double Psi_4(double phi) const;
239  const double omega;
241  static const double lambda;
243  const double coslo;
245  const double lp;
247  const double lm;
248  };
249 
258  class Kovasznay : public FlowFunction<2>
259  {
260  public:
268  Kovasznay (const double Re, bool Stokes = false);
269  virtual ~Kovasznay();
270 
271  virtual void vector_values (const std::vector<Point<2> > &points,
272  std::vector<std::vector<double> > &values) const;
273  virtual void vector_gradients (const std::vector<Point<2> > &points,
274  std::vector<std::vector<Tensor<1,2> > > &gradients) const;
275  virtual void vector_laplacians (const std::vector<Point<2> > &points,
276  std::vector<std::vector<double> > &values) const;
277 
279  double lambda () const;
280  private:
281  const double Reynolds;
282  double lbda;
283  double p_average;
284  const bool stokes;
285  };
286 
287 }
288 
289 DEAL_II_NAMESPACE_CLOSE
290 
291 #endif
std::vector< std::vector< double > > aux_values
const double lm
Auxiliary variable 1-lambda.
virtual double value(const Point< dim > &points, const unsigned int component) const
const double lp
Auxiliary variable 1+lambda.
virtual void vector_values(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const
static const double lambda
The exponent of the radius.
StokesCosine(const double viscosity=1., const double reaction=0.)
virtual void vector_gradients(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const =0
void set_parameters(const double viscosity, const double reaction)
double viscosity
The viscosity.
double Psi(double phi) const
The auxiliary function Psi.
virtual void vector_gradients(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const
virtual void vector_gradients(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const
double lambda() const
The value of lambda.
virtual void vector_values(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const
std::vector< std::vector< Tensor< 1, dim > > > aux_gradients
double reaction
The reaction parameter.
virtual void vector_value(const Point< dim > &points, Vector< double > &value) const
virtual void vector_laplacians(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const
const double omega
The angle of the reentrant corner.
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
double Psi_1(double phi) const
The derivative of Psi()
const double coslo
Cosine of lambda times omega.
double Psi_4(double phi) const
The 4th derivative of Psi()
virtual void vector_laplacians(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const =0
double Psi_3(double phi) const
The 3rd derivative of Psi()
StokesLSingularity()
Constructor setting up some data.
virtual void vector_values(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const =0
virtual void vector_laplacian_list(const std::vector< Point< dim > > &points, std::vector< Vector< double > > &values) const
virtual void vector_laplacians(const std::vector< Point< dim > > &points, std::vector< std::vector< double > > &values) const
Kovasznay(const double Re, bool Stokes=false)
double Psi_2(double phi) const
The 2nd derivative of Psi()
void pressure_adjustment(double p)