Loading canonicalNDs from the config file valid definitions:
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_canonicalND_type), | intent(out), | allocatable | :: | me(:) |
the array of canonical objects which to read in (and first allocate) |
|
type(tem_transformation_type), | intent(in) | :: | transform |
transformation for spatial object |
||
type(flu_State) | :: | conf |
lua config handle |
|||
integer, | intent(in) | :: | thandle |
table handle from which to read |
||
logical, | intent(in), | optional | :: | reqSegments |
Is true if use_get_point is true in output table |
subroutine tem_load_canonicalND_vec(me, transform, conf, thandle, reqSegments) ! -------------------------------------------------------------------------- !> the array of canonical objects which to read in (and first allocate) type(tem_canonicalND_type), allocatable, intent(out) :: me(:) !> transformation for spatial object type(tem_transformation_type), intent(in) :: transform !> lua config handle type(flu_state) :: conf !> table handle from which to read integer, intent(in) :: thandle !> Is true if use_get_point is true in output table logical, optional, intent(in) :: reqSegments ! -------------------------------------------------------------------------- ! lua handles integer :: canonicalNDs_thandle, single_thandle integer :: canonicalND_entries integer :: iCan integer :: ncanonicalNDs ! -------------------------------------------------------------------------- write(logUnit(5),*) 'Trying to read canoND shape...' call aot_table_open( L = conf, & & parent = thandle, & & thandle = canonicalNDs_thandle, & & key = 'object' ) canonicalND_entries = aot_table_length( L = conf, & & thandle = canonicalNDs_thandle ) if (canonicalNDs_thandle == 0) then write(logUnit(1),*) 'Error: object table is not defined for canoND' call tem_abort() end if ! Now check, if the first entry is a table itself: call aot_table_open( L = conf, & & parent = canonicalNDs_thandle, & & thandle = single_thandle, & & pos = 1 ) if (single_thandle == 0) then ! The entries are not tables themselves, try to interpret it as ! single canonical object. call aot_table_close( L = conf, thandle = single_thandle ) ncanonicalNDs = 1 allocate(me( ncanonicalNDs )) call tem_load_oneCanonicalND( me = me(1), & & conf = conf, & & transform = transform, & & thandle = canonicalNDs_thandle, & & reqSegments = reqSegments ) else ! First entry is a table, assume all others are also tables, ! and properly define canonicalNDs coordinates, go on reading them. call aot_table_close( L = conf, thandle = single_thandle ) ncanonicalNDs = canonicalND_entries allocate( me( ncanonicalNDs )) do iCan=1, ncanonicalNDs write(logUnit(5),"(A,I0,A)") 'Read Object: ', iCan, ' with ' call aot_table_open( L = conf, & & parent = canonicalNDs_thandle, & & thandle = single_thandle, & & pos = iCan ) call tem_load_oneCanonicalND( me = me( iCan ), & & conf = conf, & & transform = transform, & & thandle = single_thandle, & & reqSegments = reqSegments ) call aot_table_close( L = conf, thandle = single_thandle ) end do end if call aot_table_close( L = conf, thandle = canonicalNDs_thandle ) end subroutine tem_load_canonicalND_vec