|
MFC
Exascale flow solver
|
Conservative-to-primitive variable conversion, mixture property evaluation, and pressure computation. More...
Functions/Subroutines | |
| subroutine, public | s_convert_to_mixture_variables (q_vf, i, j, k, rho, gamma, pi_inf, qv, re_k, g_k, g) |
| Dispatch to the s_convert_mixture_to_mixture_variables and s_convert_species_to_mixture_variables subroutines. Replaces a procedure pointer. | |
| subroutine, public | s_compute_pressure (energy, alf, dyn_p, pi_inf, gamma, rho, qv, rhoyks, pres, t, stress, mom, g, pres_mag) |
| Compute the pressure from the appropriate equation of state. | |
| subroutine, public | s_convert_mixture_to_mixture_variables (q_vf, i, j, k, rho, gamma, pi_inf, qv) |
| Convert mixture variables to density, gamma, pi_inf, and qv for the gamma/pi_inf model. Given conservative or primitive variables, transfers the density, specific heat ratio function and the liquid stiffness function from q_vf to rho, gamma and pi_inf. | |
| subroutine, public | s_convert_species_to_mixture_variables (q_vf, k, l, r, rho, gamma, pi_inf, qv, re_k, g_k, g) |
| Convert species volume fractions and partial densities to mixture density, gamma, pi_inf, and qv. Given conservative or primitive variables, computes the density, the specific heat ratio function and the liquid stiffness function from q_vf and stores the results into rho, gamma and pi_inf. | |
| subroutine, public | s_convert_species_to_mixture_variables_acc (rho_k, gamma_k, pi_inf_k, qv_k, alpha_k, alpha_rho_k, re_k, g_k, g) |
| GPU-accelerated conversion of species volume fractions and partial densities to mixture density, gamma, pi_inf, and qv. | |
| impure subroutine, public | s_initialize_variables_conversion_module |
| Initialize the variables conversion module. | |
| subroutine, public | s_initialize_mv (qk_cons_vf, mv) |
| Initialize bubble mass-vapor values at quadrature nodes from the conserved moment statistics. | |
| subroutine, public | s_initialize_pb (qk_cons_vf, mv, pb) |
| Initialize bubble internal pressures at quadrature nodes using isothermal relations from the Preston model. | |
| subroutine, public | s_convert_conservative_to_primitive_variables (qk_cons_vf, q_t_sf, qk_prim_vf, ibounds) |
| Convert conserved variables (rho*alpha, rho*u, E, alpha) to primitives (rho, u, p, alpha). Conversion depends on model_eqns: each model has different variable sets and EOS. | |
| impure subroutine, public | s_convert_primitive_to_conservative_variables (q_prim_vf, q_cons_vf) |
| Convert primitives (rho, u, p, alpha) to conserved variables (rho*alpha, rho*u, E, alpha). | |
| subroutine, public | s_convert_primitive_to_flux_variables (qk_prim_vf, fk_vf, fk_src_vf, is1, is2, is3, s2b, s3b) |
| Convert primitive variables to Eulerian flux variables. | |
| subroutine, public | s_compute_species_fraction (q_vf, k, l, r, alpha_rho_k, alpha_k) |
| Compute partial densities and volume fractions. | |
| impure subroutine | s_finalize_variables_conversion_module () |
| Deallocate fluid property arrays and post-processing fields allocated during module initialization. | |
| subroutine | s_compute_speed_of_sound (pres, rho, gamma, pi_inf, h, adv, vel_sum, c_c, c, qv) |
| Compute the speed of sound from thermodynamic state variables, supporting multiple equation-of-state models. | |
| subroutine | s_compute_fast_magnetosonic_speed (rho, c, b, norm, c_fast, h) |
| Compute the fast magnetosonic wave speed from the sound speed, density, and magnetic field components. | |
Variables | |
| real(wp), dimension(:), allocatable, public | gammas |
| real(wp), dimension(:), allocatable, public | gs_min |
| real(wp), dimension(:), allocatable, public | pi_infs |
| real(wp), dimension(:), allocatable, public | ps_inf |
| real(wp), dimension(:), allocatable, public | cvs |
| real(wp), dimension(:), allocatable, public | qvs |
| real(wp), dimension(:), allocatable, public | qvps |
| real(wp), dimension(:), allocatable | gs_vc |
| integer, dimension(:), allocatable | bubrs_vc |
| real(wp), dimension(:,:), allocatable | res_vc |
| integer | is1b |
| integer | is2b |
| integer | is3b |
| integer | is1e |
| integer | is2e |
| integer | is3e |
| real(wp), dimension(:,:,:), allocatable, public | rho_sf |
| Scalar density function. | |
| real(wp), dimension(:,:,:), allocatable, public | gamma_sf |
| Scalar sp. heat ratio function. | |
| real(wp), dimension(:,:,:), allocatable, public | pi_inf_sf |
| Scalar liquid stiffness function. | |
| real(wp), dimension(:,:,:), allocatable, public | qv_sf |
| Scalar liquid energy reference function. | |
Conservative-to-primitive variable conversion, mixture property evaluation, and pressure computation.
|
private |
Compute the fast magnetosonic wave speed from the sound speed, density, and magnetic field components.
| [in] | h | only used for relativity |
Definition at line 2625 of file m_variables_conversion.fpp.f90.
| subroutine, public m_variables_conversion::s_compute_pressure | ( | real(stp), intent(in) | energy, |
| real(stp), intent(in) | alf, | ||
| real(wp), intent(in) | dyn_p, | ||
| real(wp), intent(in) | pi_inf, | ||
| real(wp), intent(in) | gamma, | ||
| real(wp), intent(in) | rho, | ||
| real(wp), intent(in) | qv, | ||
| real(wp), dimension(1:num_species), intent(in) | rhoyks, | ||
| real(wp), intent(out) | pres, | ||
| real(wp), intent(inout) | t, | ||
| real(stp), intent(in), optional | stress, | ||
| real(stp), intent(in), optional | mom, | ||
| real(wp), intent(in), optional | g, | ||
| real(wp), intent(in), optional | pres_mag ) |
Compute the pressure from the appropriate equation of state.
Definition at line 430 of file m_variables_conversion.fpp.f90.
| subroutine, public m_variables_conversion::s_compute_species_fraction | ( | type(scalar_field), dimension(sys_size), intent(in) | q_vf, |
| integer, intent(in) | k, | ||
| integer, intent(in) | l, | ||
| integer, intent(in) | r, | ||
| real(wp), dimension(num_fluids), intent(out) | alpha_rho_k, | ||
| real(wp), dimension(num_fluids), intent(out) | alpha_k ) |
Compute partial densities and volume fractions.
Definition at line 2315 of file m_variables_conversion.fpp.f90.
|
private |
Compute the speed of sound from thermodynamic state variables, supporting multiple equation-of-state models.
Definition at line 2541 of file m_variables_conversion.fpp.f90.
| subroutine, public m_variables_conversion::s_convert_conservative_to_primitive_variables | ( | type(scalar_field), dimension(sys_size), intent(in) | qk_cons_vf, |
| type(scalar_field), intent(inout) | q_t_sf, | ||
| type(scalar_field), dimension(sys_size), intent(inout) | qk_prim_vf, | ||
| type(int_bounds_info), dimension(1:3), intent(in) | ibounds ) |
Convert conserved variables (rho*alpha, rho*u, E, alpha) to primitives (rho, u, p, alpha). Conversion depends on model_eqns: each model has different variable sets and EOS.
Definition at line 1245 of file m_variables_conversion.fpp.f90.
| subroutine, public m_variables_conversion::s_convert_mixture_to_mixture_variables | ( | type(scalar_field), dimension(sys_size), intent(in) | q_vf, |
| integer, intent(in) | i, | ||
| integer, intent(in) | j, | ||
| integer, intent(in) | k, | ||
| real(wp), intent(out), target | rho, | ||
| real(wp), intent(out), target | gamma, | ||
| real(wp), intent(out), target | pi_inf, | ||
| real(wp), intent(out), target | qv ) |
Convert mixture variables to density, gamma, pi_inf, and qv for the gamma/pi_inf model. Given conservative or primitive variables, transfers the density, specific heat ratio function and the liquid stiffness function from q_vf to rho, gamma and pi_inf.
Definition at line 523 of file m_variables_conversion.fpp.f90.
| impure subroutine, public m_variables_conversion::s_convert_primitive_to_conservative_variables | ( | type(scalar_field), dimension(sys_size), intent(in) | q_prim_vf, |
| type(scalar_field), dimension(sys_size), intent(inout) | q_cons_vf ) |
Convert primitives (rho, u, p, alpha) to conserved variables (rho*alpha, rho*u, E, alpha).
Definition at line 1782 of file m_variables_conversion.fpp.f90.
| subroutine, public m_variables_conversion::s_convert_primitive_to_flux_variables | ( | real(wp), dimension(0:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(in) | qk_prim_vf, |
| real(wp), dimension(0:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) | fk_vf, | ||
| real(wp), dimension(0:,idwbuff(2)%beg:,idwbuff(3)%beg:,advxb:), intent(inout) | fk_src_vf, | ||
| type(int_bounds_info), intent(in) | is1, | ||
| type(int_bounds_info), intent(in) | is2, | ||
| type(int_bounds_info), intent(in) | is3, | ||
| integer, intent(in) | s2b, | ||
| integer, intent(in) | s3b ) |
Convert primitive variables to Eulerian flux variables.
Definition at line 2026 of file m_variables_conversion.fpp.f90.
| subroutine, public m_variables_conversion::s_convert_species_to_mixture_variables | ( | type(scalar_field), dimension(sys_size), intent(in) | q_vf, |
| integer, intent(in) | k, | ||
| integer, intent(in) | l, | ||
| integer, intent(in) | r, | ||
| real(wp), intent(out), target | rho, | ||
| real(wp), intent(out), target | gamma, | ||
| real(wp), intent(out), target | pi_inf, | ||
| real(wp), intent(out), target | qv, | ||
| real(wp), dimension(2), intent(out), optional | re_k, | ||
| real(wp), intent(out), optional | g_k, | ||
| real(wp), dimension(num_fluids), intent(in), optional | g ) |
Convert species volume fractions and partial densities to mixture density, gamma, pi_inf, and qv. Given conservative or primitive variables, computes the density, the specific heat ratio function and the liquid stiffness function from q_vf and stores the results into rho, gamma and pi_inf.
Definition at line 552 of file m_variables_conversion.fpp.f90.
| subroutine, public m_variables_conversion::s_convert_species_to_mixture_variables_acc | ( | real(wp), intent(out) | rho_k, |
| real(wp), intent(out) | gamma_k, | ||
| real(wp), intent(out) | pi_inf_k, | ||
| real(wp), intent(out) | qv_k, | ||
| real(wp), dimension(num_fluids), intent(inout) | alpha_k, | ||
| real(wp), dimension(num_fluids), intent(inout) | alpha_rho_k, | ||
| real(wp), dimension(2), intent(out) | re_k, | ||
| real(wp), intent(out), optional | g_k, | ||
| real(wp), dimension(num_fluids), intent(in), optional | g ) |
GPU-accelerated conversion of species volume fractions and partial densities to mixture density, gamma, pi_inf, and qv.
Definition at line 620 of file m_variables_conversion.fpp.f90.
| subroutine, public m_variables_conversion::s_convert_to_mixture_variables | ( | type(scalar_field), dimension(sys_size), intent(in) | q_vf, |
| integer, intent(in) | i, | ||
| integer, intent(in) | j, | ||
| integer, intent(in) | k, | ||
| real(wp), intent(out), target | rho, | ||
| real(wp), intent(out), target | gamma, | ||
| real(wp), intent(out), target | pi_inf, | ||
| real(wp), intent(out), target | qv, | ||
| real(wp), dimension(2), intent(out), optional | re_k, | ||
| real(wp), intent(out), optional | g_k, | ||
| real(wp), dimension(num_fluids), intent(in), optional | g ) |
Dispatch to the s_convert_mixture_to_mixture_variables and s_convert_species_to_mixture_variables subroutines. Replaces a procedure pointer.
Definition at line 412 of file m_variables_conversion.fpp.f90.
|
private |
Deallocate fluid property arrays and post-processing fields allocated during module initialization.
Definition at line 2398 of file m_variables_conversion.fpp.f90.
| subroutine, public m_variables_conversion::s_initialize_mv | ( | type(scalar_field), dimension(sys_size), intent(in) | qk_cons_vf, |
| real(stp), dimension(idwint(1)%beg:,idwint(2)%beg:,idwint(3)%beg:,1:,1:), intent(inout) | mv ) |
Initialize bubble mass-vapor values at quadrature nodes from the conserved moment statistics.
Definition at line 1160 of file m_variables_conversion.fpp.f90.
| subroutine, public m_variables_conversion::s_initialize_pb | ( | type(scalar_field), dimension(sys_size), intent(in) | qk_cons_vf, |
| real(stp), dimension(idwint(1)%beg:,idwint(2)%beg:,idwint(3)%beg:,1:,1:), intent(in) | mv, | ||
| real(stp), dimension(idwint(1)%beg:,idwint(2)%beg:,idwint(3)%beg:,1:,1:), intent(inout) | pb ) |
Initialize bubble internal pressures at quadrature nodes using isothermal relations from the Preston model.
Definition at line 1199 of file m_variables_conversion.fpp.f90.
| impure subroutine, public m_variables_conversion::s_initialize_variables_conversion_module |
Initialize the variables conversion module.
Definition at line 720 of file m_variables_conversion.fpp.f90.
|
private |
Definition at line 376 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:), allocatable, public m_variables_conversion::cvs |
Definition at line 361 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:,:,:), allocatable, public m_variables_conversion::gamma_sf |
Scalar sp. heat ratio function.
Definition at line 404 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:), allocatable, public m_variables_conversion::gammas |
Definition at line 361 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:), allocatable, public m_variables_conversion::gs_min |
Definition at line 361 of file m_variables_conversion.fpp.f90.
|
private |
Definition at line 375 of file m_variables_conversion.fpp.f90.
|
private |
Definition at line 390 of file m_variables_conversion.fpp.f90.
|
private |
Definition at line 390 of file m_variables_conversion.fpp.f90.
|
private |
Definition at line 390 of file m_variables_conversion.fpp.f90.
|
private |
Definition at line 390 of file m_variables_conversion.fpp.f90.
|
private |
Definition at line 390 of file m_variables_conversion.fpp.f90.
|
private |
Definition at line 390 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:,:,:), allocatable, public m_variables_conversion::pi_inf_sf |
Scalar liquid stiffness function.
Definition at line 405 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:), allocatable, public m_variables_conversion::pi_infs |
Definition at line 361 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:), allocatable, public m_variables_conversion::ps_inf |
Definition at line 361 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:,:,:), allocatable, public m_variables_conversion::qv_sf |
Scalar liquid energy reference function.
Definition at line 406 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:), allocatable, public m_variables_conversion::qvps |
Definition at line 361 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:), allocatable, public m_variables_conversion::qvs |
Definition at line 361 of file m_variables_conversion.fpp.f90.
|
private |
Definition at line 377 of file m_variables_conversion.fpp.f90.
| real(wp), dimension(:,:,:), allocatable, public m_variables_conversion::rho_sf |
Scalar density function.
Definition at line 403 of file m_variables_conversion.fpp.f90.