tem_load_spongeLayer_plane Subroutine

public subroutine tem_load_spongeLayer_plane(me, conf, thandle, nDim, nComp, stateName)

This subroutine load data for standard plane sponge layer Example:

 spatial = {
   -- supported options: 'spongelayer_plane', 'spongelayer_plane_1d', 
   --                    'spongelayer_plane_2d', 'viscous_spongelayer_plane'
   predefined = 'spongelayer',
   origin = {0.0,0.0,0.0},
   normal = {1.0, 0.0, 0.0},
   thickness = 0.5,
   damp_profile = 'linear', --'exponential', 'polynomial_n5', 'polynomial_n6'
   damp_factor = 0.5,
   damp_exponent = 1.0,
   target_state = {
     Default: density, velocityX, velocityY, velocityZ and pressure
     density = 1.0,
     pressure = 1.0,
     velocityX = 0.0, velocityY = 0.0, velocityZ = 0.0
 }

Arguments

Type IntentOptional Attributes Name
type(tem_spongeLayer_plane_type), intent(out) :: me

Plane spongeLayer data type

type(flu_State) :: conf

lua state type

integer, intent(in) :: thandle

aotus parent handle

integer, intent(in) :: nDim

number of Dimension for nonViscous sponges

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


Calls

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

Called by

proc~~tem_load_spongelayer_plane~~CalledByGraph proc~tem_load_spongelayer_plane tem_load_spongeLayer_plane proc~load_spatial_predefined load_spatial_predefined proc~load_spatial_predefined->proc~tem_load_spongelayer_plane 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 proc~tem_load_spacetime_single tem_load_spacetime_single proc~tem_load_spacetime_single->proc~load_spacetime_predefined

Contents


Source Code

  subroutine tem_load_spongeLayer_plane(me, conf, thandle, ndim, nComp, &
    &                                   stateName)
    ! --------------------------------------------------------------------------
    !> Plane spongeLayer data type
    type(tem_spongeLayer_plane_type), intent(out) :: me
    !> lua state type
    type(flu_State) :: conf
    !> aotus parent handle
    integer, intent(in) :: thandle
    !> number of Dimension for nonViscous sponges
    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
    ! --------------------------------------------------------------------------
    integer :: vError(3), errfatal(3)
    real(kind=rk) :: thickness
    ! --------------------------------------------------------------------------
    errfatal = aotErr_Fatal
    ! Plane_origin
    call aot_get_val( L       = conf,      &
      &               thandle = thandle,   &
      &               key     = 'origin',  &
      &               val     = me%origin, &
      &               ErrCode = vError     )
    if (any(btest( vError, errFatal )) ) then
      write(logUnit(1),*) 'ERROR reading the plane_origin of sponge layer. ' &
        &              // 'It should have 3 entries for each coordinate.'
      call tem_abort()
    end if

    write(logUnit(1),*) ' * Origin =', me%origin

    ! Plane_normal
    call aot_get_val( L       = conf,      &
      &               thandle = thandle,   &
      &               key     = 'normal',  &
      &               val     = me%normal, &
      &               ErrCode = vError     )
    if (any(btest( vError, errFatal )) ) then
      write(logUnit(1),*) 'ERROR reading the plane_normal of sponge layer. ' &
        &              // 'It should have 3 entries for each coordinate.'
      call tem_abort()
    end if
    write(logUnit(1),*) ' * Normal =', me%normal

    ! Compute thickness from normal and use it only if thickness is not
    ! defined seperately.
    thickness = sqrt(me%normal(1)**2 + me%normal(2)**2 + me%normal(3)**2)

    ! Normalize the normal
    me%normal = me%normal/thickness
    write(logUnit(1),*) ' * Normalized normal =', me%normal

    ! Load base information required for sponge layer definition like
    ! damp_factor, damp_exponent and target_state
    call load_spongeLayer( conf      = conf,                     &
      &                    thandle   = thandle,                  &
      &                    me        = me%spongeLayer_base_type, &
      &                    nDim      = nDim,                     &
      &                    nComp     = nComp,                    &
      &                    stateName = stateName,                &
      &                    thickness = thickness                 )

  end subroutine tem_load_spongeLayer_plane