tem_load_cube Subroutine

public subroutine tem_load_cube(me, conf, key, pos, parent)

This routine loads the boundCube table from config file

Arguments

Type IntentOptional Attributes Name
type(tem_cube_type), intent(out) :: me
type(flu_State) :: conf
character(len=*), intent(in), optional :: key

open cube table by given key

integer, intent(in), optional :: pos

open cube table by position

integer, intent(in), optional :: parent

if cube is to be load from pos, parent handle is required


Calls

proc~~tem_load_cube~~CallsGraph proc~tem_load_cube tem_load_cube proc~aot_table_open aot_table_open proc~tem_load_cube->proc~aot_table_open proc~aot_table_close aot_table_close proc~tem_load_cube->proc~aot_table_close proc~tem_abort tem_abort proc~tem_load_cube->proc~tem_abort proc~aot_get_val~2 aot_get_val proc~tem_load_cube->proc~aot_get_val~2 mpi_abort mpi_abort proc~tem_abort->mpi_abort

Contents

Source Code


Source Code

  subroutine tem_load_cube(me, conf, key, pos, parent)
    ! --------------------------------------------------------------------------!
    type(tem_cube_type), intent(out) :: me !< seeder cube type
    type(flu_state) :: conf !< lua state
    !> open cube table by given key
    character(len=*), optional, intent(in) :: key
    !> open cube table by position
    integer, optional, intent(in) :: pos
    !> if cube is to be load from pos, parent handle is required
    integer, optional, intent(in) :: parent
    ! --------------------------------------------------------------------------!
    integer :: thandle !< cube handle
    integer :: iError, vError(3), errFatal(3)
    real(kind=rk) :: deflen
    ! --------------------------------------------------------------------------!
    errFatal = aoterr_fatal

    me%origin = 0.0
    deflen = 1.0_rk

    if (present(key)) then
      ! open cube table from given key
      call aot_table_open(L=conf, thandle=thandle, key=trim(key))
      if (thandle == 0) then
        write(logunit(0),*) 'FATAL Error: cube definition not found with key:'
        write(logunit(0),*) (trim(key))
        call tem_abort()
      end if
    else if (present(parent) .and. present(pos)) then
      ! else open cube from given table pos
      call aot_table_open(L=conf, parent = parent, thandle=thandle, pos=pos)
      if (thandle == 0) then
        write(logunit(0),*) &
          &  'FATAL Error: cube definition not found with pos:', pos
        call tem_abort()
      end if
    else
      write(logunit(0),*) 'ERROR: Neither key nor pos is provided to load ' &
        &                 //' cube table'
      call tem_abort()
    end if

    call aot_get_val(L=conf, thandle=thandle, val=me%origin, ErrCode=vError, &
      &              key='origin', pos = 1, default=[0.0_rk, 0.0_rk, 0.0_rk] )
    if (any(btest(vError, errFatal))) then
      write(logunit(0),*) 'FATAL Error occured, while retrieving cube origin.'
      call tem_abort()
    end if

    call aot_get_val(L=conf, thandle=thandle, &
      &              val=me%extent, ErrCode=iError, &
      &              key='length', pos=2, default=deflen)
    if (btest(iError, aoterr_Fatal)) then
      write(logunit(0),*) 'FATAL Error occured, while retrieving cube length:'
      if (btest(iError, aoterr_NonExistent)) &
        &  write(logunit(0),*) 'Variable not existent!'
      if (btest(iError, aoterr_WrongType)) &
        &  write(logunit(0),*) 'Variable has wrong type!'
      call tem_abort()
    else
      if (btest(iError, aoterr_NonExistent)) &
        & write(logunit(1),*) 'Variable length not set in configuration, ' &
        &                     //'using default value!'
    end if

    me%halfwidth = 0.5_rk*me%extent
    me%center = me%origin + me%halfwidth

    write(logunit(1),*) ' Cube origin: ', me%origin
    write(logunit(1),*) ' Cube length: ', me%extent

    call aot_table_close(L=conf, thandle=thandle)

  end subroutine tem_load_cube