Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
swirl.h File Reference
#include "tracer.h"
#include "diffusion.h"
Include dependency graph for swirl.h:

Go to the source code of this file.

Functions

void event_defaults (void)
 The azimuthal velocity is zero on the axis of symmetry ( \(y=0\)).
 
void event_acceleration (void)
 The equation for \(u_y\) is solved by the centered Navier–Stokes solver combined with the axisymmetric metric, but the acceleration term \(w^2/y\) is missing.
 
void event_tracer_diffusion (void)
 The advection of \(w\) is done by the tracer solver, but we need to add diffusion.
 

Variables

scalar w []
 
scalartracers = {w}
 Here we set the gradient functions for each tracer (as defined in the user-provided tracers list).
 

Function Documentation

◆ event_acceleration()

void event_acceleration ( void  )

The equation for \(u_y\) is solved by the centered Navier–Stokes solver combined with the axisymmetric metric, but the acceleration term \(w^2/y\) is missing.

We add it here, taking care of the division by zero on the axis, and averaging \(w\) from cell center to cell face.

Event: acceleration (i++)

Definition at line 70 of file swirl.h.

References _i, a, av, sq(), w, x, vector::y, and y.

Here is the call graph for this function:

◆ event_defaults()

void event_defaults ( void  )

The azimuthal velocity is zero on the axis of symmetry ( \(y=0\)).

Event: defaults (i = 0)

We will need to add the acceleration term \(w^2/y\) in the evolution equation for \(u_y\). If the acceleration field is not allocated yet, we do so.

Event: defaults (i = 0)

The default density field is set to unity (times the metric).

We reset the multigrid parameters to their default values.

If the viscosity is non-zero, we need to allocate the face-centered viscosity field.

We also initialize the list of tracers to be advected with the VOF function \(f\) (or its complementary function).

We set limiting.

On trees, we ensure that limiting is also applied to prolongation and refinement.

We add the interface and the density to the default display.

We switch to a pure minmod limiter by default for increased robustness.

With the MUSCL scheme we use the CFL depends on the dimension of the problem.

On trees we need to replace the default bilinear refinement/prolongation with linear so that reconstructed values also use slope limiting.

The restriction/refine attributes of the charge density are those of a tracer otherwise the conservation is not guaranteed.

By default the permittivity is unity and other quantities are zero.

The (velocity) CFL is limited by the unsplit advection scheme, so is dependent on the dimension. The (gravity wave) CFL is set to 1/2 (if not already set by the user).

The gradient and prolongation/restriction functions are set for all tracer fields.

We setup the default display.

Boundary conditions

By default we set a zero Neumann boundary condition for all the components except if the bottom is an axis of symmetry.

We use (strict) minmod slope limiting for all components.

We reset the multigrid parameters to their default values.

The pressures are never dumped.

The default density field is set to unity (times the metric and the solid factors).

On trees, refinement of the face-centered velocity field needs to preserve the divergence-free condition.

When using embedded boundaries, the restriction and prolongation operators need to take the boundary into account.

We set the dimensions of the velocity field.

On trees, the refinement and restriction functions above rely on the volume fraction field f being refined/restricted before the components of velocity. To ensure this, we move f to the front of the field list (all).

We then set the refinement and restriction functions for the components of the velocity field. The boundary conditions on \(\mathbf{u}\) now depend on those on \(f\).

The default display.

Definition at line 54 of file swirl.h.

References _i, a, vector::x, and x.

◆ event_tracer_diffusion()

void event_tracer_diffusion ( void  )

The advection of \(w\) is done by the tracer solver, but we need to add diffusion.

Using the diffusion solver, we solve

\[ \theta\partial_tw = \nabla\cdot(D\nabla w) + \beta w \]

Identifying with the diffusion part of the equation for \(w\) above, we have

\[ \begin{aligned} \theta &= \rho y \\ D &= \mu y \\ \beta &= - \left(\rho u_y + \frac{\mu}{y} + \partial_y\mu \right) \end{aligned} \]

Note that the rho and mu fields (defined by the Navier–Stokes solver) already include the metric (i.e. are \(\rho y\) and \(\mu y\)), which explains the divisions by \(y\) in the code below.

Event: tracer_diffusion (i++)

Definition at line 97 of file swirl.h.

References _i, beta, diffusion(), dt, fm, mu, rho, theta, u, w, x, vector::y, and y.

Here is the call graph for this function:

Variable Documentation

◆ tracers

scalar * tracers = {w}

Here we set the gradient functions for each tracer (as defined in the user-provided tracers list).

Integral formulation for surface tension

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.

Definition at line 41 of file swirl.h.

◆ w

scalar w[]

Azimuthal velocity for axisymmetric flows

The centered Navier–Stokes solver can be combined with the axisymmetric metric but assumes zero azimuthal velocity ("swirl"). This file adds this azimuthal velocity component: the \(w\) field.

Assuming that \(x\) is the axial direction and \(y\) the radial direction, as in axi.h, the incompressible, variable-density and viscosity, axisymmetric Navier–Stokes equations (with swirl) can be written

\[ \partial_x u_x + \partial_y u_y + \frac{u_y}{y} = 0 \]

\[ \partial_t u_x + u_x \partial_x u_x + u_y \partial_y u_x = - \frac{1}{\rho} \partial_x p + \frac{1}{\rho y} \nabla \cdot (2 \mu y \nabla \mathbf{D}_x) \]

\[ \partial_t u_y + u_x \partial_x u_y + u_y \partial_y u_y - {\color{blue}\frac{w^2}{y}} = - \frac{1}{\rho} \partial_y p + \frac{1}{\rho y} \left( \nabla \cdot (2 \mu y \nabla \mathbf{D}_y) - 2 \mu \frac{u_y}{y} \right) \]

\[ {\color{blue} \partial_t w + u_x \partial_x w + u_y \partial_y w + \frac{u_y w}{y} = \frac{1}{\rho y} \left[ \nabla \cdot (\mu y \nabla w) - w \left( \frac{\mu}{y} + \partial_y \mu \right) \right] } \]

where the terms in blue are the missing "swirl" terms. We will thus need to solve an advection-diffusion equation for \(w\).

Definition at line 41 of file swirl.h.

Referenced by event_acceleration(), and event_tracer_diffusion().