|
MFC: Simulation
High-fidelity multiphase flow simulation
|
This module is used to to compute the volume-averaged bubble model. More...
Functions/Subroutines | |
| impure subroutine | s_initialize_bubbles_el_module (q_cons_vf) |
| Initializes the lagrangian subgrid bubble solver. | |
| impure subroutine | s_start_lagrange_inputs () |
| The purpose of this procedure is to start lagrange bubble parameters applying nondimensionalization if needed. | |
| impure subroutine | s_read_input_bubbles (q_cons_vf) |
| The purpose of this procedure is to obtain the initial bubbles' information. | |
| impure 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. | |
| impure 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_bubble_el_dynamics (q_prim_vf, stage) |
| Contains the bubble dynamics subroutines. | |
| subroutine | s_compute_bubbles_el_source (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_compute_cson_from_pinf (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_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. | |
| impure 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_locate_cell (pos, cell, scoord) |
| This subroutine returns the computational coordinate of the cell for the given position. | |
| impure 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. | |
| impure subroutine | s_write_lag_particles (qtime) |
| Subroutine that writes on each time step the changes of the lagrangian bubbles. | |
| impure 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. | |
| impure 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. | |
| impure subroutine | s_write_lag_bubble_stats () |
| Subroutine that writes the maximum and minimum radius of each bubble. | |
| impure subroutine | s_remove_lag_bubble (bub_id) |
| The purpose of this subroutine is to remove one specific particle if dt is too small. | |
| impure subroutine | s_finalize_lagrangian_solver () |
| The purpose of this subroutine is to deallocate variables. | |
Variables | |
| integer, dimension(:, :), allocatable | lag_id |
| Global and local IDs. | |
| real(wp), dimension(:), allocatable | bub_r0 |
| Initial bubble radius. | |
| real(wp), dimension(:), allocatable | rmax_stats |
| Maximum radius. | |
| real(wp), dimension(:), allocatable | rmin_stats |
| Minimum radius. | |
| real(wp), dimension(:), allocatable | gas_mg |
| Bubble's gas mass. | |
| real(wp), dimension(:), allocatable | gas_betat |
| heatflux model (Preston et al., 2007) | |
| real(wp), dimension(:), allocatable | gas_betac |
| massflux model (Preston et al., 2007) | |
| real(wp), dimension(:), allocatable | bub_dphidt |
| subgrid velocity potential (Maeda & Colonius, 2018) | |
| real(wp), dimension(:, :), allocatable | gas_p |
| Pressure in the bubble. | |
| real(wp), dimension(:, :), allocatable | gas_mv |
| Vapor mass in the bubble. | |
| real(wp), dimension(:, :), allocatable | intfc_rad |
| Bubble radius. | |
| real(wp), dimension(:, :), allocatable | intfc_vel |
| Velocity of the bubble interface. | |
| real(wp), dimension(:, :, :), allocatable | mtn_pos |
| Bubble's position. | |
| real(wp), dimension(:, :, :), allocatable | mtn_posprev |
| Bubble's previous position. | |
| real(wp), dimension(:, :, :), allocatable | mtn_vel |
| Bubble's velocity. | |
| real(wp), dimension(:, :, :), allocatable | mtn_s |
| Bubble's computational cell position in real format. | |
| real(wp), dimension(:, :), allocatable | intfc_draddt |
| Time derivative of bubble's radius. | |
| real(wp), dimension(:, :), allocatable | intfc_dveldt |
| Time derivative of bubble's interface velocity. | |
| real(wp), dimension(:, :), allocatable | gas_dpdt |
| Time derivative of gas pressure. | |
| real(wp), dimension(:, :), allocatable | gas_dmvdt |
| Time derivative of the vapor mass in the bubble. | |
| real(wp), dimension(:, :, :), allocatable | mtn_dposdt |
| Time derivative of the bubble's position. | |
| real(wp), dimension(:, :, :), allocatable | mtn_dveldt |
| Time derivative of the bubble's velocity. | |
| integer, private | lag_num_ts |
| Number of time stages in the time-stepping scheme. | |
| integer | nbubs |
| Number of bubbles in the local domain. | |
| real(wp) | rmax_glb |
| real(wp) | rmin_glb |
| Maximum and minimum bubbe size in the local domain Projection of the lagrangian particles in the Eulerian framework. | |
| type(scalar_field), dimension(:), allocatable | q_beta |
| integer | q_beta_idx |
| Size of the q_beta vector field. | |
This module is used to to compute the volume-averaged bubble model.
| logical function m_bubbles_el::particle_in_domain | ( | real(wp), dimension(3), intent(in) | 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).
| pos_part | Spatial coordinates of the bubble |
| logical function m_bubbles_el::particle_in_domain_physical | ( | real(wp), dimension(3), intent(in) | 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.
| pos_part | Spatial coordinates of the bubble |
| impure subroutine m_bubbles_el::s_add_bubbles | ( | real(wp), dimension(8), intent(in) | inputbubble, |
| type(scalar_field), dimension(sys_size), intent(in) | q_cons_vf, | ||
| integer, intent(in) | bub_id ) |
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 |
| subroutine m_bubbles_el::s_calculate_lag_bubble_stats |
This procedure calculates the maximum and minimum radius of each bubble.
| subroutine m_bubbles_el::s_compute_bubble_el_dynamics | ( | type(scalar_field), dimension(sys_size), intent(inout) | q_prim_vf, |
| integer, intent(in) | stage ) |
Contains the bubble dynamics subroutines.
| q_cons_vf | Conservative variables |
| q_prim_vf | Primitive variables |
| rhs_vf | Calculated change of conservative variables |
| t_step | Current time step |
| stage | Current stage in the time-stepper algorithm |
| subroutine m_bubbles_el::s_compute_bubbles_el_source | ( | type(scalar_field), dimension(sys_size), intent(inout) | q_cons_vf, |
| type(scalar_field), dimension(sys_size), intent(inout) | q_prim_vf, | ||
| type(scalar_field), dimension(sys_size), intent(inout) | 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.
| q_cons_vf | Conservative variables |
| q_prim_vf | Conservative variables |
| rhs_vf | Time derivative of the conservative variables |
| subroutine m_bubbles_el::s_compute_cson_from_pinf | ( | type(scalar_field), dimension(sys_size), intent(in) | q_prim_vf, |
| real(wp), intent(in) | pinf, | ||
| integer, dimension(3), intent(in) | cell, | ||
| real(wp), intent(in) | rhol, | ||
| real(wp), intent(in) | gamma, | ||
| real(wp), intent(in) | pi_inf, | ||
| real(wp), intent(out) | cson ) |
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 |
| impure subroutine m_bubbles_el::s_finalize_lagrangian_solver |
The purpose of this subroutine is to deallocate variables.
| subroutine m_bubbles_el::s_get_pinf | ( | integer, intent(in) | bub_id, |
| type(scalar_field), dimension(sys_size), intent(in) | q_prim_vf, | ||
| integer, intent(in) | ptype, | ||
| real(wp), intent(out) | f_pinfl, | ||
| integer, dimension(3), intent(out) | cell, | ||
| real(wp), intent(out), optional | preterm1, | ||
| real(wp), intent(out), optional | term2, | ||
| real(wp), intent(out), optional | romega ) |
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 |
| subroutine m_bubbles_el::s_gradient_dir | ( | real(stp), dimension(idwbuff(1)%beg:, idwbuff(2)%beg:, idwbuff(3)%beg:), intent(inout) | q, |
| real(stp), dimension(idwbuff(1)%beg:, idwbuff(2)%beg:, idwbuff(3)%beg:), intent(inout) | dq, | ||
| integer, intent(in) | 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.
| q | Input scalar field |
| dq | Output gradient of q |
| dir | Gradient spatial direction |
| impure subroutine m_bubbles_el::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 |
| subroutine m_bubbles_el::s_locate_cell | ( | real(wp), dimension(3), intent(in) | pos, |
| integer, dimension(3), intent(inout) | cell, | ||
| real(wp), dimension(3), intent(out) | scoord ) |
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 |
| impure subroutine m_bubbles_el::s_read_input_bubbles | ( | type(scalar_field), dimension(sys_size), intent(inout) | q_cons_vf | ) |
The purpose of this procedure is to obtain the initial bubbles' information.
| q_cons_vf | Conservative variables |
| impure subroutine m_bubbles_el::s_remove_lag_bubble | ( | integer, intent(in) | bub_id | ) |
The purpose of this subroutine is to remove one specific particle if dt is too small.
| bub_id | Particle id |
| impure subroutine m_bubbles_el::s_restart_bubbles | ( | integer, intent(inout) | bub_id, |
| integer, intent(inout) | save_count ) |
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 |
| subroutine m_bubbles_el::s_smear_voidfraction |
The purpose of this subroutine is to smear the effect of the bubbles in the Eulerian framework.
| impure subroutine m_bubbles_el::s_start_lagrange_inputs |
The purpose of this procedure is to start lagrange bubble parameters applying nondimensionalization if needed.
| impure subroutine m_bubbles_el::s_transfer_data_to_tmp |
This subroutine transfer data into the temporal variables.
| impure subroutine m_bubbles_el::s_update_lagrange_tdv_rk | ( | integer, intent(in) | 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.
| stage | Current tvd RK stage |
| impure subroutine m_bubbles_el::s_write_lag_bubble_stats |
Subroutine that writes the maximum and minimum radius of each bubble.
| impure subroutine m_bubbles_el::s_write_lag_particles | ( | real(wp), intent(in) | qtime | ) |
Subroutine that writes on each time step the changes of the lagrangian bubbles.
| q_time | Current time |
| impure subroutine m_bubbles_el::s_write_restart_lag_bubbles | ( | integer, intent(in) | t_step | ) |
Subroutine that writes the restarting files for the particles in the lagrangian solver.
| t_step | Current time step |
| impure subroutine m_bubbles_el::s_write_void_evol | ( | real(wp), intent(in) | 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.
| q_time | Current time |
| real(wp), dimension(:), allocatable m_bubbles_el::bub_dphidt |
subgrid velocity potential (Maeda & Colonius, 2018)
| real(wp), dimension(:), allocatable m_bubbles_el::bub_r0 |
Initial bubble radius.
| real(wp), dimension(:), allocatable m_bubbles_el::gas_betac |
massflux model (Preston et al., 2007)
| real(wp), dimension(:), allocatable m_bubbles_el::gas_betat |
heatflux model (Preston et al., 2007)
| real(wp), dimension(:, :), allocatable m_bubbles_el::gas_dmvdt |
Time derivative of the vapor mass in the bubble.
| real(wp), dimension(:, :), allocatable m_bubbles_el::gas_dpdt |
Time derivative of gas pressure.
| real(wp), dimension(:), allocatable m_bubbles_el::gas_mg |
Bubble's gas mass.
| real(wp), dimension(:, :), allocatable m_bubbles_el::gas_mv |
Vapor mass in the bubble.
| real(wp), dimension(:, :), allocatable m_bubbles_el::gas_p |
Pressure in the bubble.
| real(wp), dimension(:, :), allocatable m_bubbles_el::intfc_draddt |
Time derivative of bubble's radius.
| real(wp), dimension(:, :), allocatable m_bubbles_el::intfc_dveldt |
Time derivative of bubble's interface velocity.
| real(wp), dimension(:, :), allocatable m_bubbles_el::intfc_rad |
Bubble radius.
| real(wp), dimension(:, :), allocatable m_bubbles_el::intfc_vel |
Velocity of the bubble interface.
| integer, dimension(:, :), allocatable m_bubbles_el::lag_id |
Global and local IDs.
|
private |
Number of time stages in the time-stepping scheme.
| real(wp), dimension(:, :, :), allocatable m_bubbles_el::mtn_dposdt |
Time derivative of the bubble's position.
| real(wp), dimension(:, :, :), allocatable m_bubbles_el::mtn_dveldt |
Time derivative of the bubble's velocity.
| real(wp), dimension(:, :, :), allocatable m_bubbles_el::mtn_pos |
Bubble's position.
| real(wp), dimension(:, :, :), allocatable m_bubbles_el::mtn_posprev |
Bubble's previous position.
| real(wp), dimension(:, :, :), allocatable m_bubbles_el::mtn_s |
Bubble's computational cell position in real format.
| real(wp), dimension(:, :, :), allocatable m_bubbles_el::mtn_vel |
Bubble's velocity.
| integer m_bubbles_el::nbubs |
Number of bubbles in the local domain.
| type(scalar_field), dimension(:), allocatable m_bubbles_el::q_beta |
| integer m_bubbles_el::q_beta_idx |
Size of the q_beta vector field.
| real(wp) m_bubbles_el::rmax_glb |
| real(wp), dimension(:), allocatable m_bubbles_el::rmax_stats |
Maximum radius.
| real(wp) m_bubbles_el::rmin_glb |
Maximum and minimum bubbe size in the local domain Projection of the lagrangian particles in the Eulerian framework.
| real(wp), dimension(:), allocatable m_bubbles_el::rmin_stats |
Minimum radius.