|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
Go to the source code of this file.
Functions | |
| void | event_defaults (void) |
| We change the default gradient function (used for advection) to minmod-limited (rather than the centered default). | |
| void | event_pressure (void) |
At every timestep, but after all the other events for this timestep have been processed (the 'last' keyword), we update the pressure field \(p\) by solving the Poisson equation with variable coefficient \(\beta\). | |
Variables | |
| scalar | p [] |
| We allocate the pressure \(p\) and divergence field \(\zeta\). | |
| scalar | zeta [] |
| vector | beta [] |
| mgstats | mgp |
We change the default gradient function (used for advection) to minmod-limited (rather than the centered default).
Event: defaults (i = 0)
The default density field is set to unity (times the metric).
We reset the multigrid parameters to their default values.
If the viscosity is non-zero, we need to allocate the face-centered viscosity field.
We also initialize the list of tracers to be advected with the VOF function \(f\) (or its complementary function).
We set limiting.
On trees, we ensure that limiting is also applied to prolongation and refinement.
We add the interface and the density to the default display.
We switch to a pure minmod limiter by default for increased robustness.
With the MUSCL scheme we use the CFL depends on the dimension of the problem.
On trees we need to replace the default bilinear refinement/prolongation with linear so that reconstructed values also use slope limiting.
The restriction/refine attributes of the charge density are those of a tracer otherwise the conservation is not guaranteed.
By default the permittivity is unity and other quantities are zero.
Definition at line 48 of file hele-shaw.h.
References gradient, and minmod2().
At every timestep, but after all the other events for this timestep have been processed (the 'last' keyword), we update the pressure field \(p\) by solving the Poisson equation with variable coefficient \(\beta\).
Event: pressure (i++, last)
We then update the velocity field by computing the face pressure gradients.
Definition at line 60 of file hele-shaw.h.
References _i, beta, mgp, p, poisson, u, vector::x, x, and zeta.
| vector beta[] |
Definition at line 40 of file hele-shaw.h.
Referenced by event_acceleration(), event_defaults(), event_pressure(), event_tracer_diffusion(), generic_limiter(), h_relax(), and h_residual().
| mgstats mgp |
Definition at line 41 of file hele-shaw.h.
Referenced by event_pressure().
| scalar p[] |
We allocate the pressure \(p\) and divergence field \(\zeta\).
Flows dominated by friction such as Hele-Shaw flows or flows in porous media (Darcy flows) can be modelled as
\[ \mathbf{u} = \beta\nabla p \]
with \(p\) analogous to a pressure and where \(\beta\) can be a function of space and time. If the fluid is also incompressible, \(p\) needs to verify the Poisson equation
\[ \nabla\cdot(\beta\nabla p) = \zeta \]
\(\beta\) is often a function of the properties of the fluid such as its composition, temperature and/or density etc... which also requires the solution of advection–diffusion equations.
In the following we start from the advection solver and add the definition of the velocity through the Poisson equation for the pressure. This is very similar to what is done for the streamfunction–vorticity Navier–Stokes solver. The \(\beta\) coefficients need to be defined at face locations to match the locations of the face pressure gradients (and the face velocity components). These two sets of coefficients are stored in a vector field. We also allocate space to store the statistics of the Poisson solver.
Definition at line 39 of file hele-shaw.h.
Referenced by event_pressure().
| scalar zeta[] |
Definition at line 39 of file hele-shaw.h.
Referenced by event_pressure(), gotm_step(), and meanflow_updategrid().