|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
Go to the source code of this file.
Functions | |
| static void | momentum_refine (Point point, scalar u) |
| static void | momentum_restriction (Point point, scalar u) |
| void | event_defaults (void) |
| We switch-off the default advection scheme of the centered solver. | |
| void | event_stability (void) dtmax |
| We need to overload the stability event so that the CFL is taken into account (because we set stokes to true). | |
| for (int _d=0;_d< 2;_d++) static double boundary_q1_x(Point neighbor | |
| We will transport the two components of the momentum, \(q_1=f \rho_1
\mathbf{u}\) and \(q_2=(1 - f) \rho_2 \mathbf{u}\). | |
| void | event_vof (void) |
| Event: vof (i++) | |
| void | event_tracer_advection (void) |
| We set the list of interfaces back to its default value. | |
Variables | |
| Point | point |
| Point scalar | q1 |
| Point scalar bool * | data |
| Point scalar | q2 |
| static scalar * | interfaces1 = NULL |
| We overload the vof() event to transport consistently the volume fraction and the momentum of each phase. | |
We switch-off the default advection scheme of the centered solver.
Event: defaults (i = 0)
On trees, the refinement and restriction functions above rely on the volume fraction field f being refined/restricted before the components of velocity. To ensure this, we move f to the front of the field list (all).
We then set the refinement and restriction functions for the components of the velocity field. The boundary conditions on \(\mathbf{u}\) now depend on those on \(f\).
Definition at line 40 of file conserving.h.
References all, dimension, f, scalar::i, i, list_add(), momentum_refine(), momentum_restriction(), stokes, u, and x.
We need to overload the stability event so that the CFL is taken into account (because we set stokes to true).
Event: stability (i++)
We need to overload the stability event so that the CFL is taken into account (because we set stokes to true).
We need to overload the stability event so that the CFL is taken into account (because we set stokes to true).
We need to overload the stability event so that the CFL is taken into account (because we set stokes to true).
The surface tension scheme is time-explicit so the maximum timestep is the oscillation period of the smallest capillary wave.
\[ T = \sqrt{\frac{\rho_{m}\Delta_{min}^3}{\pi\sigma}} \]
with \(\rho_m=(\rho_1+\rho_2)/2.\) and \(\rho_1\), \(\rho_2\) the densities on either side of the interface.
We need to overload the stability event so that the CFL is taken into account (because we set stokes to true).
The surface tension scheme is time-explicit so the maximum timestep is the oscillation period of the smallest capillary wave.
\[ T = \sqrt{\frac{\rho_{m}\Delta_{min}^3}{\pi\sigma}} \]
with \(\rho_m=(\rho_1+\rho_2)/2.\) and \(\rho_1\), \(\rho_2\) the densities on either side of the interface.
We need to overload the stability event so that the CFL is taken into account (because we set stokes to true).
Event: stability (i++)
We first compute the minimum and maximum values of \(\alpha/f_m = 1/\rho\), as well as \(\Delta_{min}\).
The maximum timestep is set using the sum of surface tension coefficients.
We first compute the minimum and maximum values of \(\alpha/f_m = 1/\rho\), as well as \(\Delta_{min}\).
The maximum timestep is set using the sum of surface tension coefficients.
Definition at line 110 of file all-mach.h.
References _i, alpha, c, CFL, CFLa, CFLac, cube(), d, dmin, dry, dt, dtmax, dtnext(), fm, G, h, HUGE, max, min, pi, point, sigma, sound_speed(), stokes, timestep(), u, uf, vector::x, and x.
We set the list of interfaces back to its default value.
Event: tracer_advection (i++)
Definition at line 198 of file conserving.h.
References interfaces, and interfaces1.
Event: vof (i++)
We allocate two temporary vector fields to store the two components of the momentum and set the boundary conditions and prolongation functions. The boundary conditions on \(q_1\) and \(q_2\) depend on the boundary conditions on \(f\).
We split the total momentum \(q\) into its two components \(q1\) and \(q2\) associated with \(f\) and \(1 - f\) respectively.
Momentum \(q2\) is associated with \(1 - f\), so we set the inverse attribute to true. We use the same slope-limiting as for the velocity field.
We associate the transport of \(q1\) and \(q2\) with \(f\) and transport all fields consistently using the VOF scheme.
We recover the advected velocity field using the total momentum and the density
We set the list of interfaces to NULL so that the default vof() event does nothing (otherwise we would transport \(f\) twice).
Definition at line 122 of file conserving.h.
References _i, clamp(), dimension, f, free(), scalar::i, i, interfaces, interfaces1, list_add(), list_concat(), nboundary, q1, q2, rho, rho1, rho2, s, tracers, u, vof_advection(), and x.
| for | ( | ) |
We will transport the two components of the momentum, \(q_1=f \rho_1 \mathbf{u}\) and \(q_2=(1 - f) \rho_2 \mathbf{u}\).
Similarly, on trees we need prolongation functions which also follow this definition.
We will need to impose boundary conditions which match this definition. This is done using the functions below.
This file implements momentum-conserving VOF advection of the velocity components for the two-phase Navier–Stokes solver.
On trees, we first define refinement and restriction functions which guarantee conservation of each component of the total momentum. Note that these functions do not guarantee conservation of momentum for each phase.
Definition at line 16 of file conserving.h.
References cm, dimension, f, point, refine_bilinear(), rho, SEPS, u, and x.
Referenced by event_defaults().
Definition at line 26 of file conserving.h.
References cm, dimension, f, rho, SEPS, u, and x.
Referenced by event_defaults().
Definition at line 86 of file conserving.h.
We overload the vof() event to transport consistently the volume fraction and the momentum of each phase.
Definition at line 119 of file conserving.h.
Referenced by event_tracer_advection(), and event_vof().
| Point point |
Definition at line 86 of file conserving.h.
Referenced by adapt_wavelet(), advance_saint_venant(), advect(), alloc_children(), allocated(), balance(), box_boundary_level(), box_boundary_level(), box_boundary_level(), box_boundary_level_normal(), box_boundary_level_normal(), box_boundary_level_tangent(), centroids_curvature_fit(), check_two_one(), coarsen_cell(), coarsen_cell_recursive(), compile_expression(), curvature(), dirichlet(), dirichlet_gradient(), distance(), embed_area_center(), embed_flux(), embed_force(), embed_fraction_refine(), embed_geometry(), embed_gradient(), embed_vorticity(), evaluate_expression(), event_acceleration(), event_pressure(), event_properties(), event_stability(), event_viscous_term(), face_fraction(), facet_normal(), fE_refine(), for(), foreach(), foreach_boundary_dir(), foreach_cell_post_root(), foreach_child(), foreach_face_generic(), foreach_level_stencil(), foreach_mem(), foreach_neighbor(), foreach_slice_x(), foreach_slice_y(), foreach_stencil(), foreach_stencil_generic(), foreach_tree(), foreach_visible(), fraction_refine(), glvertex_normal3d(), height_curvature(), height_curvature_fit(), height_myc_normal(), height_normal(), height_position(), heights(), if(), increment_neighbors(), init_grid(), input_gfs(), interface_area(), interface_normal(), interp(), interpolate(), interpolate_array(), inverse_wavelet(), is_local_prolongation(), is_refined_check(), locals_pids(), locate(), locate(), matrix(), mg_cycle(), momentum_refine(), mpi_boundary_coarsen(), mpi_boundary_refine(), mpi_boundary_update_buffers(), msolve(), multigrid_debug(), multigrid_restriction(), neighbor(), neumann(), no_coalescence(), normal_neighbor(), OMP_PARALLEL(), output_tree(), periodic_boundary_level_x(), POINT_VARIABLES(), prolongation_elevation(), rcv_append(), rcv_pid_append(), reconstruction(), recursive(), refine(), refine_bilinear(), refine_biquadratic(), refine_cell(), refine_cm_axi(), refine_distance(), refine_face(), refine_face_solenoidal(), refine_face_x_axi(), refine_level(), refine_linear(), relax(), relax_diffusion(), relax_GN(), relax_hydro(), relax_nh(), relax_nh(), relax_viscosity(), residual(), residual_diffusion(), residual_nh(), residual_nh(), residual_nh3(), restore_mpi(), restriction_embed_linear(), restriction_face(), segment_flux(), solve(), terrain(), tree_boundary_level(), unrefine(), update_cache_f(), update_distance(), val(), VARIABLES(), vertical_diffusion(), vertical_remapping(), vertical_velocity(), vof_concentration_refine(), and wavelet().
| scalar q1 |
| scalar q2 |