|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
#include "tag.h"Go to the source code of this file.
Macros | |
| #define | EPS 1e-10 |
| This function does a fast detection of cases which may correspond to two interfaces being close to one another. | |
Functions | |
| static scalar | fclone (int i) |
| This function returns a renamed clone of the default volume fraction field f, with i its index. | |
| static bool | tracer_is_close (Point point, scalar c) |
| static bool | bubbles_are_close (Point point, scalar c, scalar b, int *b1, int *b2) |
| This is similar to the function above, but now takes into account whether the two interfaces belong to different bubbles (identified by the tag field b). | |
| trace void | no_coalescence () |
| void | event_vof (void) |
| Event: vof (i++) | |
| void | event_tracer_advection (void) |
| After VOF advection, the default volume fraction field f is updated as the sum of all VOF tracer fields. | |
| void | event_cleanup (void) |
| We need to free the list of interfaces if it has been dynamically allocated. | |
| void | event_defaults (void) |
| Event: defaults (i = 0) | |
Variables | |
| int | number_of_interfaces = 0 |
This function does a fast detection of cases which may correspond to two interfaces being close to one another.
Definition at line 75 of file no-coalescence.h.
This is similar to the function above, but now takes into account whether the two interfaces belong to different bubbles (identified by the tag field b).
If they do then the indices of the two bubbles are returned in b1 and b2.
Definition at line 97 of file no-coalescence.h.
References b, b1, b2, c, dimension, EPS, i, j, k, and x.
Referenced by no_coalescence().
We need to free the list of interfaces if it has been dynamically allocated.
Event: cleanup (i = end)
Definition at line 444 of file no-coalescence.h.
References f, free(), scalar::i, i, interfaces, and x.
Event: defaults (i = 0)
Boundary conditions for VOF-advected tracers usually depend on boundary conditions for the VOF field.
Event: defaults (i = 0)
Event: defaults (i = 0)
We set the default values.
The EHD force density, \(\mathbf{f}_e\), can be computed as the divergence of the Maxwell stress tensor \(\mathbf{M}\),
\[ M_{ij} = \varepsilon (E_i E_j - \frac{E^2}{2}\delta_{ij}) \]
where \(E_i\) is the \(i\)-component of the electric field, \(\mathbf{E}=-\nabla \phi\) and \(\delta_{ij}\) is the Kronecker delta.
We need to add the corresponding acceleration to the Navier–Stokes solver.
If the acceleration vector a (defined by the Navier–Stokes solver) is constant, we make it variable.
Event: defaults (i = 0)
Event: defaults (i = 0 )
On trees we need to ensure conservation of the tracer when refining/coarsening.
Event: defaults (i = 0)
it is an acceleration. If necessary, we allocate a new vector field to store it.
Event: defaults (i = 0)
Event: defaults (i = 0)
\[ \begin{aligned} 0 & = - \sum_k \nabla \cdot [\theta_H (hu)_k^{n + 1} + (1 - \theta_H) (hu)^n_k] \\ \frac{(hu)^{n + 1}_k - (hu)_k^n}{\Delta t} & = - \Delta tgh^{n + 1 / 2}_k (\theta_H \nabla \eta_r^{n + 1} + (1 - \theta_H) \nabla \eta_r^n) \end{aligned} \]
where \(\eta_r\) is the equivalent free-surface height (i.e. pressure) applied on the rigid lid.
Event: defaults (i = 0)
The \(w_k\) and \(\phi_k\) scalar fields are allocated and the \(w_k\) are added to the list of advected tracers.
Event: defaults (i = 0)
Event: defaults (i = 0)
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)
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.
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 464 of file no-coalescence.h.
References f, fclone(), scalar::i, interfaces, list_append(), list_copy(), list_len(), number_of_interfaces, swap, t, and x.
Event: vof (i++)
We apply the no coalescence function just before VOF advection.
Definition at line 419 of file no-coalescence.h.
References no_coalescence().
This function returns a renamed clone of the default volume fraction field f, with i its index.
When two interfaces defined by the same VOF tracer are close enough (of the order of the cell size), they automatically merge. This is one of the strength and weakness of the VOF method.
In some cases, it may be desirable to avoid coalescence entirely, for example in the case of foams, emulsions, bubble clouds etc...
A simple way to do this is to use a different VOF tracer for each bubble/droplet. When one wants to simulate more than a few bubbles, this can of course become very expensive (both in CPU and memory).
This simple idea can be improved by noting that it would be sufficient to use different VOF tracers only for bubbles which are "too close" to one another. Determining the minimum number of VOF tracers required, for a given arrangement of bubbles, is clearly a variant of the graph coloring problem. In two dimensions, the famous four color theorem states that a maximum of four VOF tracers are required. Note however that finding this optimal coloring can be very difficult (NP-complete). The important point here is that one can expect that even a non-optimal number of VOF tracers will be much smaller than the number of bubbles.
This file is typically combined with the two-phase solver. From the user point-of-view, the only thing to be aware of is that the default f volume fraction field is not transported using VOF anymore, but is the sum of all VOF tracers. Individual VOF tracers are named f0, f1, f2, ... and are stored in the interfaces list (defined by vof.h). They should be used in particular to display the actual interfaces, as displaying interfaces using f will result in coalescence artefacts.
We will need to "tag" individual bubbles. EPS is the threshold used for tagging. threshold_volume is the minimal volum to which we ignore coalesence event
Definition at line 60 of file no-coalescence.h.
References c, f, free(), i, s, scalar_clone, strdup(), and x.
Referenced by event_defaults(), and no_coalescence().
We first make a quick test to check which VOF tracers may correspond to bubbles which are too close to one another. This is essentially an optimisation which avoids calling the relatively expensive tag()* function if it is obviously not necessary.
For each VOF tracer which may be too close, we first tag the corresponding bubbles.
The next step is to build the array tc of the bubbles which are indeed too close to one another.
We need to know which tracer fields are neigboring each bubble. If the tracer of index j is neighboring the bubble of index i (in tc), then adj[i*nvar + j] is set to true. Besides we add two to nvar since 3 scalar fields might be added simulataneously in one step.
Since we are updating adj, we cannot use openMP.
We check whether bubble b[] neighbors cells containing another tracer.
If this is the first bubble we need to replace, then the only existing VOF interface is f. We create a new interface, add it to the list and remove f from the list of interfaces (since it is not advected by VOF anymore).
Array rep will contain the index of the replacement VOF tracer for each bubble.
The indices of the pair of neighboring bubbles are stored in p[i] and p[i+1]. We need to replace only one of the two bubbles. We choose to replace the bubble with the smallest number of neighboring tracers.
We check if the tags are already modified, in which case we do not want to modify them again.
We look for a replacement VOF tracer which is not already neighboring the bubble.
If we didn't find any, we create a new one.
Refresh the adj list for all pairs of bubbles which contain an index adj to p[i] or to tag_not_modified.
We perform the replacement for each bubble (which is too close).
Finally, we free the arrays and lists.
Definition at line 164 of file no-coalescence.h.
References _i, array_append(), array_free(), array_new(), assert, b, b1, b2, bubbles_are_close(), c, datasize, EPS, f, fclone(), free(), scalar::i, i, int, interfaces, j, k, l, list_append(), list_copy(), list_len(), p, point, reset, s, swap, t, tag, tracer_is_close(), and x.
Referenced by event_vof().
Definition at line 77 of file no-coalescence.h.
References c, dimension, EPS, i, j, k, and x.
Referenced by no_coalescence().
| int number_of_interfaces = 0 |
This event is called before init. It is useful in case of restoring through a dump. For this one needs to add number_of_interfaces as a command line argument. While restoring one needs to pass the argument value, which is equal to the number of VOF tracers in the interfaces list.
Definition at line 461 of file no-coalescence.h.
Referenced by event_defaults().