|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
Go to the source code of this file.
Functions | |
| void | event_defaults (void) |
| Here we set the default boundary conditions for the streamfunction. | |
| void | event_velocity (void) |
| At every timestep we update the streamfunction field \(\psi\) by solving a Poisson equation with the updated vorticity field \(\omega\) (which has just been advected/diffused). | |
Variables | |
| scalar | omega [] |
| We allocate the vorticity field \(\omega\), the streamfunction field \(\psi\) and a structure to store the statistics on the convergence of the Poisson solver. | |
| scalar | psi [] |
| mgstats | mgpsi |
| scalar * | tracers = {omega} |
Here we set the gradient functions for each tracer (as defined in the user-provided tracers list). | |
Here we set the default boundary conditions for the streamfunction.
Event: defaults (i = 0)
The default convention in Basilisk is no-flow through the boundaries of the domain, i.e. they are a streamline i.e. \(\psi=\)constant on the boundary. We set the default value for the CFL (the default in utils.h is 0.5). This is done once at the beginning of the simulation.
Event: defaults (i = 0)
The default display.
Definition at line 59 of file stream.h.
References CFL, and display().
At every timestep we update the streamfunction field \(\psi\) by solving a Poisson equation with the updated vorticity field \(\omega\) (which has just been advected/diffused).
Event: velocity (i++)
Using the new streamfunction, we can then update the components of the velocity field. Since they are defined on faces we need to average the gradients, which gives the discrete expression below (for the horizontal velocity component). The expression for the vertical velocity component is obtained by automatic permutation of the indices but requires a change of sign: this is done through the pseudo-vector f.
| mgstats mgpsi |
Definition at line 40 of file stream.h.
Referenced by event_velocity().
| scalar omega[] |
We allocate the vorticity field \(\omega\), the streamfunction field \(\psi\) and a structure to store the statistics on the convergence of the Poisson solver.
In two dimensions the incompressible, constant-density Navier–Stokes equations can be written
\[ \partial_t\omega + \mathbf{u}\cdot\nabla\omega = \nu\nabla^2\omega \]
\[ \nabla^2\psi = \omega \]
with \(\nu\) the viscosity coefficient. The vorticity \(\omega\) and streamfunction \(\psi\) are defined as
\[ \omega = \partial_x u_y - \partial_y u_x \]
\[ u_x = - \partial_y\psi \]
\[ u_y = \partial_x\psi \]
The equation for the vorticity is an advection–diffusion equation which can be solved using the flux–based advection scheme in advection.h. The equation for the streamfunction is a Poisson equation. The fields advected by the advection solver are listed in tracers.
Definition at line 39 of file stream.h.
Referenced by axistream(), event_tracer_advection(), event_velocity(), gotm_step(), and vorticity().
| scalar psi[] |
Definition at line 39 of file stream.h.
Referenced by axistream(), event_velocity(), if(), momentum(), and rectangular_source().
Here we set the gradient functions for each tracer (as defined in the user-provided tracers list).
See Al Saud et al., 2018 and Popinet & Zaleski, 1999 for details.
The surface tension field \(\sigma\) will be associated to each levelset tracer. This is done easily by adding the following field attributes.