MFC: Simulation
High-fidelity multiphase flow simulation
|
Functions/Subroutines | |
program | __m_bubbles_el_fpp_f90__ |
subroutine | s_initialize_bubbles_el_module (q_cons_vf) |
Initializes the lagrangian subgrid bubble solver. | |
subroutine | s_start_lagrange_inputs () |
The purpose of this procedure is to start lagrange bubble parameters applying nondimensionalization if needed. | |
subroutine | s_read_input_bubbles (q_cons_vf) |
The purpose of this procedure is to obtain the initial bubbles' information. | |
subroutine | s_add_bubbles (inputbubble, q_cons_vf, bub_id) |
The purpose of this procedure is to obtain the information of the bubbles when starting fresh. | |
subroutine | s_restart_bubbles (bub_id, save_count) |
The purpose of this procedure is to obtain the information of the bubbles from a restart point. | |
subroutine | s_compute_el_coupled_solver (q_cons_vf, q_prim_vf, rhs_vf, stage) |
Contains the two-way and one-way Euler-Lagrange coupled algorithm, including the bubble dynamics subroutines. | |
subroutine | s_compute_cson_from_pinf (bub_id, q_prim_vf, pinf, cell, rhol, gamma, pi_inf, cson) |
This procedure computes the speed of sound from a given driving pressure. | |
subroutine | s_smear_voidfraction () |
The purpose of this subroutine is to smear the effect of the bubbles in the Eulerian framework. | |
subroutine | s_add_rhs_sources (q_cons_vf, q_prim_vf, rhs_vf) |
The purpose of this subroutine is to obtain the bubble source terms based on Maeda and Colonius (2018) and add them to the RHS scalar field. | |
subroutine | s_get_pinf (bub_id, q_prim_vf, ptype, f_pinfl, cell, preterm1, term2, romega) |
The purpose of this procedure is obtain the bubble driving pressure p_inf. | |
subroutine | s_update_lagrange_tdv_rk (stage) |
This subroutine updates the Lagrange variables using the tvd RK time steppers. The time derivative of the bubble variables must be stored at every stage to avoid precision errors. | |
subroutine | s_update_tmp_rkck (rkstep, q_cons_ts, rhs_ts, lag_largestep) |
This subroutine updates the Euler-Lagrange temporal variables before entering to the next time-stage in the RKCK stepper. | |
subroutine | s_calculate_rkck_truncation_error (rkck_errmax) |
This subroutine calculates the maximum error between the 4th and 5th order Runge-Kutta-Cash-Karp solutions for the same time step size. If the errors are smaller than a tolerance, then the algorithm employs the 5th order solution, while if not, both eulerian/lagrangian variables are re-calculated with a smaller time step size. | |
subroutine | s_update_rkck (q_cons_ts) |
This subroutine updates the conservative fields and the lagrangian variables after accepting the performed time step. | |
subroutine | s_compute_rkck_dt (lag_largestep, restart_rkck_step, rkck_errmax) |
This subroutine computes the next time step in the adaptive RKCK stepper in the CPU. | |
subroutine | s_locate_cell (pos, cell, scoord) |
This subroutine returns the computational coordinate of the cell for the given position. | |
subroutine | s_transfer_data_to_tmp () |
This subroutine transfer data into the temporal variables. | |
logical function | particle_in_domain (pos_part) |
The purpose of this procedure is to determine if the global coordinates of the bubbles are present in the current MPI processor (including ghost cells). | |
logical function | particle_in_domain_physical (pos_part) |
The purpose of this procedure is to determine if the lagrangian bubble is located in the physical domain. The ghost cells are not part of the physical domain. | |
subroutine | s_gradient_dir (q, dq, dir) |
The purpose of this procedure is to calculate the gradient of a scalar field along the x, y and z directions following a second-order central difference considering uneven widths. | |
subroutine | s_write_lag_particles (qtime) |
Subroutine that writes on each time step the changes of the lagrangian bubbles. | |
subroutine | s_write_void_evol (qtime) |
Subroutine that writes some useful statistics related to the volume fraction of the particles (void fraction) in the computatioational domain on each time step. | |
subroutine | s_write_restart_lag_bubbles (t_step) |
Subroutine that writes the restarting files for the particles in the lagrangian solver. | |
subroutine | s_calculate_lag_bubble_stats () |
This procedure calculates the maximum and minimum radius of each bubble. | |
subroutine | s_write_lag_bubble_stats () |
Subroutine that writes the maximum and minimum radius of each bubble. | |
subroutine | s_remove_lag_bubble (bub_id) |
The purpose of this subroutine is to remove one specific particle if dt is too small. | |
subroutine | s_finalize_lagrangian_solver () |
The purpose of this subroutine is to deallocate variables. | |
program __m_bubbles_el_fpp_f90__ |
|
private |
The purpose of this procedure is to determine if the global coordinates of the bubbles are present in the current MPI processor (including ghost cells).
pos_part | Spatial coordinates of the bubble |
|
private |
The purpose of this procedure is to determine if the lagrangian bubble is located in the physical domain. The ghost cells are not part of the physical domain.
pos_part | Spatial coordinates of the bubble |
|
private |
The purpose of this procedure is to obtain the information of the bubbles when starting fresh.
inputBubble | Bubble information |
q_cons_vf | Conservative variables |
bub_id | Local id of the bubble |
|
private |
The purpose of this subroutine is to obtain the bubble source terms based on Maeda and Colonius (2018) and add them to the RHS scalar field.
q_cons_vf | Conservative variables |
q_prim_vf | Conservative variables |
rhs_vf | Time derivative of the conservative variables |
|
private |
This procedure calculates the maximum and minimum radius of each bubble.
|
private |
This subroutine calculates the maximum error between the 4th and 5th order Runge-Kutta-Cash-Karp solutions for the same time step size. If the errors are smaller than a tolerance, then the algorithm employs the 5th order solution, while if not, both eulerian/lagrangian variables are re-calculated with a smaller time step size.
rkck_errmax | Truncation error |
|
private |
This procedure computes the speed of sound from a given driving pressure.
bub_id | Bubble id |
q_prim_vf | Primitive variables |
pinf | Driving pressure |
cell | Bubble cell |
rhol | Liquid density |
gamma | Liquid specific heat ratio |
pi_inf | Liquid stiffness |
cson | Calculated speed of sound |
|
private |
Contains the two-way and one-way Euler-Lagrange coupled algorithm, including the bubble dynamics subroutines.
q_cons_vf | Conservative variables |
q_prim_vf | Primitive variables |
rhs_vf | Calculated change of conservative variables |
stage | Current stage in the time-stepper algorithm |
|
private |
This subroutine computes the next time step in the adaptive RKCK stepper in the CPU.
lag_largestep | Negative radius flag |
restart_rkck_step | Restart the current time step |
rkck_errmax | Truncation error |
|
private |
The purpose of this subroutine is to deallocate variables.
|
private |
The purpose of this procedure is obtain the bubble driving pressure p_inf.
bub_id | Particle identifier |
q_prim_vf | Primitive variables |
ptype | 1: p at infinity, 2: averaged P at the bubble location |
f_pinfl | Driving pressure |
cell | Bubble cell |
Romega | Control volume radius |
|
private |
The purpose of this procedure is to calculate the gradient of a scalar field along the x, y and z directions following a second-order central difference considering uneven widths.
q | Input scalar field |
dq | Output gradient of q |
dir | Gradient spatial direction |
subroutine __m_bubbles_el_fpp_f90__::s_initialize_bubbles_el_module | ( | type(scalar_field), dimension(sys_size), intent(inout) | q_cons_vf | ) |
Initializes the lagrangian subgrid bubble solver.
q_cons_vf | Initial conservative variables |
|
private |
This subroutine returns the computational coordinate of the cell for the given position.
pos | Input coordinates |
cell | Computational coordinate of the cell |
scoord | Calculated particle coordinates |
|
private |
The purpose of this procedure is to obtain the initial bubbles' information.
q_cons_vf | Conservative variables |
|
private |
The purpose of this subroutine is to remove one specific particle if dt is too small.
bub_id | Particle id |
|
private |
The purpose of this procedure is to obtain the information of the bubbles from a restart point.
bub_id | Local ID of the particle |
save_count | File identifier |
|
private |
The purpose of this subroutine is to smear the effect of the bubbles in the Eulerian framework.
|
private |
The purpose of this procedure is to start lagrange bubble parameters applying nondimensionalization if needed.
|
private |
This subroutine transfer data into the temporal variables.
|
private |
This subroutine updates the Lagrange variables using the tvd RK time steppers. The time derivative of the bubble variables must be stored at every stage to avoid precision errors.
stage | Current tvd RK stage |
|
private |
This subroutine updates the conservative fields and the lagrangian variables after accepting the performed time step.
q_cons_ts | Conservative variables |
|
private |
This subroutine updates the Euler-Lagrange temporal variables before entering to the next time-stage in the RKCK stepper.
RKstep | Current time step in the RKCK adaptive stepper |
q_cons_ts | Conservative variables |
rhs_ts | Time derivatives of the conservative variables |
lag_largestep | Negative radius flag |
|
private |
Subroutine that writes the maximum and minimum radius of each bubble.
|
private |
Subroutine that writes on each time step the changes of the lagrangian bubbles.
q_time | Current time |
|
private |
Subroutine that writes the restarting files for the particles in the lagrangian solver.
t_step | Current time step |
|
private |
Subroutine that writes some useful statistics related to the volume fraction of the particles (void fraction) in the computatioational domain on each time step.
q_time | Current time |