init_infinity_values Subroutine

private subroutine init_infinity_values(scheme, flow)

Set the values of the infinity variables "qp_inf"

Arguments

Type IntentOptional AttributesName
type(schemetype), intent(in) :: scheme

finite-volume Schemes: turbulence, transition model, etc

type(flowtype), intent(inout) :: flow

Information about fluid flow: freestream-speed, ref-viscosity,etc.


Calls

proc~~init_infinity_values~~CallsGraph proc~init_infinity_values init_infinity_values proc~sound_speed_inf sound_speed_inf proc~init_infinity_values->proc~sound_speed_inf debugcall debugcall proc~init_infinity_values->debugcall

Called by

proc~~init_infinity_values~~CalledByGraph proc~init_infinity_values init_infinity_values proc~setup_state setup_state proc~setup_state->proc~init_infinity_values proc~setup_solver setup_solver proc~setup_solver->proc~setup_state proc~start_run start_run proc~start_run->proc~setup_solver program~main main program~main->proc~start_run

Contents

Source Code


Source Code

        subroutine init_infinity_values(scheme, flow)
            !< Set the values of the infinity variables "qp_inf"
            !-----------------------------------------------------------

            implicit none
            type(schemetype), intent(in) :: scheme
            !< finite-volume Schemes: turbulence, transition model, etc
            type(flowtype), intent(inout) :: flow
            !< Information about fluid flow: freestream-speed, ref-viscosity,etc.
            
            DebugCall("init_infinity_values")

            flow%vel_mag = sqrt(flow%x_speed_inf**2 + flow%y_speed_inf**2 + flow%z_speed_inf**2)
            flow%MInf    = flow%vel_mag/sqrt(flow%gm*flow%pressure_inf/flow%density_inf)
            flow%Reynolds_number = flow%density_inf*flow%vel_mag*1.0/flow%mu_ref
            flow%Turb_intensity_inf = flow%tu_inf/100

            select case (trim(scheme%turbulence))
                
                case ("none")
                    continue

                case ("sst", "sst2003", "bsl")
                    flow%tk_inf = 1.5*((flow%Vel_mag*flow%Turb_Intensity_inf)**2)
                    flow%tw_inf = flow%density_inf*flow%tk_inf/(flow%mu_ref*flow%mu_ratio_inf)

                case ("kkl")
                    flow%tk_inf = 9*(1e-9)*(sound_speed_inf(flow)**2)
                    flow%tkl_inf = 1.5589*(1e-6)*(flow%mu_ref*sound_speed_inf(flow))/flow%density_inf

                case ("sa")
                     flow%tv_inf = flow%mu_ratio_inf*flow%mu_ref/flow%density_inf

                case ("saBC")
                    flow%tv_inf = 0.005*flow%mu_ratio_inf*flow%mu_ref/flow%density_inf

                case ("kw")
                    flow%tk_inf = 1.5*((flow%Vel_mag*flow%Turb_Intensity_inf)**2)
                    flow%tw_inf = flow%density_inf*flow%tk_inf/(flow%mu_ref*flow%mu_ratio_inf)

                case ("ke")
                    flow%tk_inf = 1.5*((flow%Vel_mag*flow%Turb_Intensity_inf)**2)
                    flow%tw_inf = 0.09*flow%density_inf*flow%tk_inf*flow%tk_inf/(flow%mu_ref*flow%mu_ratio_inf)

                case ("des-sst")
                    flow%tk_inf = 1.5*((flow%Vel_mag*flow%Turb_Intensity_inf)**2)
                    flow%tw_inf = flow%density_inf*flow%tk_inf/(flow%mu_ref*flow%mu_ratio_inf)

                case ("les")
                  continue
                  ! todo

                case DEFAULT
                  Fatal_error

            end select


        end subroutine init_infinity_values