Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
stream.h
Go to the documentation of this file.
1/** @file stream.h
2 */
3/**
4# A streamfunction--vorticity solver for the Navier--Stokes equations
5
6In two dimensions the incompressible, constant-density Navier--Stokes
7equations can be written
8\f[
9\partial_t\omega + \mathbf{u}\cdot\nabla\omega = \nu\nabla^2\omega
10\f]
11\f[
12\nabla^2\psi = \omega
13\f]
14with \f$\nu\f$ the viscosity coefficient. The vorticity \f$\omega\f$ and
15streamfunction \f$\psi\f$ are defined as
16\f[
17\omega = \partial_x u_y - \partial_y u_x
18\f]
19\f[
20u_x = - \partial_y\psi
21\f]
22\f[
23u_y = \partial_x\psi
24\f]
25The equation for the vorticity is an advection--diffusion equation
26which can be solved using the flux--based advection scheme in
27[`advection.h`](/src/advection.h). The equation for the streamfunction
28is a [Poisson equation](/src/poisson.h). */
29
30#include "advection.h"
31#include "poisson.h"
32
33/**
34We allocate the vorticity field \f$\omega\f$, the streamfunction field
35\f$\psi\f$ and a structure to store the statistics on the convergence of
36the Poisson solver. The fields advected by the [advection
37solver](/src/advection.h) are listed in `tracers`. */
38
42
43/**
44Here we set the default boundary conditions for the
45streamfunction. The default convention in Basilisk is no-flow through
46the boundaries of the domain, i.e. they are a streamline
47i.e. \f$\psi=\f$constant on the boundary. */
48
49/* BC: psi[right] = dirichlet(0) */
50/* BC: psi[left] = dirichlet(0) */
51/* BC: psi[top] = dirichlet(0) */
52/* BC: psi[bottom] = dirichlet(0) */
53
54/**
55We set the default value for the `CFL` (the default in `utils.h` is
560.5). This is done once at the beginning of the simulation. */
57
58/** @brief Event: defaults (i = 0) */
59void event_defaults (void) {
60 CFL = 0.8;
61
62 /**
63 The default display. */
64
65 display ("squares (color = 'omega', spread = -1);");
66}
67
68/**
69At every timestep we update the streamfunction field \f$\psi\f$ by solving
70a Poisson equation with the updated vorticity field \f$\omega\f$ (which
71has just been advected/diffused). */
72
73/** @brief Event: velocity (i++) */
74void event_velocity (void)
75{
77
78 /**
79 Using the new streamfunction, we can then update the components of
80 the velocity field. Since they are defined on faces we need to
81 average the gradients, which gives the discrete expression below
82 (for the horizontal velocity component). The expression for the
83 vertical velocity component is obtained by automatic permutation of
84 the indices but requires a change of sign: this is done through the
85 pseudo-vector `f`. */
86
87 /* trash: u */;
88 coord f = {-1.[0],1.[0]};
89 for (int _i = 0; _i < _N; _i++) /* foreach_face */
90 u.x[] = f.x*(psi[0,1] + psi[-1,1] - psi[0,-1] - psi[-1,-1])/(4.*Delta);
91}
#define u
Definition advection.h:30
int x
Definition common.h:76
void display(const char *commands, bool overwrite=false)
Definition common.h:527
scalar f[]
The primary fields are:
Definition two-phase.h:56
def _i
Definition stencils.h:405
scalar psi[]
Definition stream.h:39
scalar * tracers
Here we set the gradient functions for each tracer (as defined in the user-provided tracers list).
Definition stream.h:41
void event_defaults(void)
Here we set the default boundary conditions for the streamfunction.
Definition stream.h:59
void event_velocity(void)
At every timestep we update the streamfunction field by solving a Poisson equation with the updated ...
Definition stream.h:74
mgstats mgpsi
Definition stream.h:40
scalar omega[]
We allocate the vorticity field , the streamfunction field and a structure to store the statistics o...
Definition stream.h:39
Definition common.h:78
Information about the convergence of the solver is returned in a structure.
Definition poisson.h:112
#define poisson(...)
Finally, we overload the poisson() function called by [two-phase.h]() with our new function.
Definition thermal.h:243
double CFL
Definition utils.h:10