Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
hele-shaw.h
Go to the documentation of this file.
1/** @file hele-shaw.h
2 */
3/**
4# Hele-Shaw flow solver
5
6Flows dominated by friction such as *Hele-Shaw flows* or flows in
7porous media (*Darcy flows*) can be modelled as
8\f[
9\mathbf{u} = \beta\nabla p
10\f]
11with \f$p\f$ analogous to a pressure and where \f$\beta\f$ can be a function
12of space and time. If the fluid is also incompressible, \f$p\f$ needs to
13verify the Poisson equation
14\f[
15\nabla\cdot(\beta\nabla p) = \zeta
16\f]
17
18\f$\beta\f$ is often a function of the properties of the fluid such as its
19composition, temperature and/or density etc... which also requires the
20solution of advection--diffusion equations.
21
22In the following we start from the advection solver and add the
23definition of the velocity through the Poisson equation for the
24pressure. This is very similar to what is done for the
25[streamfunction--vorticity](navier-stokes/stream.h) Navier--Stokes
26solver. */
27
28#include "advection.h"
29#include "poisson.h"
30
31/**
32We allocate the pressure \f$p\f$ and divergence field \f$\zeta\f$. The
33\f$\beta\f$ coefficients need to be defined at face locations to
34match the locations of the face pressure gradients (and the
35face velocity components). These two sets of coefficients are
36stored in a vector field. We also allocate space to store the
37statistics of the Poisson solver. */
38
42
43/**
44We change the default gradient function (used for advection) to
45minmod-limited (rather than the centered default). */
46
47/** @brief Event: defaults (i = 0) */
48void event_defaults (void)
49{
51}
52
53/**
54At every timestep, but after all the other events for this timestep
55have been processed (the '`last`' keyword), we update the pressure
56field \f$p\f$ by solving the Poisson equation with variable coefficient
57\f$\beta\f$. */
58
59/** @brief Event: pressure (i++, last) */
60void event_pressure (void)
61{
62 mgp = poisson (p, zeta, beta);
63
64 /**
65 We then update the velocity field by computing the face pressure
66 gradients. */
67
68 /* trash: u */;
69 for (int _i = 0; _i < _N; _i++) /* foreach_face */
70 u.x[] = beta.x[]*(p[] - p[-1])/Delta;
71}
#define u
Definition advection.h:30
double(* gradient)(double, double, double)
Definition advection.h:31
int x
Definition common.h:76
void event_pressure(void)
At every timestep, but after all the other events for this timestep have been processed (the 'last' k...
Definition hele-shaw.h:60
scalar zeta[]
Definition hele-shaw.h:39
void event_defaults(void)
We change the default gradient function (used for advection) to minmod-limited (rather than the cente...
Definition hele-shaw.h:48
vector beta[]
Definition hele-shaw.h:40
scalar p[]
We allocate the pressure and divergence field .
Definition hele-shaw.h:39
mgstats mgp
Definition hele-shaw.h:41
def _i
Definition stencils.h:405
Information about the convergence of the solver is returned in a structure.
Definition poisson.h:112
scalar x
Definition common.h:47
#define poisson(...)
Finally, we overload the poisson() function called by [two-phase.h]() with our new function.
Definition thermal.h:243
double minmod2(double s0, double s1, double s2)
Definition utils.h:225