|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
Go to the source code of this file.
Macros | |
| #define | CURVATURE 1 |
Functions | |
| void | event_defaults (void) |
| Surface tension is a source term in the right-hand-side of the evolution equation for the velocity of the centered Navier–Stokes solver i.e. | |
| void | event_stability (void) |
| Event: stability (i++) | |
| static double | distance_curvature (Point point, scalar d) |
| void | event_acceleration (void) |
| Event: acceleration (i++) | |
Variables | |
| scalar * | tracers |
Here we set the gradient functions for each tracer (as defined in the user-provided tracers list). | |
| attribute | |
| #define CURVATURE 1 |
This function computes the curvature from the levelset function d using:
\[ \kappa = \frac{d^2_xd_{yy} - 2d_xd_yd_{xy} + d^2_yd_{xx}}{|\nabla d|^3} \]
Definition at line 99 of file integral.h.
Definition at line 102 of file integral.h.
References cube(), d, dx, dy, sq(), and x.
Referenced by event_acceleration().
Event: acceleration (i++)
The calculation of the acceleration is done by this event, overloaded from its definition in the centered Navier–Stokes solver.
We check whether surface tension is associated with any levelset function d.
We allocate the surface tension stress tensor "manually", because the symmetries of the default tensor allocation in Basilisk are not what we want (this is messy).
We compute the diagonal components of the tensor.
We compute the off-diagonal components of the tensor.
Finally, we add the divergence of the surface tension tensor to the acceleration and the axisymmetric term if necessary.
The axisymmetric acceleration is computed using the volumetric surface tension formulation as
\[ \frac{1}{\rho}\sigma\kappa_\theta\mathbf{n}\delta_s \]
with the principal axisymmetric curvature given by
\[ \kappa_\theta = \frac{n^r}{r} \]
and using the approximation
\[ \mathbf{n}\delta_s\approx\mathbf{\nabla}f \]
with \(f\) the volume fraction.
Definition at line 126 of file integral.h.
References _i, a, alpha, av, d, dimension, distance_curvature(), f, fm, i, kappa, n, p, point, S, SEPS, sigma, sign(), sq(), vector::x, and x.
Surface tension is a source term in the right-hand-side of the evolution equation for the velocity of the centered Navier–Stokes solver i.e.
Event: defaults (i = 0)
it is an acceleration. If necessary, we allocate a new vector field to store it.
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.
Definition at line 26 of file integral.h.
Event: stability (i++)
We need to overload the stability event so that the CFL is taken into account (because we set stokes to true).
The surface tension scheme is time-explicit so the maximum timestep is the oscillation period of the smallest capillary wave.
\[ T = \sqrt{\frac{\rho_{m}\Delta_{min}^3}{\pi\sigma}} \]
with \(\rho_m=(\rho_1+\rho_2)/2.\) and \(\rho_1\), \(\rho_2\) the densities on either side of the interface.
We first compute the minimum and maximum values of \(\alpha/f_m = 1/\rho\), as well as \(\Delta_{min}\).
The maximum timestep is set using the sum of surface tension coefficients.
Definition at line 48 of file integral.h.
References _i, alpha, cube(), d, dmin, dt, dtmax, fm, HUGE, max, min, pi, sigma, vector::x, and x.
| attribute |
Definition at line 15 of file integral.h.
|
extern |
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.
Here we set the gradient functions for each tracer (as defined in the user-provided tracers list).
\(\partial_xf_j = \partial_x(t_j/c)\) or \(\partial_xf_j = \partial_x(t_j/(1 - c))\) (for higher-order upwinding) and we need to store the computed fluxes. We first allocate the corresponding lists.