|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
#include "bcg.h"Go to the source code of this file.
Functions | |
| void | event_defaults (void) |
| On adaptive meshes, tracers need to use linear interpolation (rather than the default bilinear interpolation) to ensure conservation when refining cells. | |
| void | event_tracer_advection (void) |
| The integration is performed using the Bell-Collela-Glaz scheme. | |
| event | vof (i++, last) |
| VOF advection happens here. | |
| event | tracer_diffusion (i++, last) |
| Diffusion can be added by overloading this hook. | |
Variables | |
| scalar * | tracers |
Here we set the gradient functions for each tracer (as defined in the user-provided tracers list). | |
| vector | uf |
| We allocate the (face) velocity field. | |
| double | dt |
| These come from the multilayer solver. | |
On adaptive meshes, tracers need to use linear interpolation (rather than the default bilinear interpolation) to ensure conservation when refining cells.
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.
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.
We set limiting for q, h.
As well as default values.
On trees, we ensure that limiting is also applied to prolongation.
We setup the default display.
By default, all the layers have the same relative thickness.
We overload the default 'advance' and 'update' functions of the predictor-corrector scheme and setup the prolongation and restriction methods on trees.
On trees we make sure that slope-limiting is also used for refinement and prolongation. The prolongation/restriction functions for \(\eta\) are set and they depend on boundary conditions on \(z_b\) and \(h\).
We setup the default display.
Definition at line 28 of file tracer.h.
References refine_embed_linear(), refine_linear(), restriction_volume_average(), s, set_prolongation(), set_restriction(), and x.
|
extern |
These come from the multilayer solver.
This function approximates implicitly the EHD equation set given by the electric Poisson equation and the ohmic conduction term of the charge density conservation (the advection term is computed elsewhere using a tracer advection scheme),
\[ \nabla \cdot( \epsilon \nabla \phi^{n+1}) = -\rho_e^{n+1} \quad \mathrm{and} \quad (\rho_e^{n+1}-\rho_e^n) = \Delta t\nabla \cdot (K\nabla \phi^{n+1}) \]
where \(\rho_e\) is the charge density, and \(\phi\) is the electric potential. \(K\) and \(\epsilon\) are the conductivity and permittivity respectively.
Substituting the Poisson equation into the conservation equation gives the following time-implicit scheme,
\[ \nabla \cdot [(K \, \Delta t + \epsilon) \nabla \phi^{n+1}] = -\rho_e^{n} \]
We need the Poisson solver and the timestep dt.
The run() function below implements a generic time loop which executes events until termination.
The timestep dt can be accessed as a global variable.
Definition at line 18 of file predictor-corrector.h.
Referenced by event_tracer_advection().
|
extern |
Here we set the gradient functions for each tracer (as defined in the user-provided tracers list).
This event integrates advection equations of the form
\[ \partial_tf_i+\mathbf{u_f}\cdot\nabla f_i=0 \]
where \(\mathbf{u_f}\) is the velocity field and \(f_i\) are a list of passive tracers.
The tracers list is defined elsewhere (typically by the user), the vector field uf and the timestep dt are defined by a solver.
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.
Definition at line 60 of file hydro.h.
Referenced by event_cleanup(), event_defaults(), event_pressure(), event_tracer_advection(), for(), if(), and vof_advection().
|
extern |
We allocate the (face) velocity field.
We wish to solve the advection equations
\[ \partial_tf_i+\mathbf{u}\cdot\nabla f_i=0 \]
where \(\mathbf{u}\) is the velocity field and \(f_i\) are a list of passive tracers. This can be done with a flux-based advection scheme such as the 2nd-order, unsplit, upwind scheme of Bell-Collela-Glaz, 1989.
The main time loop is defined in [run.h](). A stable timestep needs to respect the CFL condition. For compatibility with the other solvers, we allocate it as uf and define an alias. The gradient function is used to set the type of slope-limiting required. The default is to not use any limiting (i.e. a purely centered slope estimation).
Definition at line 29 of file advection.h.
Referenced by event_tracer_advection().