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