|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
#include "poisson.h"Go to the source code of this file.
Functions | |
| static void | relax_psi (scalar *al, scalar *bl, int l, void *data) |
| static double | residual_psi (scalar *al, scalar *bl, scalar *resl, void *data) |
| mgstats | axistream (vector u, scalar psi) |
| The function axistream() takes as input the \(u\) velocity field (with \(x=z\) and \(y=r\)) and returns the corresponding axisymmetric stream function. | |
The function axistream() takes as input the \(u\) velocity field (with \(x=z\) and \(y=r\)) and returns the corresponding axisymmetric stream function.
Definition at line 84 of file axistream.h.
References _i, mg_solve(), omega, psi, relax_psi(), residual_psi(), u, and y.
The axisymmetric or Stokes stream function \(\psi\) verifies
\[ u_r = -\frac{1}{r}\partial_z\psi \]
\[ u_z = +\frac{1}{r}\partial_r\psi \]
with \(u_r\) and \(u_z\) the radial and longitudinal velocity components of an incompressible axisymmetric flow.
Given the definition of the vorticity \(\omega\)
\[ \omega = \partial_zu_r - \partial_ru_z \]
\(\psi\) verifies the variable-coefficient Poisson equation
\[ \partial_z\left(\frac{1}{r}\partial_z\psi\right) + \partial_r\left(\frac{1}{r}\partial_r\psi\right) = - \omega \]
This is not well-conditioned due to the divergence of the coefficients at \(r = 0\) and can be rewritten instead as
\[ \frac{\partial^2\psi}{\partial z^2} + \frac{\partial^2\psi}{\partial r^2} - \frac{1}{r}\partial_r\psi = - \omega r \]
which is better behaved.
This equation can be inverted using the multigrid solver combined with the relaxation and residual functions defined below.
Definition at line 41 of file axistream.h.
References _i, a, b, l, sq(), and y.
Referenced by axistream().