Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
vof.h
Go to the documentation of this file.
1/** @file vof.h
2 */
3/**
4# Volume-Of-Fluid advection
5
6We want to approximate the solution of the advection equations
7\f[
8\partial_tc_i + \mathbf{u}_f\cdot\nabla c_i = 0
9\f]
10where \f$c_i\f$ are volume fraction fields describing sharp interfaces.
11
12This can be done using a conservative, non-diffusive geometric VOF
13scheme.
14
15We also add the option to transport diffusive tracers (aka "VOF
16concentrations") confined to one side of the interface i.e. solve the
17equations
18\f[
19\partial_tt_{i,j} + \nabla\cdot(\mathbf{u}_ft_{i,j}) = 0
20\f]
21with \f$t_{i,j} = c_if_j\f$ (or \f$t_{i,j} = (1 - c_i)f_j\f$) and \f$f_j\f$ is a
22volumetric tracer concentration (see [Lopez-Herrera et al, 2015](#lopez2015)).
23
24The list of tracers associated with the volume fraction is stored in
25the *tracers* attribute. For each tracer, the "side" of the interface
26(i.e. either \f$c\f$ or \f$1 - c\f$) is controlled by the *inverse*
27attribute). */
28
30 scalar * tracers, c;
31 bool inverse;
32}
33
34/**
35We will need basic functions for volume fraction computations. */
36
37#include "fractions.h"
38
39/**
40The list of volume fraction fields `interfaces`, will be provided by
41the user.
42
43The face velocity field `uf` will be defined by a solver as well
44as the timestep. */
45
46extern scalar * interfaces;
47extern vector uf;
48extern double dt;
49
50/**
51The gradient of a VOF-concentration `t` is computed using a standard
52three-point scheme if we are far enough from the interface (as
53controlled by *cmin*), otherwise a two-point scheme biased away from
54the interface is used. */
55
56for (int _d = 0; _d < dimension; _d++)
58{
59 static const double cmin = 0.5;
60 double cl = c[-1], cc = c[], cr = c[1];
61 if (t.inverse)
62 cl = 1. - cl, cc = 1. - cc, cr = 1. - cr;
63 if (cc >= cmin && t.gradient != zero) {
64 if (cr >= cmin) {
65 if (cl >= cmin) {
66 if (t.gradient)
67 return t.gradient (t[-1]/cl, t[]/cc, t[1]/cr)/Delta;
68 else
69 return (t[1]/cr - t[-1]/cl)/(2.*Delta);
70 }
71 else
72 return (t[1]/cr - t[]/cc)/Delta;
73 }
74 else if (cl >= cmin)
75 return (t[]/cc - t[-1]/cl)/Delta;
76 }
77 return 0.;
78}
79
80/**
81On trees, VOF concentrations need to be refined properly i.e. using
82volume-fraction-weighted linear interpolation of the concentration. */
83
84#if TREE
86{
87 scalar f = s.c;
88 if (cm[] == 0. || (!s.inverse && f[] <= 0.) || (s.inverse && f[] >= 1.))
89 for (int _c = 0; _c < 4; _c++) /* foreach_child */
90 s[] = 0.;
91 else {
92 coord g;
93 for (int _d = 0; _d < dimension; _d++)
95 double sc = s.inverse ? s[]/(1. - f[]) : s[]/f[], cmc = 4.*cm[];
96 for (int _c = 0; _c < 4; _c++) /* foreach_child */ {
97 s[] = sc;
98 for (int _d = 0; _d < dimension; _d++)
99 s[] += child.x*g.x*cm[-child.x]/cmc;
100 s[] *= s.inverse ? 1. - f[] : f[];
101 }
102 }
103}
104
105/**
106On trees, we need to setup the appropriate prolongation and
107refinement functions for the volume fraction fields. */
108
109/** @brief Event: defaults (i = 0) */
110void event_defaults (void)
111{
112 for (int _c = 0; _c < 1; _c++) /* scalar in interfaces */ {
113 c.refine = fraction_refine;
115 scalar * tracers = c.tracers;
116 for (int _t = 0; _t < 1; _t++) /* scalar in tracers */ {
120 t.c = c;
121 }
122 }
123}
124#endif // TREE
125
126/**
127Boundary conditions for VOF-advected tracers usually depend on
128boundary conditions for the VOF field. */
129
130/** @brief Event: defaults (i = 0) */
131void event_defaults (void)
132{
133 for (int _c = 0; _c < 1; _c++) /* scalar in interfaces */ {
134 scalar * tracers = c.tracers;
135 for (int _t = 0; _t < 1; _t++) /* scalar in tracers */
136 t.depends = list_add (t.depends, c);
137 }
138}
139
140/**
141We need to make sure that the CFL is smaller than 0.5 to ensure
142stability of the VOF scheme. */
143
144/** @brief Event: stability (i++) */
145void event_stability (void) {
146 if (CFL > 0.5)
147 CFL = 0.5;
148}
149
150/**
151## One-dimensional advection
152
153The simplest way to implement a multi-dimensional VOF advection scheme
154is to use dimension-splitting i.e. advect the field along each
155dimension successively using a one-dimensional scheme.
156
157We implement the one-dimensional scheme along the x-dimension and use
158the [for (int _d = 0; _d < dimension; _d++)](/Basilisk C#foreach_dimension) operator to
159automatically derive the corresponding functions along the other
160dimensions. */
161
162for (int _d = 0; _d < dimension; _d++)
163static void sweep_x (scalar c, scalar cc, scalar * tcl)
164{
165 vector n[];
167 double cfl = 0.;
168
169 /**
170 If we are also transporting tracers associated with \f$c\f$, we need to
171 compute their gradient i.e. \f$\partial_xf_j = \partial_x(t_j/c)\f$ or
172 \f$\partial_xf_j = \partial_x(t_j/(1 - c))\f$ (for higher-order
173 upwinding) and we need to store the computed fluxes. We first
174 allocate the corresponding lists. */
175
176 scalar * tracers = c.tracers, * gfl = NULL, * tfluxl = NULL;
177 if (tracers) {
178 for (int _t = 0; _t < 1; _t++) /* scalar in tracers */ {
179 scalar gf = {0} /* new scalar */, flux = {0} /* new scalar */;
180 gfl = list_append (gfl, gf);
182 }
183
184 /**
185 The gradient is computed using the "interface-biased" scheme above. */
186
187 for (int _i = 0; _i < _N; _i++) /* foreach */ {
188 scalar t, gf;
189 for (t,gf in tracers,gfl)
191 }
192 }
193
194 /**
195 We reconstruct the interface normal \f$\mathbf{n}\f$ and the intercept
196 \f$\alpha\f$ for each cell. Then we go through each (vertical) face of
197 the grid. */
198
200 for (int _i = 0; _i < _N; _i++) /* foreach_face */) {
201
202 /**
203 To compute the volume fraction flux, we check the sign of the velocity
204 component normal to the face and compute the index `i` of the
205 corresponding *upwind* cell (either 0 or -1). */
206
207 double un = uf.x[]*dt/(Delta*fm.x[] + SEPS), s = sign(un);
208 int i = -(s + 1.)/2.;
209
210 /**
211 We also check that we are not violating the CFL condition. */
212
213#if EMBED
214 if (cs[] >= 1.)
215#endif
216 if (un*fm.x[]*s/(cm[] + SEPS) > cfl)
217 cfl = un*fm.x[]*s/(cm[] + SEPS);
218
219 /**
220 If we assume that `un` is negative i.e. `s` is -1 and `i` is 0, the
221 volume fraction flux through the face of the cell is given by the dark
222 area in the figure below. The corresponding volume fraction can be
223 computed using the `rectangle_fraction()` function.
224
225 ![Volume fraction flux](figures/flux.svg)
226
227 When the upwind cell is entirely full or empty we can avoid this
228 computation. */
229
230 double cf = (c[i] <= 0. || c[i] >= 1.) ? c[i] :
231 rectangle_fraction ((coord){-s*n.x[i], n.y[i], n.z[i]}, alpha[i],
232 (coord){-0.5, -0.5, -0.5},
233 (coord){s*un - 0.5, 0.5, 0.5});
234
235 /**
236 Once we have the upwind volume fraction *cf*, the volume fraction
237 flux through the face is simply: */
238
239 flux[] = cf*uf.x[];
240
241 /**
242 If we are transporting tracers, we compute their flux using the
243 upwind volume fraction *cf* and a tracer value upwinded using the
244 Bell--Collela--Glaz scheme and the gradient computed above. */
245
246 scalar t, gf, tflux;
247 for (t,gf,tflux in tracers,gfl,tfluxl) {
248 double cf1 = cf, ci = c[i];
249 if (t.inverse)
250 cf1 = 1. - cf1, ci = 1. - ci;
251 if (ci > 1e-10) {
252 double ff = t[i]/ci + s*min(1., 1. - s*un)*gf[i]*Delta/2.;
253 tflux[] = ff*cf1*uf.x[];
254 }
255 else
256 tflux[] = 0.;
257 }
258 }
259 delete (gfl); free (gfl);
260
261 /**
262 We warn the user if the CFL condition has been violated. */
263
264 if (cfl > 0.5 + 1e-6)
265 fprintf (ferr,
266 "src/vof.h:%d: warning: CFL must be <= 0.5 for VOF (cfl - 0.5 = %g)\n",
267 LINENO, cfl - 0.5), fflush (ferr);
268
269 /**
270 Once we have computed the fluxes on all faces, we can update the
271 volume fraction field according to the one-dimensional advection
272 equation
273 \f[
274 \partial_tc = -\nabla_x\cdot(\mathbf{u}_f c) + c\nabla_x\cdot\mathbf{u}_f
275 \f]
276 The first term is computed using the fluxes. The second term -- which is
277 non-zero for the one-dimensional velocity field -- is approximated using
278 a centered volume fraction field `cc` which will be defined below.
279
280 For tracers, the corresponding one-dimensional update is
281 \f[
282 \partial_tt_j = -\nabla_x\cdot(\mathbf{u}_f t_j) + t_j\nabla_x\cdot\mathbf{u}_f
283 \f]
284 The compressive second-term can be removed by defining the
285 NO_1D_COMPRESSION macro. */
286
287#if !EMBED
288 for (int _i = 0; _i < _N; _i++) /* foreach */ {
289 c[] += dt*(flux[] - flux[1] + cc[]*(uf.x[1] - uf.x[]))/(cm[]*Delta);
290#if NO_1D_COMPRESSION
291 scalar t, tflux;
292 for (t, tflux in tracers, tfluxl)
293 t[] += dt*(tflux[] - tflux[1])/(cm[]*Delta);
294#else // !NO_1D_COMPRESSION
295 scalar t, tc, tflux;
296 for (t, tc, tflux in tracers, tcl, tfluxl)
297 t[] += dt*(tflux[] - tflux[1] + tc[]*(uf.x[1] - uf.x[]))/(cm[]*Delta);
298#endif // !NO_1D_COMPRESSION
299 }
300#else // EMBED
301 /**
302 When dealing with embedded boundaries, we simply ignore the fraction
303 occupied by the solid. This is a simple approximation which has the
304 advantage of ensuring boundedness of the volume fraction and
305 conservation of the total tracer mass (if it is computed also
306 ignoring the volume occupied by the solid in partial cells). */
307
308 for (int _i = 0; _i < _N; _i++) /* foreach */
309 if (cs[] > 0.) {
310 c[] += dt*cs[]*(flux[] - flux[1] + cc[]*(uf.x[1] - uf.x[]))/(cm[]*Delta);
311#if NO_1D_COMPRESSION
312 for (t, tflux in tracers, tfluxl)
313 t[] += dt*cs[]*(tflux[] - tflux[1])/(cm[]*Delta);
314#else // !NO_1D_COMPRESSION
315 scalar t, tc, tflux;
316 for (t, tc, tflux in tracers, tcl, tfluxl)
317 t[] += dt*cs[]*(tflux[] - tflux[1] + tc[]*(uf.x[1] - uf.x[]))/(cm[]*Delta);
318#endif // !NO_1D_COMPRESSION
319 }
320#endif // EMBED
321
322 delete (tfluxl); free (tfluxl);
323}
324
325/**
326## Multi-dimensional advection
327
328The multi-dimensional advection is performed by the event below. */
329
331{
332 for (int _c = 0; _c < 1; _c++) /* scalar in interfaces */ {
333
334 /**
335 We first define the volume fraction field used to compute the
336 divergent term in the one-dimensional advection equation above. We
337 follow [Weymouth & Yue, 2010](/src/references.bib#weymouth2010) and use a
338 step function which guarantees exact mass conservation for the
339 multi-dimensional advection scheme (provided the advection velocity
340 field is exactly non-divergent). */
341
342 scalar cc[], * tcl = NULL, * tracers = c.tracers;
343 for (int _t = 0; _t < 1; _t++) /* scalar in tracers */ {
344#if !NO_1D_COMPRESSION
345 scalar tc = {0} /* new scalar */;
346 tcl = list_append (tcl, tc);
347#endif // !NO_1D_COMPRESSION
348#if TREE
349 if (t.refine != vof_concentration_refine) {
353 t.c = c;
354 }
355#endif // TREE
356 }
357 for (int _i = 0; _i < _N; _i++) /* foreach */ {
358 cc[] = (c[] > 0.5);
359#if !NO_1D_COMPRESSION
360 scalar t, tc;
361 for (t, tc in tracers, tcl) {
362 if (t.inverse)
363 tc[] = c[] < 0.5 ? t[]/(1. - c[]) : 0.;
364 else
365 tc[] = c[] > 0.5 ? t[]/c[] : 0.;
366 }
367#endif // !NO_1D_COMPRESSION
368 }
369
370 /**
371 We then apply the one-dimensional advection scheme along each
372 dimension. To try to minimise phase errors, we alternate dimensions
373 according to the parity of the iteration index `i`. */
374
376 int d = 0;
377 for (int _d = 0; _d < dimension; _d++)
378 sweep[d++] = sweep_x;
379 for (d = 0; d < dimension; d++)
380 sweep[(i + d) % dimension] (c, cc, tcl);
381 delete (tcl), free (tcl);
382 }
383}
384
385/** @brief Event: vof (i++) */
386void event_vof (void)
388
389/**
390## References
391
392~~~bib
393@Article{lopez2015,
394 title = {A VOF numerical study on the electrokinetic effects in the
395 breakup of electrified jets},
396 author = {J. M. Lopez-Herrera and A. M. Ganan-Calvo and S. Popinet and
397 M. A. Herrada},
398 journal = {International Journal of Multiphase Flows},
399 pages = {14-22},
400 volume = {71},
401 year = {2015},
402 doi = {doi.org/10.1016/j.ijmultiphaseflow.2014.12.005},
403 url = {https://gerris.dalembert.upmc.fr/papers/lopez2015.pdf}
404}
405~~~
406*/
vector g[]
We store the combined pressure gradient and acceleration field in g*.
Definition all-mach.h:65
vector un[]
Definition atmosphere.h:5
int min
Definition balance.h:192
#define dimension
Definition bitree.h:3
int x
Definition common.h:76
double zero(double s0, double s1, double s2)
Definition common.h:116
const scalar cm[]
Definition common.h:397
scalar * list_add(scalar *list, scalar s)
Definition common.h:208
static int sign(number x)
Definition common.h:13
scalar * list_append(scalar *list, scalar s)
Definition common.h:188
const vector fm[]
Definition common.h:396
#define SEPS
Definition common.h:402
scalar f[]
The primary fields are:
Definition two-phase.h:56
Point point
Definition conserving.h:86
else return n
Definition curvature.h:101
scalar s
Definition embed-tree.h:56
scalar cs[]
The volume and area fractions are stored in these fields.
Definition embed.h:21
double d[2]
Definition embed.h:383
void fraction_refine(Point point, scalar c)
Definition fractions.h:42
double rectangle_fraction(coord n, double alpha, coord a, coord b)
VOF algorithms require the computation of volume fractions on (rectangular) parts of the initial squa...
Definition geometry.h:269
static void restriction_volume_average(Point point, scalar s)
void set_prolongation(scalar s, void(*prolongation)(Point, scalar))
void set_restriction(scalar s, void(*restriction)(Point, scalar))
def _i
Definition stencils.h:405
Definition linear.h:21
Definition common.h:78
scalar x
Definition common.h:47
double CFL
Definition utils.h:10
double dt
These come from the multilayer solver.
scalar * interfaces
We will need basic functions for volume fraction computations.
Definition two-phase.h:56
void event_vof(void) vof_advection(interfaces
Event: vof (i++)
void event_stability(void)
We need to make sure that the CFL is smaller than 0.5 to ensure stability of the VOF scheme.
Definition vof.h:145
void vof_advection(scalar *interfaces, int i)
Definition vof.h:330
double cr
Definition vof.h:60
bool inverse
Definition vof.h:31
free(gfl)
void i
Definition vof.h:387
scalar flux[]
Definition vof.h:166
scalar * tracers
If we are also transporting tracers associated with , we need to compute their gradient i....
Definition vof.h:176
scalar * gfl
Definition vof.h:176
double cfl
Definition vof.h:167
double cc
Definition vof.h:60
scalar c
Definition vof.h:57
attribute
Definition vof.h:29
scalar scalar t
Definition vof.h:58
src vof fflush(ferr)
void event_defaults(void)
On trees, we need to setup the appropriate prolongation and refinement functions for the volume fract...
Definition vof.h:110
src vof LINENO
Definition vof.h:267
reconstruction(c, n, alpha)
We reconstruct the interface normal and the intercept for each cell.
static void vof_concentration_refine(Point point, scalar s)
On trees, VOF concentrations need to be refined properly i.e.
Definition vof.h:85
scalar alpha[]
Definition vof.h:166
double cl
Definition vof.h:60
scalar * tfluxl
Definition vof.h:176
vector uf
We allocate the (face) velocity field.
Definition advection.h:29
scalar scalar * tcl
Definition vof.h:164