init_state_with_infinity_values Subroutine

private subroutine init_state_with_infinity_values(qp, scheme, flow, dims)

Initialize the state based on the infinity values

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(inout), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 1:dims%n_var), target:: qp

Store primitive variable at cell center

type(schemetype), intent(in) :: scheme

finite-volume Schemes

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

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

type(extent), intent(in) :: dims

Extent of the domain:imx,jmx,kmx


Calls

proc~~init_state_with_infinity_values~~CallsGraph proc~init_state_with_infinity_values init_state_with_infinity_values debugcall debugcall proc~init_state_with_infinity_values->debugcall

Called by

proc~~init_state_with_infinity_values~~CalledByGraph proc~init_state_with_infinity_values init_state_with_infinity_values proc~initstate initstate proc~initstate->proc~init_state_with_infinity_values proc~setup_state setup_state proc~setup_state->proc~initstate 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

        subroutine init_state_with_infinity_values(qp, scheme, flow, dims)
            !< Initialize the state based on the infinity values
            !-----------------------------------------------------------
            
            implicit none
            type(schemetype), intent(in) :: scheme
            !< finite-volume Schemes
            type(flowtype), intent(in) :: flow
            !< Information about fluid flow: freestream-speed, ref-viscosity,etc.
            type(extent), intent(in) :: dims
            !< Extent of the domain:imx,jmx,kmx
            real(wp), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 1:dims%n_var), intent(inout), target :: qp
            !< Store primitive variable at cell center
            
            DebugCall("init_state_with_infinity_values")

            !density = density_inf
            qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 1) = flow%density_inf
            !x_speed = x_speed_inf
            qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 2) = flow%x_speed_inf
            !y_speed = y_speed_inf
            qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 3) = flow%y_speed_inf
            !z_speed = z_speed_inf
            qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 4) = flow%z_speed_inf
            !pressure = pressure_inf
            qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 5)= flow%pressure_inf

            select case (trim(scheme%turbulence))

                case ("none")
                    !include nothing
                    continue
                
                case ("sst", "sst2003", "bsl", "des-sst", "kw")
                    !tk = tk_inf
                    qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 6) = flow%tk_inf
                    !tw = tw_inf
                    qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 7) = flow%tw_inf

                case ("kkl")
                    !tk = tk_inf
                    qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 6) = flow%tk_inf
                    !tkl = tkl_inf
                    qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 7) = flow%tkl_inf

                case ("sa", "saBC")
                    !tv = tv_inf
                    qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 6) = flow%tv_inf

                case ("ke")
                    !tk = tk_inf
                    qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 6) = flow%tk_inf
                    !te = te_inf
                    qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 7) = flow%te_inf

                case ("les")
                  continue
                  ! todo

                case DEFAULT
                  Fatal_error

            end select


            ! Transition modeling
            select case(trim(scheme%transition))

              case('lctm2015')
                !tgm = tgm_inf
                qp(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 8) = flow%tgm_inf

              case('bc', 'none')
                !do nothing
                continue

              case DEFAULT
                Fatal_error

            end Select
            
        end subroutine init_state_with_infinity_values