|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
Go to the source code of this file.
Macros | |
| #define | wmin(h, hmin) (h >= hmin ? 0. : sq((hmin - h)/hmin)) |
| #define | wmax(h, hmax) (h <= hmax ? 0. : sq((h - hmax)/hmax)) |
Functions | |
| void | event_viscous_term (void) |
| Event: viscous_term (i++) | |
Variables | |
| double | omr |
| double | Cm |
| double * | hmin |
| The maximum and minimum layer thicknesses ( \(h_k^+\) and \(h_k^-\) respectively in H&H, 2000) are provided by the user and used to define the weighing factors appearing in the definition of the "diapycnal
mixing velocities" \(\omega_k^+\) and \(\omega_k^-\) in H&H, 2000. | |
| double * | hmax |
| double | dt |
| These come from the multilayer solver. | |
| double | dry |
Definition at line 25 of file entrainment.h.
Definition at line 24 of file entrainment.h.
Event: viscous_term (i++)
If the average entrainment velocity is negative or null, we do nothing.
We compute the volumed-averaged entrainment/detrainment for each layer (bottom layer excepted).
The local entrainment/detrainment is then applied as the sum of a global contribution and a local contribution which depends on the local thickness.
These are the entrainment and detrainment terms for the velocity components in eqs. 1) and 2) of H&H, 2000.
We then update h and u.
Definition at line 33 of file entrainment.h.
References _i, _layer, Cm, dimension, dry, dt, dv, flux, foreach_layer, h, hmax, hmin, hn, l, nl, o, omr, point, u, un, w, wmax, wmin, vector::x, and x.
| double Cm |
Definition at line 15 of file entrainment.h.
Referenced by event_viscous_term().
| double dry |
Definition at line 30 of file entrainment.h.
Referenced by area_z2(), bflux(), box_matrix(), event_acceleration(), event_face_fields(), event_init(), event_pressure(), event_viscous_term(), fault(), geostrophic_velocity(), horizontal_diffusion(), prolongation_elevation(), refine_elevation(), relax_GN(), residual_nh(), restriction_elevation(), update_green_naghdi(), and vertical_remapping().
|
extern |
These come from the multilayer solver.
This function approximates implicitly the EHD equation set given by the electric Poisson equation and the ohmic conduction term of the charge density conservation (the advection term is computed elsewhere using a tracer advection scheme),
\[ \nabla \cdot( \epsilon \nabla \phi^{n+1}) = -\rho_e^{n+1} \quad \mathrm{and} \quad (\rho_e^{n+1}-\rho_e^n) = \Delta t\nabla \cdot (K\nabla \phi^{n+1}) \]
where \(\rho_e\) is the charge density, and \(\phi\) is the electric potential. \(K\) and \(\epsilon\) are the conductivity and permittivity respectively.
Substituting the Poisson equation into the conservation equation gives the following time-implicit scheme,
\[ \nabla \cdot [(K \, \Delta t + \epsilon) \nabla \phi^{n+1}] = -\rho_e^{n} \]
We need the Poisson solver and the timestep dt.
The run() function below implements a generic time loop which executes events until termination.
The timestep dt can be accessed as a global variable.
Definition at line 18 of file predictor-corrector.h.
Referenced by event_viscous_term().
| double * hmax |
Definition at line 23 of file entrainment.h.
Referenced by event_viscous_term().
|
extern |
The maximum and minimum layer thicknesses ( \(h_k^+\) and \(h_k^-\) respectively in H&H, 2000) are provided by the user and used to define the weighing factors appearing in the definition of the "diapycnal mixing velocities" \(\omega_k^+\) and \(\omega_k^-\) in H&H, 2000.
Referenced by event_viscous_term().
|
extern |
This implements the entrainment/detrainment between isopycnal layers briefly described in Hurlburt & Hogan, 2000, section 2, equations 1), 2) and 3).
omr is the average entrainment velocity between layers ( \(\tilde{\omega}_k\) in H&H, 2000) and Cm is the coefficient of additional interfacial friction associated with entrainment ( \(C_M\) in H&H, 2000). They are parameters provided by the user.
Referenced by event_viscous_term().