MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_derived_variables Module Reference

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)
 Derive the specific heat ratio from the specific heat ratio function gamma_sf. The latter is stored in the derived flow quantity storage variable, q_sf.
subroutine, public s_derive_liquid_stiffness (q_sf)
 Compute the liquid stiffness from the specific heat ratio function gamma_sf and the liquid stiffness function 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)
 Compute the speed of sound from the primitive variables, density, 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)
 Derive 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, public s_derive_vorticity_component (i, q_prim_vf, q_sf)
 Compute the specified component of the vorticity from 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)
 Compute the Q_M criterion from the primitive variables. 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)
 Compute 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)
 Compute 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

type(fd_context), public fd
 Finite-difference state: density gradient magnitude and centered FD coefficients in x-, y-, and z-directions.

Detailed Description

Computes derived flow quantities (sound speed, vorticity, Schlieren, etc.) from conservative and primitive variables.

Function/Subroutine Documentation

◆ s_derive_flux_limiter()

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 )

Derive 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.

Definition at line 128 of file m_derived_variables.fpp.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ s_derive_liquid_stiffness()

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)

Compute the liquid stiffness from the specific heat ratio function gamma_sf and the liquid stiffness function 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.

Definition at line 73 of file m_derived_variables.fpp.f90.

Here is the caller graph for this function:

◆ s_derive_liutex()

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 )

Compute the Liutex vector and its magnitude based on Xu et al. (2019).

Parameters
[out]liutex_magLiutex magnitude
[out]liutex_axisLiutex rigid rotation axis

Definition at line 322 of file m_derived_variables.fpp.f90.

Here is the caller graph for this function:

◆ s_derive_numerical_schlieren_function()

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 )

Compute the values of the numerical Schlieren function, which are subsequently stored in the derived flow quantity storage variable, q_sf.

Definition at line 427 of file m_derived_variables.fpp.f90.

Here is the caller graph for this function:

◆ s_derive_qm()

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 )

Compute the Q_M criterion from the primitive variables. The Q_M function, which are subsequently stored in the derived flow quantity storage variable, q_sf.

Definition at line 264 of file m_derived_variables.fpp.f90.

Here is the caller graph for this function:

◆ s_derive_sound_speed()

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 )

Compute the speed of sound from the primitive variables, density, 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.

Definition at line 92 of file m_derived_variables.fpp.f90.

◆ s_derive_specific_heat_ratio()

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)

Derive the specific heat ratio from the specific heat ratio function gamma_sf. The latter is stored in the derived flow quantity storage variable, q_sf.

Definition at line 54 of file m_derived_variables.fpp.f90.

Here is the caller graph for this function:

◆ s_derive_vorticity_component()

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 )

Compute the specified component of the vorticity from 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.

Definition at line 199 of file m_derived_variables.fpp.f90.

Here is the caller graph for this function:

◆ s_finalize_derived_variables_module()

impure subroutine, public m_derived_variables::s_finalize_derived_variables_module

Deallocation procedures for the module.

Definition at line 525 of file m_derived_variables.fpp.f90.

◆ s_initialize_derived_variables_module()

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 29 of file m_derived_variables.fpp.f90.

Here is the caller graph for this function:

Variable Documentation

◆ fd

type(fd_context), public m_derived_variables::fd

Finite-difference state: density gradient magnitude and centered FD coefficients in x-, y-, and z-directions.

Definition at line 24 of file m_derived_variables.fpp.f90.