Reference documentation for deal.II version 9.1.0-pre
flow_function.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2007 - 2017 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/flow_function.h>
17 #include <deal.II/base/point.h>
18 #include <deal.II/base/tensor.h>
19 
20 #include <deal.II/lac/vector.h>
21 
22 #include <cmath>
23 
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 
28 namespace Functions
29 {
30  template <int dim>
32  : Function<dim>(dim + 1)
33  , mean_pressure(0)
34  , aux_values(dim + 1)
35  , aux_gradients(dim + 1)
36  {}
37 
38 
39 
40  template <int dim>
41  void
43  {
44  mean_pressure = p;
45  }
46 
47 
48  template <int dim>
49  void
51  const std::vector<Point<dim>> &points,
52  std::vector<Vector<double>> & values) const
53  {
54  const unsigned int n_points = points.size();
55  Assert(values.size() == n_points,
56  ExcDimensionMismatch(values.size(), n_points));
57 
58  // guard access to the aux_*
59  // variables in multithread mode
61 
62  for (unsigned int d = 0; d < dim + 1; ++d)
63  aux_values[d].resize(n_points);
64  vector_values(points, aux_values);
65 
66  for (unsigned int k = 0; k < n_points; ++k)
67  {
68  Assert(values[k].size() == dim + 1,
69  ExcDimensionMismatch(values[k].size(), dim + 1));
70  for (unsigned int d = 0; d < dim + 1; ++d)
71  values[k](d) = aux_values[d][k];
72  }
73  }
74 
75 
76  template <int dim>
77  void
79  Vector<double> & value) const
80  {
81  Assert(value.size() == dim + 1,
82  ExcDimensionMismatch(value.size(), dim + 1));
83 
84  const unsigned int n_points = 1;
85  std::vector<Point<dim>> points(1);
86  points[0] = point;
87 
88  // guard access to the aux_*
89  // variables in multithread mode
91 
92  for (unsigned int d = 0; d < dim + 1; ++d)
93  aux_values[d].resize(n_points);
94  vector_values(points, aux_values);
95 
96  for (unsigned int d = 0; d < dim + 1; ++d)
97  value(d) = aux_values[d][0];
98  }
99 
100 
101  template <int dim>
102  double
104  const unsigned int comp) const
105  {
106  Assert(comp < dim + 1, ExcIndexRange(comp, 0, dim + 1));
107  const unsigned int n_points = 1;
108  std::vector<Point<dim>> points(1);
109  points[0] = point;
110 
111  // guard access to the aux_*
112  // variables in multithread mode
114 
115  for (unsigned int d = 0; d < dim + 1; ++d)
116  aux_values[d].resize(n_points);
117  vector_values(points, aux_values);
118 
119  return aux_values[comp][0];
120  }
121 
122 
123  template <int dim>
124  void
126  const std::vector<Point<dim>> & points,
127  std::vector<std::vector<Tensor<1, dim>>> &values) const
128  {
129  const unsigned int n_points = points.size();
130  Assert(values.size() == n_points,
131  ExcDimensionMismatch(values.size(), n_points));
132 
133  // guard access to the aux_*
134  // variables in multithread mode
136 
137  for (unsigned int d = 0; d < dim + 1; ++d)
138  aux_gradients[d].resize(n_points);
140 
141  for (unsigned int k = 0; k < n_points; ++k)
142  {
143  Assert(values[k].size() == dim + 1,
144  ExcDimensionMismatch(values[k].size(), dim + 1));
145  for (unsigned int d = 0; d < dim + 1; ++d)
146  values[k][d] = aux_gradients[d][k];
147  }
148  }
149 
150 
151  template <int dim>
152  void
154  const std::vector<Point<dim>> &points,
155  std::vector<Vector<double>> & values) const
156  {
157  const unsigned int n_points = points.size();
158  Assert(values.size() == n_points,
159  ExcDimensionMismatch(values.size(), n_points));
160 
161  // guard access to the aux_*
162  // variables in multithread mode
164 
165  for (unsigned int d = 0; d < dim + 1; ++d)
166  aux_values[d].resize(n_points);
167  vector_laplacians(points, aux_values);
168 
169  for (unsigned int k = 0; k < n_points; ++k)
170  {
171  Assert(values[k].size() == dim + 1,
172  ExcDimensionMismatch(values[k].size(), dim + 1));
173  for (unsigned int d = 0; d < dim + 1; ++d)
174  values[k](d) = aux_values[d][k];
175  }
176  }
177 
178 
179  template <int dim>
180  std::size_t
182  {
183  Assert(false, ExcNotImplemented());
184  return 0;
185  }
186 
187 
188  //----------------------------------------------------------------------//
189 
190  template <int dim>
191  PoisseuilleFlow<dim>::PoisseuilleFlow(const double r, const double Re)
192  : radius(r)
193  , Reynolds(Re)
194  {
195  Assert(Reynolds != 0., ExcMessage("Reynolds number cannot be zero"));
196  }
197 
198 
199 
200  template <int dim>
201  void
203  const std::vector<Point<dim>> & points,
204  std::vector<std::vector<double>> &values) const
205  {
206  unsigned int n = points.size();
207  double stretch = 1. / radius;
208 
209  Assert(values.size() == dim + 1,
210  ExcDimensionMismatch(values.size(), dim + 1));
211  for (unsigned int d = 0; d < dim + 1; ++d)
212  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
213 
214  for (unsigned int k = 0; k < n; ++k)
215  {
216  const Point<dim> &p = points[k];
217  // First, compute the
218  // square of the distance to
219  // the x-axis divided by the
220  // radius.
221  double r2 = 0;
222  for (unsigned int d = 1; d < dim; ++d)
223  r2 += p(d) * p(d) * stretch * stretch;
224 
225  // x-velocity
226  values[0][k] = 1. - r2;
227  // other velocities
228  for (unsigned int d = 1; d < dim; ++d)
229  values[d][k] = 0.;
230  // pressure
231  values[dim][k] = -2 * (dim - 1) * stretch * stretch * p(0) / Reynolds +
232  this->mean_pressure;
233  }
234  }
235 
236 
237 
238  template <int dim>
239  void
241  const std::vector<Point<dim>> & points,
242  std::vector<std::vector<Tensor<1, dim>>> &values) const
243  {
244  unsigned int n = points.size();
245  double stretch = 1. / radius;
246 
247  Assert(values.size() == dim + 1,
248  ExcDimensionMismatch(values.size(), dim + 1));
249  for (unsigned int d = 0; d < dim + 1; ++d)
250  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
251 
252  for (unsigned int k = 0; k < n; ++k)
253  {
254  const Point<dim> &p = points[k];
255  // x-velocity
256  values[0][k][0] = 0.;
257  for (unsigned int d = 1; d < dim; ++d)
258  values[0][k][d] = -2. * p(d) * stretch * stretch;
259  // other velocities
260  for (unsigned int d = 1; d < dim; ++d)
261  values[d][k] = 0.;
262  // pressure
263  values[dim][k][0] = -2 * (dim - 1) * stretch * stretch / Reynolds;
264  for (unsigned int d = 1; d < dim; ++d)
265  values[dim][k][d] = 0.;
266  }
267  }
268 
269 
270 
271  template <int dim>
272  void
274  const std::vector<Point<dim>> & points,
275  std::vector<std::vector<double>> &values) const
276  {
277  unsigned int n = points.size();
278  (void)n;
279  Assert(values.size() == dim + 1,
280  ExcDimensionMismatch(values.size(), dim + 1));
281  for (unsigned int d = 0; d < dim + 1; ++d)
282  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
283 
284  for (unsigned int d = 0; d < values.size(); ++d)
285  for (unsigned int k = 0; k < values[d].size(); ++k)
286  values[d][k] = 0.;
287  }
288 
289  //----------------------------------------------------------------------//
290 
291  template <int dim>
292  StokesCosine<dim>::StokesCosine(const double nu, const double r)
293  : viscosity(nu)
294  , reaction(r)
295  {}
296 
297 
298 
299  template <int dim>
300  void
301  StokesCosine<dim>::set_parameters(const double nu, const double r)
302  {
303  viscosity = nu;
304  reaction = r;
305  }
306 
307 
308  template <int dim>
309  void
311  const std::vector<Point<dim>> & points,
312  std::vector<std::vector<double>> &values) const
313  {
314  unsigned int n = points.size();
315 
316  Assert(values.size() == dim + 1,
317  ExcDimensionMismatch(values.size(), dim + 1));
318  for (unsigned int d = 0; d < dim + 1; ++d)
319  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
320 
321  for (unsigned int k = 0; k < n; ++k)
322  {
323  const Point<dim> &p = points[k];
324  const double x = numbers::PI / 2. * p(0);
325  const double y = numbers::PI / 2. * p(1);
326  const double cx = cos(x);
327  const double cy = cos(y);
328  const double sx = sin(x);
329  const double sy = sin(y);
330 
331  if (dim == 2)
332  {
333  values[0][k] = cx * cx * cy * sy;
334  values[1][k] = -cx * sx * cy * cy;
335  values[2][k] = cx * sx * cy * sy + this->mean_pressure;
336  }
337  else if (dim == 3)
338  {
339  const double z = numbers::PI / 2. * p(2);
340  const double cz = cos(z);
341  const double sz = sin(z);
342 
343  values[0][k] = cx * cx * cy * sy * cz * sz;
344  values[1][k] = cx * sx * cy * cy * cz * sz;
345  values[2][k] = -2. * cx * sx * cy * sy * cz * cz;
346  values[3][k] = cx * sx * cy * sy * cz * sz + this->mean_pressure;
347  }
348  else
349  {
350  Assert(false, ExcNotImplemented());
351  }
352  }
353  }
354 
355 
356 
357  template <int dim>
358  void
360  const std::vector<Point<dim>> & points,
361  std::vector<std::vector<Tensor<1, dim>>> &values) const
362  {
363  unsigned int n = points.size();
364 
365  Assert(values.size() == dim + 1,
366  ExcDimensionMismatch(values.size(), dim + 1));
367  for (unsigned int d = 0; d < dim + 1; ++d)
368  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
369 
370  for (unsigned int k = 0; k < n; ++k)
371  {
372  const Point<dim> &p = points[k];
373  const double x = numbers::PI / 2. * p(0);
374  const double y = numbers::PI / 2. * p(1);
375  const double c2x = cos(2 * x);
376  const double c2y = cos(2 * y);
377  const double s2x = sin(2 * x);
378  const double s2y = sin(2 * y);
379  const double cx2 = .5 + .5 * c2x; // cos^2 x
380  const double cy2 = .5 + .5 * c2y; // cos^2 y
381 
382  if (dim == 2)
383  {
384  values[0][k][0] = -.25 * numbers::PI * s2x * s2y;
385  values[0][k][1] = .5 * numbers::PI * cx2 * c2y;
386  values[1][k][0] = -.5 * numbers::PI * c2x * cy2;
387  values[1][k][1] = .25 * numbers::PI * s2x * s2y;
388  values[2][k][0] = .25 * numbers::PI * c2x * s2y;
389  values[2][k][1] = .25 * numbers::PI * s2x * c2y;
390  }
391  else if (dim == 3)
392  {
393  const double z = numbers::PI / 2. * p(2);
394  const double c2z = cos(2 * z);
395  const double s2z = sin(2 * z);
396  const double cz2 = .5 + .5 * c2z; // cos^2 z
397 
398  values[0][k][0] = -.125 * numbers::PI * s2x * s2y * s2z;
399  values[0][k][1] = .25 * numbers::PI * cx2 * c2y * s2z;
400  values[0][k][2] = .25 * numbers::PI * cx2 * s2y * c2z;
401 
402  values[1][k][0] = .25 * numbers::PI * c2x * cy2 * s2z;
403  values[1][k][1] = -.125 * numbers::PI * s2x * s2y * s2z;
404  values[1][k][2] = .25 * numbers::PI * s2x * cy2 * c2z;
405 
406  values[2][k][0] = -.5 * numbers::PI * c2x * s2y * cz2;
407  values[2][k][1] = -.5 * numbers::PI * s2x * c2y * cz2;
408  values[2][k][2] = .25 * numbers::PI * s2x * s2y * s2z;
409 
410  values[3][k][0] = .125 * numbers::PI * c2x * s2y * s2z;
411  values[3][k][1] = .125 * numbers::PI * s2x * c2y * s2z;
412  values[3][k][2] = .125 * numbers::PI * s2x * s2y * c2z;
413  }
414  else
415  {
416  Assert(false, ExcNotImplemented());
417  }
418  }
419  }
420 
421 
422 
423  template <int dim>
424  void
426  const std::vector<Point<dim>> & points,
427  std::vector<std::vector<double>> &values) const
428  {
429  unsigned int n = points.size();
430 
431  Assert(values.size() == dim + 1,
432  ExcDimensionMismatch(values.size(), dim + 1));
433  for (unsigned int d = 0; d < dim + 1; ++d)
434  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
435 
436  if (reaction != 0.)
437  {
438  vector_values(points, values);
439  for (unsigned int d = 0; d < dim; ++d)
440  for (unsigned int k = 0; k < values[d].size(); ++k)
441  values[d][k] *= -reaction;
442  }
443  else
444  {
445  for (unsigned int d = 0; d < dim; ++d)
446  for (unsigned int k = 0; k < values[d].size(); ++k)
447  values[d][k] = 0.;
448  }
449 
450 
451  for (unsigned int k = 0; k < n; ++k)
452  {
453  const Point<dim> &p = points[k];
454  const double x = numbers::PI / 2. * p(0);
455  const double y = numbers::PI / 2. * p(1);
456  const double c2x = cos(2 * x);
457  const double c2y = cos(2 * y);
458  const double s2x = sin(2 * x);
459  const double s2y = sin(2 * y);
460  const double pi2 = .25 * numbers::PI * numbers::PI;
461 
462  if (dim == 2)
463  {
464  values[0][k] += -viscosity * pi2 * (1. + 2. * c2x) * s2y -
465  numbers::PI / 4. * c2x * s2y;
466  values[1][k] += viscosity * pi2 * s2x * (1. + 2. * c2y) -
467  numbers::PI / 4. * s2x * c2y;
468  values[2][k] = 0.;
469  }
470  else if (dim == 3)
471  {
472  const double z = numbers::PI * p(2);
473  const double c2z = cos(2 * z);
474  const double s2z = sin(2 * z);
475 
476  values[0][k] +=
477  -.5 * viscosity * pi2 * (1. + 2. * c2x) * s2y * s2z -
478  numbers::PI / 8. * c2x * s2y * s2z;
479  values[1][k] += .5 * viscosity * pi2 * s2x * (1. + 2. * c2y) * s2z -
480  numbers::PI / 8. * s2x * c2y * s2z;
481  values[2][k] +=
482  -.5 * viscosity * pi2 * s2x * s2y * (1. + 2. * c2z) -
483  numbers::PI / 8. * s2x * s2y * c2z;
484  values[3][k] = 0.;
485  }
486  else
487  {
488  Assert(false, ExcNotImplemented());
489  }
490  }
491  }
492 
493 
494  //----------------------------------------------------------------------//
495 
496  const double StokesLSingularity::lambda = 0.54448373678246;
497 
499  : omega(3. / 2. * numbers::PI)
500  , coslo(cos(lambda * omega))
501  , lp(1. + lambda)
502  , lm(1. - lambda)
503  {}
504 
505 
506  inline double
507  StokesLSingularity::Psi(double phi) const
508  {
509  return coslo * (sin(lp * phi) / lp - sin(lm * phi) / lm) - cos(lp * phi) +
510  cos(lm * phi);
511  }
512 
513 
514  inline double
515  StokesLSingularity::Psi_1(double phi) const
516  {
517  return coslo * (cos(lp * phi) - cos(lm * phi)) + lp * sin(lp * phi) -
518  lm * sin(lm * phi);
519  }
520 
521 
522  inline double
523  StokesLSingularity::Psi_2(double phi) const
524  {
525  return coslo * (lm * sin(lm * phi) - lp * sin(lp * phi)) +
526  lp * lp * cos(lp * phi) - lm * lm * cos(lm * phi);
527  }
528 
529 
530  inline double
531  StokesLSingularity::Psi_3(double phi) const
532  {
533  return coslo * (lm * lm * cos(lm * phi) - lp * lp * cos(lp * phi)) +
534  lm * lm * lm * sin(lm * phi) - lp * lp * lp * sin(lp * phi);
535  }
536 
537 
538  inline double
539  StokesLSingularity::Psi_4(double phi) const
540  {
541  return coslo *
542  (lp * lp * lp * sin(lp * phi) - lm * lm * lm * sin(lm * phi)) +
543  lm * lm * lm * lm * cos(lm * phi) -
544  lp * lp * lp * lp * cos(lp * phi);
545  }
546 
547 
548  void
549  StokesLSingularity::vector_values(
550  const std::vector<Point<2>> & points,
551  std::vector<std::vector<double>> &values) const
552  {
553  unsigned int n = points.size();
554 
555  Assert(values.size() == 2 + 1, ExcDimensionMismatch(values.size(), 2 + 1));
556  for (unsigned int d = 0; d < 2 + 1; ++d)
557  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
558 
559  for (unsigned int k = 0; k < n; ++k)
560  {
561  const Point<2> &p = points[k];
562  const double x = p(0);
563  const double y = p(1);
564 
565  if ((x < 0) || (y < 0))
566  {
567  const double phi = std::atan2(y, -x) + numbers::PI;
568  const double r2 = x * x + y * y;
569  const double rl = pow(r2, lambda / 2.);
570  const double rl1 = pow(r2, lambda / 2. - .5);
571  values[0][k] =
572  rl * (lp * sin(phi) * Psi(phi) + cos(phi) * Psi_1(phi));
573  values[1][k] =
574  rl * (lp * cos(phi) * Psi(phi) - sin(phi) * Psi_1(phi));
575  values[2][k] = -rl1 * (lp * lp * Psi_1(phi) + Psi_3(phi)) / lm +
576  this->mean_pressure;
577  }
578  else
579  {
580  for (unsigned int d = 0; d < 3; ++d)
581  values[d][k] = 0.;
582  }
583  }
584  }
585 
586 
587 
588  void
589  StokesLSingularity::vector_gradients(
590  const std::vector<Point<2>> & points,
591  std::vector<std::vector<Tensor<1, 2>>> &values) const
592  {
593  unsigned int n = points.size();
594 
595  Assert(values.size() == 2 + 1, ExcDimensionMismatch(values.size(), 2 + 1));
596  for (unsigned int d = 0; d < 2 + 1; ++d)
597  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
598 
599  for (unsigned int k = 0; k < n; ++k)
600  {
601  const Point<2> &p = points[k];
602  const double x = p(0);
603  const double y = p(1);
604 
605  if ((x < 0) || (y < 0))
606  {
607  const double phi = std::atan2(y, -x) + numbers::PI;
608  const double r2 = x * x + y * y;
609  const double r = sqrt(r2);
610  const double rl = pow(r2, lambda / 2.);
611  const double rl1 = pow(r2, lambda / 2. - .5);
612  const double rl2 = pow(r2, lambda / 2. - 1.);
613  const double psi = Psi(phi);
614  const double psi1 = Psi_1(phi);
615  const double psi2 = Psi_2(phi);
616  const double cosp = cos(phi);
617  const double sinp = sin(phi);
618 
619  // Derivatives of u with respect to r, phi
620  const double udr = lambda * rl1 * (lp * sinp * psi + cosp * psi1);
621  const double udp = rl * (lp * cosp * psi + lp * sinp * psi1 -
622  sinp * psi1 + cosp * psi2);
623  // Derivatives of v with respect to r, phi
624  const double vdr = lambda * rl1 * (lp * cosp * psi - sinp * psi1);
625  const double vdp = rl * (lp * (cosp * psi1 - sinp * psi) -
626  cosp * psi1 - sinp * psi2);
627  // Derivatives of p with respect to r, phi
628  const double pdr =
629  -(lambda - 1.) * rl2 * (lp * lp * psi1 + Psi_3(phi)) / lm;
630  const double pdp = -rl1 * (lp * lp * psi2 + Psi_4(phi)) / lm;
631  values[0][k][0] = cosp * udr - sinp / r * udp;
632  values[0][k][1] = -sinp * udr - cosp / r * udp;
633  values[1][k][0] = cosp * vdr - sinp / r * vdp;
634  values[1][k][1] = -sinp * vdr - cosp / r * vdp;
635  values[2][k][0] = cosp * pdr - sinp / r * pdp;
636  values[2][k][1] = -sinp * pdr - cosp / r * pdp;
637  }
638  else
639  {
640  for (unsigned int d = 0; d < 3; ++d)
641  values[d][k] = 0.;
642  }
643  }
644  }
645 
646 
647 
648  void
649  StokesLSingularity::vector_laplacians(
650  const std::vector<Point<2>> & points,
651  std::vector<std::vector<double>> &values) const
652  {
653  unsigned int n = points.size();
654  (void)n;
655  Assert(values.size() == 2 + 1, ExcDimensionMismatch(values.size(), 2 + 1));
656  for (unsigned int d = 0; d < 2 + 1; ++d)
657  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
658 
659  for (unsigned int d = 0; d < values.size(); ++d)
660  for (unsigned int k = 0; k < values[d].size(); ++k)
661  values[d][k] = 0.;
662  }
663 
664 
665  //----------------------------------------------------------------------//
666 
667  Kovasznay::Kovasznay(double Re, bool stokes)
668  : Reynolds(Re)
669  , stokes(stokes)
670  {
671  long double r2 = Reynolds / 2.;
672  long double b = 4 * numbers::PI * numbers::PI;
673  long double l = -b / (r2 + std::sqrt(r2 * r2 + b));
674  lbda = l;
675  // mean pressure for a domain
676  // spreading from -.5 to 1.5 in
677  // x-direction
678  p_average = 1 / (8 * l) * (std::exp(3. * l) - std::exp(-l));
679  }
680 
681 
682 
683  void
684  Kovasznay::vector_values(const std::vector<Point<2>> & points,
685  std::vector<std::vector<double>> &values) const
686  {
687  unsigned int n = points.size();
688 
689  Assert(values.size() == 2 + 1, ExcDimensionMismatch(values.size(), 2 + 1));
690  for (unsigned int d = 0; d < 2 + 1; ++d)
691  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
692 
693  for (unsigned int k = 0; k < n; ++k)
694  {
695  const Point<2> &p = points[k];
696  const double x = p(0);
697  const double y = 2. * numbers::PI * p(1);
698  const double elx = std::exp(lbda * x);
699 
700  values[0][k] = 1. - elx * cos(y);
701  values[1][k] = .5 / numbers::PI * lbda * elx * sin(y);
702  values[2][k] = -.5 * elx * elx + p_average + this->mean_pressure;
703  }
704  }
705 
706 
707  void
708  Kovasznay::vector_gradients(
709  const std::vector<Point<2>> & points,
710  std::vector<std::vector<Tensor<1, 2>>> &gradients) const
711  {
712  unsigned int n = points.size();
713 
714  Assert(gradients.size() == 3, ExcDimensionMismatch(gradients.size(), 3));
715  Assert(gradients[0].size() == n,
716  ExcDimensionMismatch(gradients[0].size(), n));
717 
718  for (unsigned int i = 0; i < n; ++i)
719  {
720  const double x = points[i](0);
721  const double y = points[i](1);
722 
723  const double elx = std::exp(lbda * x);
724  const double cy = cos(2 * numbers::PI * y);
725  const double sy = sin(2 * numbers::PI * y);
726 
727  // u
728  gradients[0][i][0] = -lbda * elx * cy;
729  gradients[0][i][1] = 2. * numbers::PI * elx * sy;
730  gradients[1][i][0] = lbda * lbda / (2 * numbers::PI) * elx * sy;
731  gradients[1][i][1] = lbda * elx * cy;
732  // p
733  gradients[2][i][0] = -lbda * elx * elx;
734  gradients[2][i][1] = 0.;
735  }
736  }
737 
738 
739 
740  void
741  Kovasznay::vector_laplacians(const std::vector<Point<2>> & points,
742  std::vector<std::vector<double>> &values) const
743  {
744  unsigned int n = points.size();
745  Assert(values.size() == 2 + 1, ExcDimensionMismatch(values.size(), 2 + 1));
746  for (unsigned int d = 0; d < 2 + 1; ++d)
747  Assert(values[d].size() == n, ExcDimensionMismatch(values[d].size(), n));
748 
749  if (stokes)
750  {
751  const double zp = 2. * numbers::PI;
752  for (unsigned int k = 0; k < n; ++k)
753  {
754  const Point<2> &p = points[k];
755  const double x = p(0);
756  const double y = zp * p(1);
757  const double elx = std::exp(lbda * x);
758  const double u = 1. - elx * cos(y);
759  const double ux = -lbda * elx * cos(y);
760  const double uy = elx * zp * sin(y);
761  const double v = lbda / zp * elx * sin(y);
762  const double vx = lbda * lbda / zp * elx * sin(y);
763  const double vy = zp * lbda / zp * elx * cos(y);
764 
765  values[0][k] = u * ux + v * uy;
766  values[1][k] = u * vx + v * vy;
767  values[2][k] = 0.;
768  }
769  }
770  else
771  {
772  for (unsigned int d = 0; d < values.size(); ++d)
773  for (unsigned int k = 0; k < values[d].size(); ++k)
774  values[d][k] = 0.;
775  }
776  }
777 
778  double
780  {
781  return lbda;
782  }
783 
784 
785 
786  template class FlowFunction<2>;
787  template class FlowFunction<3>;
788  template class PoisseuilleFlow<2>;
789  template class PoisseuilleFlow<3>;
790  template class StokesCosine<2>;
791  template class StokesCosine<3>;
792 } // namespace Functions
793 
794 
795 
796 DEAL_II_NAMESPACE_CLOSE
std::vector< std::vector< double > > aux_values
virtual void vector_laplacians(const std::vector< Point< dim >> &points, std::vector< std::vector< double >> &values) const override
const double lm
Auxiliary variable 1-lambda.
virtual void vector_laplacians(const std::vector< Point< dim >> &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
const double lp
Auxiliary variable 1+lambda.
virtual void vector_laplacians(const std::vector< Point< dim >> &points, std::vector< std::vector< double >> &values) const =0
virtual void vector_value_list(const std::vector< Point< dim >> &points, std::vector< Vector< double >> &values) const override
static const double lambda
StokesCosine(const double viscosity=1., const double reaction=0.)
size_type size() const
virtual void vector_values(const std::vector< Point< dim >> &points, std::vector< std::vector< double >> &values) const override=0
virtual void vector_gradients(const std::vector< Point< dim >> &points, std::vector< std::vector< Tensor< 1, dim >>> &gradients) const override=0
void set_parameters(const double viscosity, const double reaction)
double viscosity
The viscosity.
static::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
double Psi_3(double phi) const
The 3rd derivative of Psi()
static const double PI
Definition: numbers.h:143
static::ExceptionBase & ExcMessage(std::string arg1)
virtual void vector_values(const std::vector< Point< dim >> &points, std::vector< std::vector< double >> &values) const override
#define Assert(cond, exc)
Definition: exceptions.h:1227
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
double reaction
The reaction parameter.
std::vector< std::vector< Tensor< 1, dim > > > aux_gradients
PoisseuilleFlow(const double r, const double Re)
virtual void vector_value(const Point< dim > &points, Vector< double > &value) const override
virtual void vector_values(const std::vector< Point< dim >> &points, std::vector< std::vector< double >> &values) const override
double Psi_2(double phi) const
The 2nd derivative of Psi()
virtual void vector_laplacian_list(const std::vector< Point< dim >> &points, std::vector< Vector< double >> &values) const override
const double coslo
Cosine of lambda times omega.
virtual double value(const Point< dim > &points, const unsigned int component) const override
StokesLSingularity()
Constructor setting up some data.
static::ExceptionBase & ExcNotImplemented()
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 override
Kovasznay(const double Re, bool Stokes=false)
void pressure_adjustment(double p)
double Psi_1(double phi) const
The derivative of Psi()
double Psi_4(double phi) const
The 4th derivative of Psi()