tem_subres_prop_load Subroutine

public subroutine tem_subres_prop_load(me, tree, coloring)

Load the subresolution property from disk.

Before this can be done, the coloring information has to have been loaded.

Arguments

Type IntentOptional Attributes Name
type(tem_subres_prop_type), intent(out) :: me

Color definitions to load.

type(treelmesh_type), intent(in) :: tree

Tree to build the polynomial subresolution information for

type(tem_color_prop_type), intent(in) :: coloring

Information on the colors in the mesh.


Calls

proc~~tem_subres_prop_load~~CallsGraph proc~tem_subres_prop_load tem_subres_prop_load proc~tem_open tem_open proc~tem_subres_prop_load->proc~tem_open mpi_exscan mpi_exscan proc~tem_subres_prop_load->mpi_exscan proc~upper_to_lower upper_to_lower proc~tem_open->proc~upper_to_lower proc~tem_abort tem_abort proc~tem_open->proc~tem_abort proc~newunit newunit proc~tem_open->proc~newunit mpi_abort mpi_abort proc~tem_abort->mpi_abort

Contents

Source Code


Source Code

  subroutine tem_subres_prop_load( me, tree, coloring )
    ! -----------------------------------------------------------------------!
    !> Color definitions to load.
    type(tem_subres_prop_type), intent(out) :: me

    !> Tree to build the polynomial subresolution information for
    type(treelmesh_type), intent(in) :: tree

    !> Information on the colors in the mesh.
    type(tem_color_prop_type), intent(in) :: coloring
    ! -----------------------------------------------------------------------!
    integer :: rl
    integer :: fUnit
    character(len=pathLen) :: datafile
    integer :: iColor
    integer :: iProp
    integer :: colChar, colBit
    integer(kind=long_k), allocatable :: long_counts(:)
    integer :: iError
    integer :: iElem
    integer :: ice
    ! -----------------------------------------------------------------------!

    prp_loop: do iprop=1, tree%global%nProperties
      if (tree%global%Property(iprop)%bitpos == prp_hasPolynomial) then
        me%header => tree%global%Property(iprop)
        me%property => tree%property(iprop)

        datafile   = trim(tree%global%dirname)//'subres.ascii'

        allocate(me%subresolved_colors(coloring%nChars, me%property%nElems))
        allocate(me%nElems(coloring%nColors))
        allocate(me%Offset(coloring%nColors))
        me%Offset = 0_long_k

        ! If there are actually subresolved elements on the local process,
        ! read them now.
        if (me%property%nElems > 0) then

          allocate(me%elem(coloring%nColors))
          inquire(iolength=rl) me%subresolved_colors(:,1)
          call tem_open( newunit = fUnit,         &
            &            file    = datafile,      &
            &            action  = 'read',        &
            &            access  = 'stream',      &
            &            form    = 'unformatted', &
            &            status  = 'old'          )
          read(fUnit, pos=me%property%offset+1) me%subresolved_colors
          close(fUnit)

          do iColor=1,coloring%nColors
            colChar = (iColor-1)/colors_per_char + 1
            colBit = mod(iColor-1, colors_per_char)
            me%nElems(iColor)                                            &
              &  = count( btest(ichar(me%subresolved_colors(ColChar,:)), &
              &                 ColBit)                                  )
            allocate( me%elem(iColor)%id(me%nElems(iColor)) )

            ! Store the element link for each color.
            if (me%nElems(iColor) > 0) then
              ice = 0
              do iElem=1,me%property%nElems
                if ( btest(ichar(me%subresolved_colors(ColChar,iElem)), &
                  &        colBit) ) then
                  ice = ice + 1
                  me%elem(iColor)%id(ice) = me%property%elemid(iElem)
                end if
              end do
            end if

          end do

        else
          ! No subresolved elements on the local partition:
          me%nElems = 0
        end if

        allocate( long_counts(coloring%nColors) )
        long_counts = me%nElems
        call MPI_Exscan( long_counts, me%Offset, coloring%nColors,       &
          &              MPI_INTEGER8, MPI_SUM, tree%global%comm, iError )

        EXIT prp_loop

      end if
    end do prp_loop

  end subroutine tem_subres_prop_load