Load the color property from disk.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_color_prop_type), | intent(out) | :: | me |
Color definitions to load. |
||
type(treelmesh_type), | intent(in) | :: | tree |
Tree to build the polynomial subresolution information for |
||
integer, | intent(in) | :: | myPart |
Partition to load |
||
integer, | intent(in) | :: | comm |
Communicator to use |
subroutine tem_color_prop_load( me, tree, myPart, comm )
! --------------------------------------------------------------------------!
!> Color definitions to load.
type(tem_color_prop_type), intent(out) :: me
!> Tree to build the polynomial subresolution information for
type(treelmesh_type), intent(in) :: tree
!> Partition to load
integer, intent(in) :: myPart
!> Communicator to use
integer, intent(in) :: comm
! --------------------------------------------------------------------------!
integer, parameter :: root = 0
type(flu_State) :: conf
integer :: iError
integer :: rl
integer :: thandle
integer :: fUnit
integer :: i
integer :: iProp
character(len=pathLen) :: headerfile
character(len=pathLen) :: datafile
! --------------------------------------------------------------------------!
me%nColors = 0
prp_loop: do iprop=1, tree%global%nProperties
if (tree%global%Property(iprop)%bitpos == prp_isColored) then
me%header => tree%global%Property(iprop)
me%property => tree%property(iprop)
headerfile = trim(tree%global%dirname)//'colors.lua'
datafile = trim(tree%global%dirname)//'colors.ascii'
if (myPart == root) then
! Read the header only on the root process, broadcast to all others
! open mesh header file
call open_config_file( L = conf, filename = headerfile )
call aot_get_val( L = conf, &
& key = 'nColors', &
& val = me%nColors, &
& ErrCode = iError )
end if
call MPI_Bcast(me%nColors, 1, MPI_INTEGER, root, comm, iError)
! The number of colors that can be stored per character is fixed, thus
! the number of characters required by a given number of colors is
! immediatly known.
me%nChars = ceiling(real(me%nColors)/real(colors_per_char))
allocate(me%color_label(me%nColors))
allocate(me%color_fill(me%nColors))
allocate(me%color_void(me%nColors))
allocate(me%colored_bit(me%nChars, me%property%nElems))
if (myPart == root) then
! Now read the color labels on the root process.
call aot_table_open( L = conf, thandle = thandle, &
& key = 'color_label' )
do i=1,me%nColors
call aot_get_val( L = conf, &
& thandle = thandle, &
& pos = i, &
& val = me%color_label(i), &
& ErrCode = iError )
end do
call aot_table_close( L = conf, thandle = thandle )
! Now read the color fill values on the root process.
call aot_table_open( L = conf, thandle = thandle, &
& key = 'color_fill' )
do i=1,me%nColors
call aot_get_val( L = conf, &
& thandle = thandle, &
& pos = i, &
& val = me%color_fill(i), &
& ErrCode = iError )
end do
call aot_table_close( L = conf, thandle = thandle )
! Now read the color void values on the root process.
call aot_table_open( L = conf, thandle = thandle, &
& key = 'color_void' )
do i=1,me%nColors
call aot_get_val( L = conf, &
& thandle = thandle, &
& pos = i, &
& val = me%color_void(i), &
& ErrCode = iError )
end do
call aot_table_close( L = conf, thandle = thandle )
call close_config(conf)
end if
call MPI_Bcast( me%color_label, LabelLen*me%nColors, MPI_CHARACTER, &
& root, comm, iError )
call MPI_Bcast( me%color_fill, me%nColors, rk_mpi, &
& root, comm, iError )
call MPI_Bcast( me%color_void, me%nColors, rk_mpi, &
& root, comm, iError )
! If there are actually colored elements on the local process,
! read them now.
if (me%property%nElems > 0) then
inquire(iolength=rl) me%colored_bit(:,1)
call tem_open( newunit = fUnit, &
& file = datafile, &
& action = 'read', &
& access = 'stream', &
& form = 'unformatted', &
& status = 'old' )
read(fUnit, pos=me%property%offset+1) me%colored_bit
close(fUnit)
end if
EXIT prp_loop
end if
end do prp_loop
end subroutine tem_color_prop_load