Read one canonical object definition into a tem_canonicalND_type from a lua table.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_canonicalND_type), | intent(out) | :: | me |
contains canonicalND data |
||
type(tem_transformation_type), | intent(in) | :: | transform |
transformation for spatial object |
||
type(flu_State) | :: | conf |
lua state |
|||
integer, | intent(in) | :: | thandle |
lua table identification |
||
logical, | intent(in), | optional | :: | reqSegments |
Is true if use_get_point is true in output table |
subroutine tem_load_oneCanonicalND(me, transform, conf, thandle, reqSegments)
! --------------------------------------------------------------------------
!> contains canonicalND data
type(tem_canonicalND_type), intent(out) :: me
!> transformation for spatial object
type(tem_transformation_type), intent(in) :: transform
!> lua state
type(flu_state) :: conf
!> lua table identification
integer, intent(in) :: thandle
!> Is true if use_get_point is true in output table
logical, optional, intent(in) :: reqSegments
!> output messages?
! --------------------------------------------------------------------------
! error identifiers
integer :: iError
integer :: vError(3), errfatal(3)
! counter and total number of vectors
integer :: iVec, nVecs
! lua handles
integer :: vec_thandle, halfvec_thandle, single_thandle, halfsingle_thandle
real(kind=rk), allocatable :: local_vecs(:,:)
real(kind=rk) :: local_length
real(kind=rk) :: vecAbsSum
character(len=labelLen) :: buffer
logical :: reqSegments_loc
! --------------------------------------------------------------------------
! Load segments only if reqSegments is true i.e. use_get_point in output
! table is true
if (present(reqSegments)) then
reqSegments_loc = reqSegments
else
reqSegments_loc = .false.
end if
errfatal = aotErr_Fatal
! initialize the array of vectors and the array of segments
me%vec = 0._rk
! open the vec and halvec table and get its size
call aot_table_open( L = conf, &
& parent = thandle, &
& thandle = vec_thandle, &
& key = 'vec' )
call aot_table_open( L = conf, &
& parent = thandle, &
& thandle = halfvec_thandle, &
& key = 'halfvec' )
! check if the vec table is given
if( vec_thandle .ne. 0 ) then
nVecs = aot_table_length( L = conf, thandle = vec_thandle )
! Now check, if the first entry is a table itself:
call aot_table_open( L = conf, &
& parent = vec_thandle, &
& thandle = single_thandle, &
& pos = 1 )
if (single_thandle == 0) then
! The entries are not tables themselves, try to interpret it as
! single vector.
call aot_table_close( L = conf, thandle = single_thandle )
nVecs = 1
allocate(local_vecs( 3, nVecs ))
local_vecs = 0._rk
call aot_get_val( L = conf, &
& thandle = thandle, &
& val = local_vecs(:,1), &
& key = 'vec', &
& ErrCode = vError )
if( any( btest( vError, errFatal ))) then
write(logUnit(1),*)' FATAL ERROR: in loading vec'
write(logUnit(1),*)' You should provide 3 numbers for this vector.'
call tem_abort()
end if
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 )
allocate( local_vecs( 3, nVecs ))
local_vecs = 0._rk
do iVec=1, nVecs
call aot_get_val( L = conf, &
& thandle = vec_thandle, &
& val = local_vecs(:,iVec), &
& pos = iVec, &
& ErrCode = vError )
end do
end if
if( any( btest( vError, errFatal ))) then
write(logUnit(1),*) ' FATAL ERROR: in loading vec'
call tem_abort()
endif
! if vec table is not set try to read halfvec
else if( halfvec_thandle /= 0 ) then
call aot_table_open( L = conf, &
& parent = halfvec_thandle, &
& thandle = halfsingle_thandle, &
& pos = 1 )
nVecs = aot_table_length( L = conf, thandle = halfvec_thandle )
if( halfsingle_thandle == 0 ) then
! The entries are not tables themselves, try to interpret it as
! single vector.
call aot_table_close( L = conf, thandle = halfsingle_thandle )
nVecs = 1
allocate( local_vecs( 3, nVecs ))
local_vecs = 0._rk
call aot_get_val( L = conf, &
& thandle = thandle, &
& val = local_vecs(:,1), &
& key = 'halfvec', &
& ErrCode = vError )
local_vecs = local_vecs * 2._rk
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 = halfsingle_thandle )
allocate( local_vecs( 3, nVecs ))
local_vecs = 0._rk
do iVec=1, nVecs
call aot_get_val( L = conf, &
& thandle = halfvec_thandle, &
& val = local_vecs(:,iVec), &
& pos = iVec, &
& ErrCode = vError )
end do
local_vecs = local_vecs * 2._rk
end if
if( any( btest( vError, errFatal ))) then
write(logUnit(1),*) ' FATAL ERROR: in loading halfvec'
call tem_abort()
endif
! if vec and halfvec table is not set, try to read length (for a cube)
else
call aot_get_val( L = conf, &
& thandle = thandle, &
& val = local_length, &
& key = 'length', &
& default = 0.0_rk, &
& ErrCode = iError )
! set the axis parallel entries of local_vecs to local_length
nVecs = 3
allocate( local_vecs( 3, nVecs ))
local_vecs = 0._rk
local_vecs(1,1) = local_length
local_vecs(2,2) = local_length
local_vecs(3,3) = local_length
end if
call aot_table_close( L = conf, thandle = vec_thandle )
call aot_table_close( L = conf, thandle = halfvec_thandle )
! copy back the information from the local vector array
me%vec(:,:nVecs) = local_vecs
me%vec(:,nVecs+1:) = 0.0_rk
! read origin of the canonicalND
call aot_get_val( L = conf, &
& thandle = thandle, &
& val = me%origin, &
& key = 'origin', &
& ErrCode = vError )
if( any( btest( vError, errFatal ))) then
call aot_get_val( L = conf, &
& thandle = thandle, &
& val = me%origin, &
& key = 'center', &
& ErrCode = vError )
if (any(btest( vError, errFatal )) ) then
write(logUnit(1),*) 'ERROR reading the canonical objects, '// &
& 'origin/center is not well defined. STOPPING'
call tem_abort()
endif
me%origin = me%origin - me%vec(:,1)/2.0_rk - me%vec(:,2)/2.0_rk - &
& me%vec(:,3)/2.0_rk
end if
call aot_get_val( L = conf, &
& thandle = thandle, &
& val = me%distribution, &
& ErrCode = iError, &
& key = 'distribution', &
& default = 'equal' )
me%active = .false.
me%nDim = 0
! check the kind and the active vectors
do iVec = 1, 3
vecAbsSum = abs(me%vec(1,iVec)) + abs(me%vec(2,iVec)) &
& + abs(me%vec(3,iVec))
if (vecAbsSum .fne. 0._rk) then
me%active(iVec) = .true.
me%nDim = me%nDim + 1
end if
end do
if(me%nDim == 0) me%kind = 'point'
if(me%nDim == 1) me%kind = 'line'
if(me%nDim == 2) me%kind = 'plane'
if(me%nDim == 3) me%kind = 'box'
! Load segments for objects other than point if reqSegments is true
if (trim(me%kind) == 'point') then
! set default segments to 1
me%segments = 1
else
if (reqSegments_loc) then
! set segments to 1 in all dimensions as default
me%segments = 1
! try loading segments as integer if fails try
! to loading segments as vector integer
call aot_get_val( L = conf, &
& thandle = thandle, &
& val = me%segments(1), &
& key = 'segments', &
& ErrCode = iError )
if ( btest(iError, aoterr_Fatal) ) then
call aot_get_val( L = conf, &
& thandle = thandle, &
& val = me%segments(:nVecs), &
& key = 'segments', &
& ErrCode = vError(:nVecs) )
if ( any( btest(vError(:nVecs), errFatal(:nVecs)) ) ) then
write(logUnit(1),*) ' FATAL ERROR: in loading segments'
write(logUnit(1), '(A)') ' "segments" are not defined in ' &
& // 'shape.object table.'
call tem_abort('Segments are required for output.use_get_point.')
end if
end if
else
! segments are not required setting it to -1 for error detection
me%segments = -1
end if
end if
! if canonical object is box then check for only_surface
if( me%kind == 'box' ) then
call aot_get_val( L = conf, &
& thandle = thandle, &
& val = me%only_surface, &
& ErrCode = iError, &
& key = 'only_surface', &
& pos = 4, &
& default = .false. )
if (btest(iError, aoterr_WrongType)) then
write(logUnit(1),*) 'Error occured, while retrieving canonical box '// &
& 'only_surface'
write(logUnit(1),*) 'Variable has wrong type!'
write(logUnit(1),*) 'Should be a LOGICAL!'
call tem_abort()
endif
endif
! if tranformation is defined then tranform canoND before converting them
! into their respective shapes
call tem_transformCanoND(canoND = me, transform = transform)
! Convert canonical shape to respective kind
select case(upper_to_lower(trim(me%kind)))
case('point')
me%point%coord = me%origin
case('line')
me%line%origin = me%origin
me%line%vec = me%vec(:,1)
case('plane')
call tem_createPlane(me = me%plane, &
& origin = me%origin, &
& vecA = me%vec(:,1), &
& vecB = me%vec(:,2) )
case('box')
call tem_createBox(me = me%box, &
& origin = me%origin, &
& vecA = me%vec(:,1), &
& vecB = me%vec(:,2), &
& vecC = me%vec(:,3), &
& only_surface = me%only_surface)
case default
call tem_abort('Unknown canonical kind')
end select
! @todo: DEB: Update this when converting arrays to strings is available
write(logUnit(6),"(A)") ' kind: '//trim( me%kind )
write(logUnit(6),*) ' origin:', me%origin(1:3)
do iVec=1,me%nDim
write(buffer,"(A,I0,A)") ' vec ', iVec, ': '
write(logUnit(6),*) trim(buffer), me%vec(:,iVec)
enddo
if (trim(me%kind) == 'plane') then
write(logUnit(6),*) ' unit normal: ' &
& // trim(tem_toStr(me%plane%unitNormal(1))) &
& // ', ' // trim(tem_toStr(me%plane%unitNormal(2))) &
& // ', ' // trim(tem_toStr(me%plane%unitNormal(3)))
end if
write(logUnit(6),*) ' only_surface: ', me%only_surface
do iVec=1,me%nDim
write(buffer,"(A,I0,A)") 'segments ', iVec, ': '
write(logUnit(6),"(A,I0)") trim(buffer), me%segments(iVec)
enddo
write(logUnit(6),"(A)") ' distribution: '//trim( me%distribution )
end subroutine tem_load_oneCanonicalND