Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
tracer.h
Go to the documentation of this file.
1/** @file tracer.h
2 */
3/**
4# Tracer advection event
5
6This event integrates advection equations of the form
7\f[
8\partial_tf_i+\mathbf{u_f}\cdot\nabla f_i=0
9\f]
10where \f$\mathbf{u_f}\f$ is the velocity field and \f$f_i\f$ are a list of
11passive tracers.
12
13The `tracers` list is defined elsewhere (typically by the user), the
14vector field `uf` and the timestep `dt` are defined by a
15solver. */
16
17extern scalar * tracers;
18extern vector uf;
19extern double dt;
20
21/**
22On adaptive meshes, tracers need to use linear interpolation (rather
23than the default bilinear interpolation) to ensure conservation when
24refining cells. */
25
26#if TREE
27/** @brief Event: defaults (i = 0) */
28void event_defaults (void) {
29 for (int _s = 0; _s < 1; _s++) /* scalar in tracers */ {
30#if EMBED
31 s.refine = refine_embed_linear;
33#else
34 s.refine = refine_linear;
35#endif
37 }
38}
39#endif
40
41/**
42The integration is performed using the Bell-Collela-Glaz scheme. */
43
44#include "bcg.h"
45
46/** @brief Event: tracer_advection (i++,last) */
49}
50
51/**
52VOF advection happens here. */
53
54event vof (i++, last);
55
56/**
57Diffusion can be added by overloading this hook. */
58
59event tracer_diffusion (i++,last);
trace void advection(scalar *tracers, vector u, double dt, scalar *src=NULL)
The function below uses the tracer_fluxes function to integrate the advection equation,...
Definition bcg.h:71
int x
Definition common.h:76
static void refine_embed_linear(Point point, scalar s)
Definition embed-tree.h:270
scalar s
Definition embed-tree.h:56
scalar int i
Definition embed.h:74
static void restriction_volume_average(Point point, scalar s)
static void refine_linear(Point point, scalar s)
void set_prolongation(scalar s, void(*prolongation)(Point, scalar))
void set_restriction(scalar s, void(*restriction)(Point, scalar))
double dt
These come from the multilayer solver.
event tracer_diffusion(i++, last)
Diffusion can be added by overloading this hook.
scalar * tracers
Here we set the gradient functions for each tracer (as defined in the user-provided tracers list).
Definition hydro.h:60
event vof(i++, last)
VOF advection happens here.
void event_defaults(void)
On adaptive meshes, tracers need to use linear interpolation (rather than the default bilinear interp...
Definition tracer.h:28
void event_tracer_advection(void)
The integration is performed using the Bell-Collela-Glaz scheme.
Definition tracer.h:47
vector uf
We allocate the (face) velocity field.
Definition advection.h:29