|
MFC
Exascale flow solver
|
Computes derived flow quantities (sound speed, vorticity, Schlieren, etc.) from conservative and primitive variables. More...
Functions/Subroutines | |
| impure subroutine, public | s_initialize_derived_variables_module |
| Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the module. | |
| subroutine, public | s_derive_specific_heat_ratio (q_sf) |
| This subroutine receives as input the specific heat ratio function, gamma_sf, and derives from it the specific heat ratio. The latter is stored in the derived flow quantity storage variable, q_sf. | |
| subroutine, public | s_derive_liquid_stiffness (q_sf) |
| This subroutine admits as inputs the specific heat ratio function and the liquid stiffness function, gamma_sf and pi_inf_sf, respectively. These are used to calculate the values of the liquid stiffness, which are stored in the derived flow quantity storage variable, q_sf. | |
| subroutine, public | s_derive_sound_speed (q_prim_vf, q_sf) |
| This subroutine admits as inputs the primitive variables, the density, the specific heat ratio function and liquid stiffness function. It then computes from those variables the values of the speed of sound, which are stored in the derived flow quantity storage variable, q_sf. | |
| subroutine, public | s_derive_flux_limiter (i, q_prim_vf, q_sf) |
| This subroutine derives the flux_limiter at cell boundary i+1/2. This is an approximation because the velocity used to determine the upwind direction is the velocity at the cell center i instead of the contact velocity at the cell boundary from the Riemann solver. | |
| subroutine | s_solve_linear_system (a, b, sol, ndim) |
| Computes the solution to the linear system Ax=b w/ sol = x. | |
| subroutine, public | s_derive_vorticity_component (i, q_prim_vf, q_sf) |
| This subroutine receives as inputs the indicator of the component of the vorticity that should be outputted and the primitive variables. From those inputs, it proceeds to calculate values of the desired vorticity component, which are subsequently stored in derived flow quantity storage variable, q_sf. | |
| subroutine, public | s_derive_qm (q_prim_vf, q_sf) |
| This subroutine gets as inputs the primitive variables. From those inputs, it proceeds to calculate the value of the Q_M function, which are subsequently stored in the derived flow quantity storage variable, q_sf. | |
| impure subroutine, public | s_derive_liutex (q_prim_vf, liutex_mag, liutex_axis) |
| This subroutine gets as inputs the primitive variables. From those inputs, it proceeds to calculate the Liutex vector and its magnitude based on Xu et al. (2019). | |
| impure subroutine, public | s_derive_numerical_schlieren_function (q_cons_vf, q_sf) |
| This subroutine gets as inputs the conservative variables and density. From those inputs, it proceeds to calculate the values of the numerical Schlieren function, which are subsequently stored in the derived flow quantity storage variable, q_sf. | |
| impure subroutine, public | s_finalize_derived_variables_module |
| Deallocation procedures for the module. | |
Variables | |
| real(wp), dimension(:, :, :), allocatable | gm_rho_sf |
| Gradient magnitude (gm) of the density for each cell of the computational sub-domain. This variable is employed in the calculation of the numerical Schlieren function. | |
| integer, private | flg |
| Flagging (flg) variable used to annotate the dimensionality of the dataset that is undergoing the post-process. A flag value of 1 indicates that the dataset is 3D, while a flag value of 0 indicates that it is not. This flg variable is necessary to avoid cycling through the third dimension of the flow variable(s) when the simulation is not 3D and the size of the buffer is non-zero. Note that a similar procedure does not have to be applied to the second dimension since in 1D, the buffer size is always zero. | |
Finite-difference (fd) coefficients in x-, y- and z-coordinate directions. | |
Note that because sufficient boundary information is available for all the active coordinate directions, the centered family of the finite-difference schemes is used. | |
| real(wp), dimension(:, :), allocatable, public | fd_coeff_x |
| real(wp), dimension(:, :), allocatable, public | fd_coeff_y |
| real(wp), dimension(:, :), allocatable, public | fd_coeff_z |
Computes derived flow quantities (sound speed, vorticity, Schlieren, etc.) from conservative and primitive variables.
| subroutine, public m_derived_variables::s_derive_flux_limiter | ( | integer, intent(in) | i, |
| type(scalar_field), dimension(sys_size), intent(in) | q_prim_vf, | ||
| real(wp), dimension(-offset_x%beg:m + offset_x%end, -offset_y%beg:n + offset_y%end, -offset_z%beg:p + offset_z%end), intent(inout) | q_sf ) |
This subroutine derives the flux_limiter at cell boundary i+1/2. This is an approximation because the velocity used to determine the upwind direction is the velocity at the cell center i instead of the contact velocity at the cell boundary from the Riemann solver.
| i | Component indicator |
| q_prim_vf | Primitive variables |
| q_sf | Flux limiter |
Definition at line 231 of file m_derived_variables.fpp.f90.
| subroutine, public m_derived_variables::s_derive_liquid_stiffness | ( | real(wp), dimension(-offset_x%beg:m + offset_x%end, -offset_y%beg:n + offset_y%end, -offset_z%beg:p + offset_z%end), intent(inout) | q_sf | ) |
This subroutine admits as inputs the specific heat ratio function and the liquid stiffness function, gamma_sf and pi_inf_sf, respectively. These are used to calculate the values of the liquid stiffness, which are stored in the derived flow quantity storage variable, q_sf.
| q_sf | Liquid stiffness |
Definition at line 145 of file m_derived_variables.fpp.f90.
| impure subroutine, public m_derived_variables::s_derive_liutex | ( | type(scalar_field), dimension(sys_size), intent(in) | q_prim_vf, |
| real(wp), dimension(-offset_x%beg:m + offset_x%end, -offset_y%beg:n + offset_y%end, -offset_z%beg:p + offset_z%end), intent(out) | liutex_mag, | ||
| real(wp), dimension(-offset_x%beg:m + offset_x%end, -offset_y%beg:n + offset_y%end, -offset_z%beg:p + offset_z%end, nm), intent(out) | liutex_axis ) |
This subroutine gets as inputs the primitive variables. From those inputs, it proceeds to calculate the Liutex vector and its magnitude based on Xu et al. (2019).
| q_prim_vf | Primitive variables | |
| [out] | liutex_mag | Liutex magnitude |
| [out] | liutex_axis | Liutex rigid rotation axis |
Definition at line 561 of file m_derived_variables.fpp.f90.
| impure subroutine, public m_derived_variables::s_derive_numerical_schlieren_function | ( | type(scalar_field), dimension(sys_size), intent(in) | q_cons_vf, |
| real(wp), dimension(-offset_x%beg:m + offset_x%end, -offset_y%beg:n + offset_y%end, -offset_z%beg:p + offset_z%end), intent(inout) | q_sf ) |
This subroutine gets as inputs the conservative variables and density. From those inputs, it proceeds to calculate the values of the numerical Schlieren function, which are subsequently stored in the derived flow quantity storage variable, q_sf.
| q_cons_vf | Conservative variables |
| q_sf | Numerical Schlieren function |
Definition at line 691 of file m_derived_variables.fpp.f90.
| subroutine, public m_derived_variables::s_derive_qm | ( | type(scalar_field), dimension(sys_size), intent(in) | q_prim_vf, |
| real(wp), dimension(-offset_x%beg:m + offset_x%end, -offset_y%beg:n + offset_y%end, -offset_z%beg:p + offset_z%end), intent(inout) | q_sf ) |
This subroutine gets as inputs the primitive variables. From those inputs, it proceeds to calculate the value of the Q_M function, which are subsequently stored in the derived flow quantity storage variable, q_sf.
| q_prim_vf | Primitive variables |
| q_sf | Q_M |
Definition at line 477 of file m_derived_variables.fpp.f90.
| subroutine, public m_derived_variables::s_derive_sound_speed | ( | type(scalar_field), dimension(sys_size), intent(in) | q_prim_vf, |
| real(wp), dimension(-offset_x%beg:m + offset_x%end, -offset_y%beg:n + offset_y%end, -offset_z%beg:p + offset_z%end), intent(inout) | q_sf ) |
This subroutine admits as inputs the primitive variables, the density, the specific heat ratio function and liquid stiffness function. It then computes from those variables the values of the speed of sound, which are stored in the derived flow quantity storage variable, q_sf.
| q_prim_vf | Primitive variables |
| q_sf | Speed of sound |
Definition at line 174 of file m_derived_variables.fpp.f90.
| subroutine, public m_derived_variables::s_derive_specific_heat_ratio | ( | real(wp), dimension(-offset_x%beg:m + offset_x%end, -offset_y%beg:n + offset_y%end, -offset_z%beg:p + offset_z%end), intent(inout) | q_sf | ) |
This subroutine receives as input the specific heat ratio function, gamma_sf, and derives from it the specific heat ratio. The latter is stored in the derived flow quantity storage variable, q_sf.
| q_sf | Specific heat ratio |
Definition at line 118 of file m_derived_variables.fpp.f90.
| subroutine, public m_derived_variables::s_derive_vorticity_component | ( | integer, intent(in) | i, |
| type(scalar_field), dimension(sys_size), intent(in) | q_prim_vf, | ||
| real(wp), dimension(-offset_x%beg:m + offset_x%end, -offset_y%beg:n + offset_y%end, -offset_z%beg:p + offset_z%end), intent(inout) | q_sf ) |
This subroutine receives as inputs the indicator of the component of the vorticity that should be outputted and the primitive variables. From those inputs, it proceeds to calculate values of the desired vorticity component, which are subsequently stored in derived flow quantity storage variable, q_sf.
| i | Vorticity component indicator |
| q_prim_vf | Primitive variables |
| q_sf | Vorticity component |
Definition at line 375 of file m_derived_variables.fpp.f90.
| impure subroutine, public m_derived_variables::s_finalize_derived_variables_module |
Deallocation procedures for the module.
Definition at line 816 of file m_derived_variables.fpp.f90.
| impure subroutine, public m_derived_variables::s_initialize_derived_variables_module |
Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the module.
Definition at line 62 of file m_derived_variables.fpp.f90.
| subroutine m_derived_variables::s_solve_linear_system | ( | real(wp), dimension(ndim, ndim), intent(inout) | a, |
| real(wp), dimension(ndim), intent(inout) | b, | ||
| real(wp), dimension(ndim), intent(out) | sol, | ||
| integer, intent(in) | ndim ) |
Computes the solution to the linear system Ax=b w/ sol = x.
| A | Input matrix |
| b | right-hane-side |
| sol | Solution |
| ndim | Problem size |
Definition at line 325 of file m_derived_variables.fpp.f90.
| real(wp), dimension(:, :), allocatable, public m_derived_variables::fd_coeff_x |
Definition at line 44 of file m_derived_variables.fpp.f90.
| real(wp), dimension(:, :), allocatable, public m_derived_variables::fd_coeff_y |
Definition at line 45 of file m_derived_variables.fpp.f90.
| real(wp), dimension(:, :), allocatable, public m_derived_variables::fd_coeff_z |
Definition at line 46 of file m_derived_variables.fpp.f90.
|
private |
Flagging (flg) variable used to annotate the dimensionality of the dataset that is undergoing the post-process. A flag value of 1 indicates that the dataset is 3D, while a flag value of 0 indicates that it is not. This flg variable is necessary to avoid cycling through the third dimension of the flow variable(s) when the simulation is not 3D and the size of the buffer is non-zero. Note that a similar procedure does not have to be applied to the second dimension since in 1D, the buffer size is always zero.
Definition at line 49 of file m_derived_variables.fpp.f90.
| real(wp), dimension(:, :, :), allocatable m_derived_variables::gm_rho_sf |
Gradient magnitude (gm) of the density for each cell of the computational sub-domain. This variable is employed in the calculation of the numerical Schlieren function.
Definition at line 34 of file m_derived_variables.fpp.f90.