verify_read_control Subroutine

public subroutine verify_read_control(control, scheme)

Verify all the variable being asked to read in the output file. This is a fail-safe subroutine which do not allow to read the incorrect input variable. Based on previous flow type some varible might be skipped

Arguments

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

Calls

proc~~verify_read_control~~CallsGraph proc~verify_read_control verify_read_control proc~lcase lcase proc~verify_read_control->proc~lcase

Called by

proc~~verify_read_control~~CalledByGraph proc~verify_read_control verify_read_control proc~read_file~3 read_file proc~read_file~3->proc~verify_read_control proc~initstate initstate proc~initstate->proc~read_file~3 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


Source Code

    subroutine verify_read_control(control, scheme)
      !< Verify all the variable being asked to read in the output file. 
      !< This is a fail-safe subroutine which do not allow to read the incorrect input variable. 
      !< Based on previous flow type some varible might be skipped
      implicit none
      type(controltype), intent(inout) :: control
      type(schemetype) , intent(in) :: scheme
      integer :: n
      character(len=*), parameter :: err="Control Error: can't read variable - "

      do n = 1,control%r_count

        select case (trim(lcase(control%r_list(n))))
        
          case('velocity','vel','speed','u','v')
            control%r_list(n) = "Velocity"

          case('density','rho')
            control%r_list(n) = "Density"
          
          case('pressure','presssure','p')
            control%r_list(n) = "Pressure"

          case('mu','viscosity','mu_l','laminar_viscosity','muv','mu_v')
            control%r_list(n) = "do not read"
           ! if (flow%mu_ref/=0.0) then
           !   control%r_list(n) = "Mu"
           ! else
           !   print*, err//trim(control%r_list(n))//" from file"
           !   control%r_list(n) = "do not read"
           ! end if
            
          case('mu_t','turbulent_viscosity','mut')
            control%r_list(n) = "do not read"
            !if (scheme%turbulence/='none') then
            !  control%r_list(n) = "Mu_t"
            !else
            !  print*, err//trim(control%r_list(n))//" from file"
            !  control%r_list(n) = "do not read"
            !end if
            
          case('tke','tk','turbulent_kinetic_enrgy','k')
            select case (trim(scheme%turbulence))
              case('sst','sst2003','kw','bsl','kkl','ke','des-sst')
                select case (trim(control%previous_flow_type))
                  case('sst','sst2003','kw','bsl','kkl','ke','des-sst')
                    control%r_list(n) = "TKE"
                end select
              case DEFAULT
                print*, err//trim(control%w_list(n))//" from file"
                control%r_list(n) = "do not read"
            end select

          case('omega','tw')
            select case (trim(scheme%turbulence))
              case('sst','sst2003','kw','bsl','des-sst')
                select case (trim(control%previous_flow_type))
                  case('sst','sst2003','kw','bsl','des-sst')
                    control%r_list(n) = "Omega"
                  case DEFAULT
                    print*, err//trim(control%w_list(n))//" from file"
                    control%r_list(n) = "do not read"
                end select
              case DEFAULT
                print*, err//trim(control%w_list(n))//" from file"
                control%r_list(n) = "do not read"
            end select

          case('dissipation','te','teps','eps')
            select case (trim(scheme%turbulence))
              case('ke')
                select case (trim(control%previous_flow_type))
                  case('ke')
                    control%r_list(n) = "Dissipation"
                  case DEFAULT
                    print*, err//trim(control%w_list(n))//" to file"
                    control%r_list(n) = "do not write"
                end select
              case DEFAULT
                print*, err//trim(control%w_list(n))//" from file"
                control%r_list(n) = "do not read"
            end select

          case('kl')
            select case (trim(scheme%turbulence))
              case('kkl')
                select case (trim(control%previous_flow_type))
                  case('kkl')
                    control%r_list(n) = "Kl"
                  case DEFAULT
                    print*, err//trim(control%w_list(n))//" to file"
                    control%r_list(n) = "do not write"
                end select
              case DEFAULT
                print*, err//trim(control%w_list(n))//" from file"
                control%r_list(n) = "do not read"
            end select

          case('tv')
            select case (trim(scheme%turbulence))
              case('sa', 'saBC')
                select case (trim(control%previous_flow_type))
                  case('sa', 'saBC')
                    control%r_list(n) = "tv"
                  case DEFAULT
                    print*, err//trim(control%w_list(n))//" to file"
                    control%r_list(n) = "do not write"
                end select
              case DEFAULT
                print*, err//trim(control%w_list(n))//" from file"
                control%r_list(n) = "do not read"
            end select

          case('wall_distance', 'dist', 'wall_dist', 'wdist')
            control%r_list(n) = "do not read"
            !if(scheme%turbulence/="none") then
            !  control%r_list(n) = "Wall_distance"
            !else
            !  print*, err//trim(control%r_list(n))//" from file"
            !  control%r_list(n) = "do not read"
            !end if

          case('intermittency')
            select case (trim(scheme%turbulence))
              case('saBC')
                select case(trim(control%previous_flow_type))
                  case('saBC')
                    control%r_list(n) = "Intermittency"
                  case DEFAULT
                    print*, err//trim(control%r_list(n))//" to file"
                    control%r_list(n) = "do not read"
                end select
              case DEFAULT
                print*, err//trim(control%r_list(n))//" to file"
                control%r_list(n) = "do not read"
             end select

          case('tgm')
            select case (trim(scheme%transition))
              case('lctm2015')
                select case(trim(control%previous_flow_type))
                  case('sst','sst2003')
                    control%r_list(n) = "tgm"
                  case DEFAULT
                    print*, err//trim(control%r_list(n))//" to file"
                    control%r_list(n) = "do not read"
                end select
              case DEFAULT
                print*, err//trim(control%r_list(n))//" to file"
                control%r_list(n) = "do not read"
             end select

          case('resnorm')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Resnorm"

          case('tke_residue')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "TKE_residue"

          case('omega_residue')
            control%r_list(n) = "do not read"
            !control%w_list(n) = "Omega_residue"

          case('f1')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "F1"

          case('dudx')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dudx"

          case('dudy')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dudy"

          case('dudz')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dudz"

          case('dvdx')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dvdx"

          case('dvdy')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dvdy"

          case('dvdz')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dvdz"

          case('dwdx')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dwdx"

          case('dwdy')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dwdy"

          case('dwdz')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dwdz"

          case('dTdx')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "DTdx"

          case('dTdy')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "DTdy"

          case('dTdz')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "DTdz"

          case('dtkdx')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dtkdx"

          case('dtkdy')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dtkdy"

          case('dtkdz')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dtkdz"

          case('dtwdx')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dtwdx"

          case('dtwdy')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dtwdy"

          case('dtwdz')
            control%r_list(n) = "do not read"
            !control%r_list(n) = "Dtwdz"

          case('extravar1','extravar2', 'extravar3', 'extravar4', 'extravar5')
            control%r_list(n) = "do not read"
            !control%r_list(n) = trim(lcase(control%w_list(n)))

          case Default
            print*, err//trim(control%r_list(n))//" from file"
            control%r_list(n) = "do not read"

        end select
      end do

    end subroutine verify_read_control