load_spongeLayer Subroutine

private subroutine load_spongeLayer(conf, thandle, me, ndim, nComp, stateName, thickness)

This routine load base info for sponge layer

Arguments

Type IntentOptional Attributes Name
type(flu_State) :: conf

lua state type

integer, intent(in) :: thandle

aotus parent handle

type(spongeLayer_base_type), intent(out) :: me

base spongeLayer data type

integer, intent(in) :: ndim

number of spatial dimensions

integer, intent(in) :: nComp

Number of component of St-Fun variable under which this spatial function is defined

character(len=*), intent(in), optional :: stateName

Load stateName from target_state table

real(kind=rk), intent(in), optional :: thickness

Thickness computed from sponge layer plane normal. Use this thickness If thickness is not defined.


Calls

proc~~load_spongelayer~~CallsGraph proc~load_spongelayer load_spongeLayer proc~load_defaulttargetstate load_defaultTargetState proc~load_spongelayer->proc~load_defaulttargetstate proc~aot_table_open aot_table_open proc~load_spongelayer->proc~aot_table_open proc~aot_table_close aot_table_close proc~load_spongelayer->proc~aot_table_close proc~aot_get_val aot_get_val proc~load_spongelayer->proc~aot_get_val proc~tem_abort tem_abort proc~load_spongelayer->proc~tem_abort proc~load_defaulttargetstate->proc~aot_table_open proc~load_defaulttargetstate->proc~aot_table_close proc~load_defaulttargetstate->proc~aot_get_val proc~load_defaulttargetstate->proc~tem_abort mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~load_spongelayer~~CalledByGraph proc~load_spongelayer load_spongeLayer proc~tem_load_spongelayer_plane tem_load_spongeLayer_plane proc~tem_load_spongelayer_plane->proc~load_spongelayer proc~tem_load_spongelayer_box tem_load_spongeLayer_box proc~tem_load_spongelayer_box->proc~load_spongelayer proc~tem_load_spongelayer_radial tem_load_spongeLayer_radial proc~tem_load_spongelayer_radial->proc~load_spongelayer proc~load_spatial_predefined load_spatial_predefined proc~load_spatial_predefined->proc~tem_load_spongelayer_plane proc~load_spatial_predefined->proc~tem_load_spongelayer_box proc~load_spatial_predefined->proc~tem_load_spongelayer_radial proc~tem_load_spatial tem_load_spatial proc~tem_load_spatial->proc~load_spatial_predefined proc~load_spacetime_predefined load_spacetime_predefined proc~load_spacetime_predefined->proc~tem_load_spatial proc~tem_load_ic tem_load_ic proc~tem_load_ic->proc~tem_load_spatial

Contents

Source Code


Source Code

  subroutine load_spongeLayer(conf, thandle, me, ndim, nComp, stateName, &
    &                         thickness)
    ! --------------------------------------------------------------------------
    !> lua state type
    type(flu_State) :: conf
    !> aotus parent handle
    integer, intent(in) :: thandle
    !> base spongeLayer data type
    type(spongeLayer_base_type), intent(out) :: me
    !> number of spatial dimensions
    integer, intent(in) :: ndim
    !> Number of component of St-Fun variable under which this spatial function
    !! is defined
    integer, intent(in) :: nComp
    !> Load stateName from target_state table
    character(len=*), intent(in), optional :: stateName
    !> Thickness computed from sponge layer plane normal. Use this thickness
    !! If thickness is not defined.
    real(kind=rk), intent(in), optional :: thickness
    ! --------------------------------------------------------------------------
    integer :: iError, ts_handle
    ! --------------------------------------------------------------------------
    ! Thickness
    call aot_get_val( L       = conf,         &
      &               thandle = thandle,      &
      &               key     = 'thickness',  &
      &               val     = me%thickness, &
      &               ErrCode = iError        )
    if (btest( iError, aotErr_Fatal )) then
      if (present(thickness)) then
        me%thickness = thickness
      else
        write(logUnit(1),*) 'ERROR reading the thickness of sponge layer.'
        write(logUnit(1),*) 'Thickness is required to calculate sponge end.'
        call tem_abort()
      end if
    end if

    write(logUnit(1),*) ' * Thickness =', me%thickness

    !damp_factor
    call aot_get_val( L       = conf,          &
      &               thandle = thandle,       &
      &               key     = 'damp_factor', &
      &               val     = me%dampFactor, &
      &               ErrCode = iError         )
    if (btest( iError, aotErr_Fatal )) then
      write(logUnit(1),*) 'ERROR reading the damp_factor of sponge layer.'
      call tem_abort()
    end if
    write(logUnit(1),*) ' * Damp_factor =', me%dampFactor

    !damp_profile
    call aot_get_val( L       = conf,           &
      &               thandle = thandle,        &
      &               key     = 'damp_profile', &
      &               val     = me%dampProfile, &
      &               ErrCode = iError          )

    ! Viscous sponge works only with exponential profile so no need load 
    ! damp_profile
    if (present(stateName)) then
      if (trim(stateName) == 'viscosity') then
        me%dampProfile = 'exponential'
        iError = 0
      end if
    end if

    if (btest( iError, aotErr_Fatal )) then
      write(logUnit(1),*) 'ERROR reading the damp_profile of sponge layer.'
      call tem_abort()
    end if
    write(logUnit(1),*) ' * Damp_profile =', me%dampProfile

    select case (trim(me%dampProfile))
    case ('linear')
      me%dampExponent = 1.0
    case ('exponential')
      !damp_exponent
      call aot_get_val( L       = conf,            &
        &               thandle = thandle,         &
        &               key     = 'damp_exponent', &
        &               val     = me%dampExponent, &
        &               ErrCode = iError,          &
        &               default = 1.0_rk           )

      if (btest( iError, aotErr_Fatal )) then
        write(logUnit(1),*) 'ERROR reading the damp_exponent of exponential ' &
          &                 //'sponge layer.'
        call tem_abort()
      end if
    case ('polynomial_n5', 'polynomial_n6')
      me%dampExponent = 1.0
    case default
      write(logUnit(1),*) 'ERROR unknown damp_profile for sponge layer.'
      write(logUnit(1),*) 'Supported options: '
      write(logUnit(1),*) '  linear, exponential, polynomial_n5, polynomial_n6'
      call tem_abort()
    end select

    write(logUnit(1),*) ' * Damp_exponent =', me%dampExponent


    ! Load stateName provided by caller function
    if ( present(stateName) .and. nComp == 1) then
      allocate(me%targetState(1))
      call aot_table_open( L       = conf,          &
        &                  parent  = thandle,       &
        &                  thandle = ts_handle,     &
        &                  key     = 'target_state' )

      call aot_get_val( L       = conf,               &
        &               thandle = ts_handle,          &
        &               key     = trim(stateName),    &
        &               val     = me%targetState(1),  &
        &               ErrCode = iError              )

      if (btest(iError, aotErr_Fatal)) then
        write(logUnit(1),*) 'ERROR reading the target state: '//trim(stateName)
        call tem_abort()
      end if
      call aot_table_close(conf, thandle)

      write(logUnit(1),*) ' * Target state:'
      write(logUnit(1),*) '   '//trim(stateName)//' =', me%targetState(1)
    else if (nComp > 1) then
      call load_defaultTargetState( conf        = conf,          &
        &                           parent      = thandle,       &
        &                           nDim        = nDim,          &
        &                           nComp       = nComp,         &
        &                           targetState = me%targetState )
    else
      write(logUnit(1),*) 'WARNING: nComp = 1 so no target states are loaded'
    end if

  end subroutine load_spongeLayer