|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
#include "multilayer.h"#include "predictor-corrector.h"#include "riemann.h"#include "elevation.h"#include "gauges.h"Go to the source code of this file.
Macros | |
| #define | radiation(ref) (sqrt (G*max(h[],0.)) - sqrt(G*max((ref) - zb[], 0.))) |
Functions | |
| static trace void | advance_saint_venant (scalar *soutput, scalar *sinput, scalar *updates, double dt) |
| We need to overload the default advance function of the predictor-corrector scheme, because the evolving variables ( \(h\) and \(\mathbf{u}\)) are not the conserved variables \(h\) and \(h\mathbf{u}\). | |
| static void | refine_eta (Point point, scalar eta) |
| When using an adaptive discretisation (i.e. | |
| static void | restriction_eta (Point point, scalar eta) |
| trace double | update_saint_venant (scalar *evolving, scalar *updates, double dtmax) |
| void | event_defaults (void) |
| Event: defaults (i = 0) | |
| void | event_init (void) |
| The event below will happen after all the other initial events to take into account user-defined field initialisations. | |
| void | event_cleanup (void) |
| At the end of the simulation, we free the memory allocated in the defaults* event. | |
Variables | |
| scalar | zb [] |
| scalar | h [] |
| scalar | eta [] |
| vector | u [] |
| double | G = 1. |
| The only physical parameter is the acceleration of gravity G. | |
| double | dry = 1e-10 |
| int | nl = 1 |
| By default there is only a single layer i.e. | |
| scalar * | evolving = NULL |
| The generic time-integration scheme in predictor-corrector.h needs to know which fields are updated. | |
This can be used to implement open boundary conditions at low Froude numbers. The idea is to set the velocity normal to the boundary so that the water level relaxes towards its desired value (ref).
Definition at line 435 of file saint-venant.h.
|
static |
We need to overload the default advance function of the predictor-corrector scheme, because the evolving variables ( \(h\) and \(\mathbf{u}\)) are not the conserved variables \(h\) and \(h\mathbf{u}\).
In the case of multiple layers we add the viscous friction between layers.
Definition at line 97 of file saint-venant.h.
References _i, avector, dh, dimension, dry, dt, eta, l, nl, point, vertical_viscosity(), x, and zb.
Referenced by event_defaults().
Event: defaults (i = 0)
We use the main time loop (in the predictor-corrector scheme) to setup the initial defaults.
By default, all the layers have the same relative thickness.
We overload the default 'advance' and 'update' functions of the predictor-corrector scheme and setup the prolongation and restriction methods on trees.
On trees we make sure that slope-limiting is also used for refinement and prolongation. The prolongation/restriction functions for \(\eta\) are set and they depend on boundary conditions on \(z_b\) and \(h\).
We setup the default display.
Definition at line 343 of file saint-venant.h.
References _i, advance, advance_saint_venant(), assert, dimension, display(), eta, evolving, h, l, layer, list_append(), list_concat(), list_copy(), nl, qmalloc, refine_eta(), refine_linear(), restriction_eta(), restriction_volume_average(), s, set_prolongation(), set_restriction(), u, ul, update, update_saint_venant(), vectors_append(), w, wl, vector::x, x, and zb.
When using an adaptive discretisation (i.e.
a tree)., we need to make sure that \(\eta\) is maintained as \(z_b + h\) whenever cells are refined or restricted.
Definition at line 140 of file saint-venant.h.
Referenced by event_defaults().
Definition at line 146 of file saint-venant.h.
Referenced by event_defaults().
Various approximate Riemann solvers are defined in [riemann.h]().
We first recover the currently evolving height and velocity (as set by the predictor-corrector scheme).
Fh* and Fq will contain the fluxes for \(h\) and \(h\mathbf{u}\) respectively and S is necessary to store the asymmetric topographic source term.
The gradients are stored in locally-allocated fields. First-order reconstruction is used for the gradient fields.
We go through each layer.
We recover the velocity field for the current layer and compute its gradient (for the first layer the gradient has already been computed above).
The faces which are "wet" on at least one side are traversed.
The gradients computed above are used to reconstruct the left and right states of the primary fields \(h\), \(\mathbf{u}\), \(z_b\). The "interface" topography \(z_{lr}\) is reconstructed using the hydrostatic reconstruction of Audusse et al, 2004
We can now call one of the approximate Riemann solvers to get the fluxes.
In the case of adaptive refinement, care must be taken to ensure well-balancing at coarse/fine faces (see [notes/balanced.tm]()).
We store the divergence of the fluxes in the update fields. Note that these are updates for \(h\) and \(h\mathbf{u}\) (not \(\mathbf{u}\)).
For multiple layers we need to store the divergence in each layer.
We also need to add the metric terms. They can be written (see eq. (8) of Popinet, 2011)
\[ S_g = h \left(\begin{array}{c} 0\\ \frac{g}{2} h \partial_{\lambda} m_{\theta} + f_G u_y\\ \frac{g}{2} h \partial_{\theta} m_{\lambda} - f_G u_x \end{array}\right) \]
with
\[ f_G = u_y \partial_{\lambda} m_{\theta} - u_x \partial_{\theta} m_{\lambda} \]
For multiple layers we need to add fluxes between layers.
Definition at line 160 of file saint-venant.h.
References _i, avector, cell(), cm, coarse(), dh, dimension, dry, dtmax, dx, eta, evolving, fm, G, gradients(), h, hn, is_prolongation, kurganov(), l, layer, max, neighbor(), nl, refine_linear(), s, S, sq(), u, vertical_fluxes(), wl, vector::x, x, vector::y, zb, zero(), and zn.
Referenced by event_defaults(), and update_green_naghdi().
Definition at line 59 of file saint-venant.h.
Referenced by advance_saint_venant(), and update_saint_venant().
| scalar eta[] |
Definition at line 50 of file saint-venant.h.
Referenced by advance_saint_venant(), event_defaults(), event_init(), refine_eta(), restriction_eta(), and update_saint_venant().
The generic time-integration scheme in predictor-corrector.h needs to know which fields are updated.
The generic time-integration scheme in predictor-corrector.h needs to know which fields are updated i.e.
Time integration will be done with a generic predictor-corrector scheme. The list will be constructed in the defaults* event below.
Definition at line 88 of file saint-venant.h.
Referenced by event_cleanup(), event_defaults(), run(), and update_saint_venant().
| double G = 1. |
The only physical parameter is the acceleration of gravity G.
Cells are considered "dry" when the water depth is less than the dry parameter (this should not require tweaking).
Definition at line 58 of file saint-venant.h.
Referenced by update_saint_venant().
| scalar h[] |
Definition at line 50 of file saint-venant.h.
Referenced by event_defaults(), event_init(), refine_eta(), restriction_eta(), and update_saint_venant().
| int nl = 1 |
By default there is only a single layer i.e.
this is the classical Saint-Venant system above. This can be changed by setting nl to a different value. The extension of the Saint-Venant system to multiple layers is implemented in [multilayer.h]().
Definition at line 68 of file saint-venant.h.
Referenced by advance_saint_venant(), event_defaults(), and update_saint_venant().
| vector u[] |
Definition at line 51 of file saint-venant.h.
Referenced by event_defaults(), and update_saint_venant().
| scalar zb[] |
The Saint-Venant equations can be written in integral form as the hyperbolic system of conservation laws
\[ \partial_t \int_{\Omega} \mathbf{q} d \Omega = \int_{\partial \Omega} \mathbf{f} ( \mathbf{q}) \cdot \mathbf{n}d \partial \Omega - \int_{\Omega} hg \nabla z_b \]
where \(\Omega\) is a given subset of space, \(\partial \Omega\) its boundary and \(\mathbf{n}\) the unit normal vector on this boundary. For conservation of mass and momentum in the shallow-water context, \(\Omega\) is a subset of bidimensional space and \(\mathbf{q}\) and \(\mathbf{f}\) are written
\[ \mathbf{q} = \left(\begin{array}{c} h\\ hu_x\\ hu_y \end{array}\right), \;\;\;\;\;\; \mathbf{f} (\mathbf{q}) = \left(\begin{array}{cc} hu_x & hu_y\\ hu_x^2 + \frac{1}{2} gh^2 & hu_xu_y\\ hu_xu_y & hu_y^2 + \frac{1}{2} gh^2 \end{array}\right) \]
where \(\mathbf{u}\) is the velocity vector, \(h\) the water depth and \(z_b\) the height of the topography. See also Popinet, 2011 for a more detailed introduction.
The primary fields are the water depth \(h\), the bathymetry \(z_b\) and the flow speed \(\mathbf{u}\). \(\eta\) is the water level i.e. \(z_b + h\). Note that the order of the declarations is important as \(z_b\) needs to be refined before \(h\) and \(h\) before \(\eta\).
Definition at line 50 of file saint-venant.h.
Referenced by advance_saint_venant(), event_defaults(), event_init(), refine_eta(), restriction_eta(), and update_saint_venant().