|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
Go to the source code of this file.
Data Structures | |
| struct | HDiffusion |
Functions | |
| void | event_defaults (void) |
| Event: defaults (i = 0) | |
| void | event_vof (void) |
| Event: vof (i++) | |
| static void | h_relax (scalar *al, scalar *bl, int l, void *data) |
| static double | h_residual (scalar *al, scalar *bl, scalar *resl, void *data) |
| void | event_tracer_diffusion (void) |
| Event: tracer_diffusion (i++) | |
Variables | |
| attribute | |
| scalar | phi1 |
| scalar | phi2 |
| scalar * | stracers |
| static scalar * | phi_tracers = NULL |
Event: defaults (i = 0)
On trees we need to ensure conservation of the tracer when refining/coarsening.
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 49 of file henry.h.
References refine_embed_linear(), refine_linear(), restriction_volume_average(), s, set_prolongation(), set_restriction(), and x.
Event: tracer_diffusion (i++)
The advected concentration is computed from \(\phi_1\) and \(\phi_2\) as
\[ c = \phi_1 + \phi_2 \]
and these fields are then discarded.
The diffusion equation for \(c\) is then solved using the multigrid solver and the residual and relaxation functions defined above.
Definition at line 178 of file henry.h.
References _i, beta, c, cm, D, dt, f, fm, free(), h_relax(), h_residual(), mg_solve(), phi1, phi2, phi_tracers, q, restriction, vector::x, and x.
Event: vof (i++)
\(\phi_1\) and \(\phi_2\) are computed from \(c\) as
\[ \phi_1 = c \frac{\alpha f}{\alpha f + (1 - f)} \]
\[ \phi_2 = c \frac{1 - f}{\alpha f + (1 - f)} \]
Definition at line 74 of file henry.h.
References _i, a, c, f, list_append(), phi1, phi2, phi_tracers, scalar_clone, and x.
| attribute |
We consider the transport and diffusion of a tracer \(c\) with different solubilities in the two-phases described by two-phase.h.
The diffusion coefficients in the two phases are \(D_1\) and \(D_2\) and the jump in concentration at the interface is given by
\[ c_1 = \alpha c_2 \]
The advection/diffusion equation for \(c\) can then be written
\[ \partial_t c + \nabla\cdot(\mathbf{u} c) = \nabla\cdot\left(D\nabla c - D \frac{c(\alpha - 1)} {\alpha f + (1 - f)}\nabla f\right) \]
with \(f\) the volume fraction and
\[ D = \frac{D_1 D_2}{D_2 f + D_1 (1 - f)} \]
the harmonic mean of the diffusion coefficients (see Haroun et al., 2010).
The diffusion coefficients and solubility are attributes of each tracer.
The stracers list of soluble tracers must be defined by the calling code.
| scalar phi1 |
Definition at line 36 of file henry.h.
Referenced by event_tracer_diffusion(), and event_vof().
| scalar phi2 |
Definition at line 36 of file henry.h.
Referenced by event_tracer_diffusion(), and event_vof().
To avoid numerical diffusion through the interface we use the VOF tracer transport scheme for the temporary fields \(\phi_1\) and \(\phi_2\), see section 3.2 of Farsoiya et al., 2021.
Definition at line 71 of file henry.h.
Referenced by event_tracer_diffusion(), and event_vof().
|
extern |