Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=long_k), | intent(out) | :: | val(:) | Vector read from the Lua table. |
||
integer, | intent(out) | :: | ErrCode(:) | Error code describing problems encountered in each of the components. This array has to have the same length as val. |
||
type(flu_State) | :: | L | Handle to the lua script |
|||
integer(kind=long_k), | intent(in), | optional | :: | default(:) | A default vector to use, if no proper definition is found. Components will be filled with the help of this default definition. |
subroutine get_top_long_v(val, ErrCode, L, default)
type(flu_State) :: L !! Handle to the lua script
!> Vector read from the Lua table.
integer(kind=long_k), intent(out) :: val(:)
!> Error code describing problems encountered in each of the components.
!! This array has to have the same length as val.
integer, intent(out) :: ErrCode(:)
!> A default vector to use, if no proper definition is found.
!! Components will be filled with the help of this default definition.
integer(kind=long_k), intent(in), optional :: default(:)
integer :: vect_handle
integer :: table_len, vect_len, def_len, val_len
integer :: vect_lb
integer :: iComp
ErrCode = 0
! Try to interpret it as table.
vect_handle = aot_table_top(L=L)
table_len = aot_table_length(L=L, thandle=vect_handle)
val_len = size(val)
vect_len = min(table_len, val_len)
! Find the length of the default value, if it is not provided, its 0.
def_len = 0
if (present(default)) def_len = size(default)
! Now parse the table with the vector entries.
if (aot_table_first(L, vect_handle).and.(vect_len > 0)) then
! Up to the length of the default value, provide the default settings.
do iComp=1,def_len
call aot_top_get_val(val(iComp), ErrCode(iComp), L, &
& default(iComp))
if (.not. flu_next(L, vect_handle)) exit
end do
vect_lb = def_len+1
! After def_len entries no default values for the components are
! available anymore, proceed without a default setting for the rest.
do iComp=vect_lb,vect_len
call aot_top_get_val(val(iComp), ErrCode(iComp), L)
if (.not. flu_next(L, vect_handle)) exit
end do
! If the table in the Lua script is not long enough, fill the remaining
! components with the default components, as far as they are defined.
do iComp=vect_len+1,def_len
ErrCode(iComp) = ibSet(ErrCode(iComp), aoterr_NonExistent)
val(iComp) = default(iComp)
end do
vect_lb = max(vect_len+1, def_len+1)
do iComp=vect_lb,val_len
ErrCode(iComp) = ibSet(ErrCode(iComp), aoterr_Fatal)
end do
else
! No vector definition found in the Lua script, use the default.
ErrCode = ibSet(ErrCode, aoterr_NonExistent)
if (present(default)) then
val(:def_len) = default(:def_len)
if (def_len < val_len) then
ErrCode(def_len+1:) = ibSet(ErrCode(def_len+1:), aoterr_Fatal)
end if
else
! No vector definition in the Lua script and no default provided.
ErrCode = ibSet(ErrCode, aoterr_Fatal)
end if
end if
call aot_table_close(L, vect_handle)
end subroutine get_top_long_v