load_spatial_asConst Subroutine

private subroutine load_spatial_asConst(const, conf, errCode, parent, key, nComp)

Load spatial as constant

Arguments

Type IntentOptional Attributes Name
real(kind=rk), allocatable :: const(:)

value to be filled

type(flu_State) :: conf

lua state type

integer, intent(out) :: errCode

errCode = -1, if spatial function is not defined as constant

integer, intent(in) :: parent

aotus parent handle

character(len=*), intent(in) :: key

key is either "local_key" or "const"

integer, intent(in) :: nComp

number of components of the variable


Calls

proc~~load_spatial_asconst~~CallsGraph proc~load_spatial_asconst load_spatial_asConst proc~aot_get_val~2 aot_get_val proc~load_spatial_asconst->proc~aot_get_val~2 proc~tem_abort tem_abort proc~load_spatial_asconst->proc~tem_abort mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~load_spatial_asconst~~CalledByGraph proc~load_spatial_asconst load_spatial_asConst proc~tem_load_spatial tem_load_spatial proc~tem_load_spatial->proc~load_spatial_asconst proc~load_spatial_predefined load_spatial_predefined proc~tem_load_spatial->proc~load_spatial_predefined proc~load_spatial_parabol load_spatial_parabol proc~load_spatial_parabol->proc~load_spatial_asconst proc~load_spacetime_predefined load_spacetime_predefined proc~load_spacetime_predefined->proc~tem_load_spatial proc~load_spatial_predefined->proc~load_spatial_parabol 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 proc~tem_load_spacetime_single->proc~tem_load_spacetime_single proc~tem_load_spacetime_table tem_load_spacetime_table proc~tem_load_spacetime_table->proc~tem_load_spacetime_single interface~tem_load_spacetime tem_load_spacetime interface~tem_load_spacetime->proc~tem_load_spacetime_single

Contents

Source Code


Source Code

  subroutine load_spatial_asConst( const, conf, errCode, parent, key, nComp )
    ! -------------------------------------------------------------------- !
    !> value to be filled
    real(kind=rk), allocatable :: const(:)
    !> lua state type
    type(flu_State) :: conf
    !> errCode = -1, if spatial function is not defined as constant
    integer, intent(out) :: errCode
    !> aotus parent handle
    integer, intent(in) :: parent
    !> key is either "local_key" or "const"
    character(len=*), intent(in) :: key
    !> number of components of the variable
    integer, intent(in) :: nComp
    ! -------------------------------------------------------------------- !
    integer, allocatable :: vError(:), vErr_NonExistent(:), vErr_WrongType(:)
    ! -------------------------------------------------------------------- !

    errCode = 0

    if (nComp == 1) then
      write(logUnit(5),*) 'Load spatial as single constant'
      allocate(const(nComp), vError(nComp))
      call aot_get_val( L       = conf,     &
        &               thandle = parent,   &
        &               val     = const(1), &
        &               key     = key,      &
        &               ErrCode = vError(1) )

      if (btest(vError(1), aoterr_NonExistent)) then
        ! if previous checks fails return errCode =-1
        ! and set default kind = none and value = 1
        write(logUnit(3),*) 'ERROR: Fatal error occured in loading ' &
          & // 'scalar constants'
        errCode = -1
        return
      end if
    else !nComp > 1
      write(logUnit(3),"(A,I0)") &
        & 'Load spatial as array of constant with nComp ', nComp
      !Try to load constant of nComp
      call aot_get_val( L         = conf,   &
        &               thandle   = parent, &
        &               val       = const,  &
        &               maxLength = nComp,  &
        &               key       = key,    &
        &               ErrCode   = vError  )
      if (size(const) == 0) then
        write(logUnit(3),*) 'Error loading array constant from key "'&
          & // trim(key) // '"'
        errCode = -1
        return
      else
        ! if any error index is fatal then return errCode - 1
        allocate(vErr_NonExistent(size(vError)))
        vErr_NonExistent = aoterr_NonExistent
        if(any(btest(vError, vErr_NonExistent))) then
          write(logUnit(3),*) 'ERROR: Fatal error occured in loading ' &
            & // 'array of constants with nComp:', nComp
          errCode = -1
          return
        else if (size(const) /= nComp) then
          ! no fatal error but number of constant defined /= nComp
          write(logUnit(1),*) 'Spatial as constant. table size /= nComp', &
            &                 nComp
          call tem_abort()
        end if
      end if
    end if

    allocate(vErr_WrongType(size(vError)))
    vErr_WrongType = aoterr_WrongType
    if (any(btest(vError, vErr_WrongType))) then
      write(logUnit(1),*)'FATAL Error occured in definition of ' &
        & // 'the spatial'
      write(logUnit(1),*)'while retrieving spatial constant:'
      write(logUnit(1),*)'Variable has wrong type (should be a ' &
        & // 'real number)!'
      call tem_abort()
    end if

  end subroutine load_spatial_asConst