Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
advection.h
Go to the documentation of this file.
1/** @file advection.h
2 */
3/**
4# An advection solver
5
6We wish to solve the advection equations
7\f[
8\partial_tf_i+\mathbf{u}\cdot\nabla f_i=0
9\f]
10where \f$\mathbf{u}\f$ is the velocity field and \f$f_i\f$ are a list of
11passive tracers. This can be done with a flux-based advection scheme
12such as the 2nd-order, unsplit, upwind scheme of [Bell-Collela-Glaz,
131989](references.bib#bell89).
14
15The main time loop is defined in [run.h](). A stable timestep needs to
16respect the [CFL
17condition](http://en.wikipedia.org/wiki/Courant%E2%80%93Friedrichs%E2%80%93Lewy_condition). */
18
19#include "run.h"
20#include "timestep.h"
21
22/**
23We allocate the (face) velocity field. For compatibility with the
24other solvers, we allocate it as `uf` and define an alias. The
25`gradient` function is used to set the type of slope-limiting
26required. The default is to not use any limiting (i.e. a purely
27centered slope estimation). */
28
30#define u uf
31double (* gradient) (double, double, double) = NULL;
32
33/**
34Here we set the gradient functions for each tracer (as defined in the
35user-provided `tracers` list). */
36
37extern scalar * tracers;
38
39/** @brief Event: defaults (i = 0) */
40void event_defaults (void) {
41 for (int _f = 0; _f < 1; _f++) /* scalar in tracers */
42 f.gradient = gradient;
43}
44
45/**
46User initialisation happens here. */
47
48event init (i = 0);
49
50/**
51The timestep is set using the velocity field and the CFL
52criterion. The integration itself is performed in the events of
53[tracer.h](). */
54
55/** @brief Event: velocity (i++,last) */
56void event_velocity (void) {
57 dt = dtnext (timestep (u, DT));
58}
59
60#include "tracer.h"
#define u
Definition advection.h:30
vector uf[]
We allocate the (face) velocity field.
Definition advection.h:29
scalar * tracers
Here we set the gradient functions for each tracer (as defined in the user-provided tracers list).
Definition hydro.h:60
void event_defaults(void)
Event: defaults (i = 0)
Definition advection.h:40
event init(i=0)
User initialisation happens here.
double(* gradient)(double, double, double)
Definition advection.h:31
void event_velocity(void)
The timestep is set using the velocity field and the CFL criterion.
Definition advection.h:56
trace double timestep(void)
Definition atmosphere.h:37
scalar f[]
The primary fields are:
Definition two-phase.h:56
double dt
scalar int i
Definition embed.h:74
double dtnext(double dt)
Definition events.h:276
double DT
Definition utils.h:10