bc_set_Read Subroutine

private impure subroutine bc_set_Read(this, iter, time, name)

Reads boundary conditions from file.

Type Bound

bc_set

Arguments

Type IntentOptional Attributes Name
class(bc_set), intent(inout) :: this

Boundary conditions utility

integer, intent(out) :: iter

Iteration read from file

real(kind=wp), intent(out) :: time

Time read from file

character(len=*), intent(in), optional :: name

Name of file to write


Calls

proc~~bc_set_read~~CallsGraph proc~bc_set_read bc_set%bc_set_Read none~readattributes hdf5_obj%ReadAttributes proc~bc_set_read->none~readattributes proc~bc_set_add bc_set%bc_set_Add proc~bc_set_read->proc~bc_set_add proc~bc_set_getbcpointer bc_set%bc_set_GetBCPointer proc~bc_set_read->proc~bc_set_getbcpointer proc~bc_set_getextents bc_set%bc_set_GetExtents proc~bc_set_read->proc~bc_set_getextents proc~hdf5_obj_final hdf5_obj%hdf5_obj_Final proc~bc_set_read->proc~hdf5_obj_final proc~hdf5_obj_init hdf5_obj%hdf5_obj_Init proc~bc_set_read->proc~hdf5_obj_init proc~hdf5_obj_readgroupnames hdf5_obj%hdf5_obj_ReadGroupNames proc~bc_set_read->proc~hdf5_obj_readgroupnames proc~hdf5_obj_readattributes0d hdf5_obj%hdf5_obj_ReadAttributes0D none~readattributes->proc~hdf5_obj_readattributes0d proc~hdf5_obj_readattributes1d hdf5_obj%hdf5_obj_ReadAttributes1D none~readattributes->proc~hdf5_obj_readattributes1d proc~bc_set_checkbounds bc_set%bc_set_CheckBounds proc~bc_set_add->proc~bc_set_checkbounds proc~bc_set_expand bc_set%bc_set_Expand proc~bc_set_add->proc~bc_set_expand proc~bc_set_getsidedirbynormal bc_set%bc_set_GetSideDirByNormal proc~bc_set_add->proc~bc_set_getsidedirbynormal proc~bc_set_updateextents bc_set%bc_set_UpdateExtents proc~bc_set_add->proc~bc_set_updateextents proc~hashtbl_obj_hashstring hashtbl_obj%hashtbl_obj_HashString proc~bc_set_add->proc~hashtbl_obj_hashstring proc~hashtbl_obj_put hashtbl_obj%hashtbl_obj_Put proc~bc_set_add->proc~hashtbl_obj_put proc~region_obj_init region_obj%region_obj_Init proc~bc_set_add->proc~region_obj_init proc~bc_set_getregionindex bc_set%bc_set_GetRegionIndex proc~bc_set_getbcpointer->proc~bc_set_getregionindex proc~region_obj_getbcindex region_obj%region_obj_GetBCIndex proc~bc_set_getbcpointer->proc~region_obj_getbcindex proc~bc_set_getextents->proc~bc_set_getregionindex h5garbage_collect_f h5garbage_collect_f proc~hdf5_obj_final->h5garbage_collect_f proc~hashtbl_obj_final hashtbl_obj%hashtbl_obj_Final proc~hdf5_obj_final->proc~hashtbl_obj_final h5open_f h5open_f proc~hdf5_obj_init->h5open_f proc~hashtbl_obj_init hashtbl_obj%hashtbl_obj_Init proc~hdf5_obj_init->proc~hashtbl_obj_init h5gn_members_f h5gn_members_f proc~hdf5_obj_readgroupnames->h5gn_members_f h5iget_name_f h5iget_name_f proc~hdf5_obj_readgroupnames->h5iget_name_f h5oclose_f h5oclose_f proc~hdf5_obj_readgroupnames->h5oclose_f h5oget_info_f h5oget_info_f proc~hdf5_obj_readgroupnames->h5oget_info_f h5oopen_by_idx_f h5oopen_by_idx_f proc~hdf5_obj_readgroupnames->h5oopen_by_idx_f proc~hdf5_obj_fixgroupname hdf5_obj%hdf5_obj_FixGroupName proc~hdf5_obj_readgroupnames->proc~hdf5_obj_fixgroupname proc~bc_set_getregionindex->proc~hashtbl_obj_hashstring none~get~3 hashtbl_obj%Get proc~bc_set_getregionindex->none~get~3 proc~bc_set_updateextents->proc~bc_set_getregionindex proc~sllist_obj_put sllist_obj%sllist_obj_Put proc~hashtbl_obj_put->proc~sllist_obj_put h5aclose_f h5aclose_f proc~hdf5_obj_readattributes0d->h5aclose_f h5aget_space_f h5aget_space_f proc~hdf5_obj_readattributes0d->h5aget_space_f h5aget_type_f h5aget_type_f proc~hdf5_obj_readattributes0d->h5aget_type_f h5aopen_f h5aopen_f proc~hdf5_obj_readattributes0d->h5aopen_f h5aread_f h5aread_f proc~hdf5_obj_readattributes0d->h5aread_f h5sclose_f h5sclose_f proc~hdf5_obj_readattributes0d->h5sclose_f h5sget_simple_extent_dims_f h5sget_simple_extent_dims_f proc~hdf5_obj_readattributes0d->h5sget_simple_extent_dims_f proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_readattributes0d->proc~hdf5_obj_closegroup proc~hdf5_obj_getgroupobject hdf5_obj%hdf5_obj_GetGroupObject proc~hdf5_obj_readattributes0d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_opengroup hdf5_obj%hdf5_obj_OpenGroup proc~hdf5_obj_readattributes0d->proc~hdf5_obj_opengroup proc~hdf5_obj_readattributes1d->h5aclose_f proc~hdf5_obj_readattributes1d->h5aget_space_f proc~hdf5_obj_readattributes1d->h5aget_type_f proc~hdf5_obj_readattributes1d->h5aopen_f proc~hdf5_obj_readattributes1d->h5aread_f proc~hdf5_obj_readattributes1d->h5sclose_f proc~hdf5_obj_readattributes1d->h5sget_simple_extent_dims_f proc~hdf5_obj_readattributes1d->proc~hdf5_obj_closegroup proc~hdf5_obj_readattributes1d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_readattributes1d->proc~hdf5_obj_opengroup proc~region_obj_getbcindex->proc~hashtbl_obj_hashstring proc~region_obj_getbcindex->none~get~3 proc~region_obj_init->proc~hashtbl_obj_init proc~hashtbl_obj_get_int4 hashtbl_obj%hashtbl_obj_Get_int4 none~get~3->proc~hashtbl_obj_get_int4 proc~hashtbl_obj_get_int8 hashtbl_obj%hashtbl_obj_Get_int8 none~get~3->proc~hashtbl_obj_get_int8 proc~hashtbl_obj_get_real_dp hashtbl_obj%hashtbl_obj_Get_real_dp none~get~3->proc~hashtbl_obj_get_real_dp proc~hashtbl_obj_get_real_sp hashtbl_obj%hashtbl_obj_Get_real_sp none~get~3->proc~hashtbl_obj_get_real_sp proc~hdf5_obj_closegroup->proc~hashtbl_obj_hashstring proc~hdf5_obj_closegroup->proc~hdf5_obj_fixgroupname proc~hdf5_obj_closegroup->proc~hdf5_obj_getgroupobject h5gclose_f h5gclose_f proc~hdf5_obj_closegroup->h5gclose_f proc~hashtbl_obj_remove hashtbl_obj%hashtbl_obj_Remove proc~hdf5_obj_closegroup->proc~hashtbl_obj_remove proc~hdf5_obj_getgroupobject->proc~hashtbl_obj_hashstring proc~hdf5_obj_getgroupobject->proc~hdf5_obj_fixgroupname proc~hdf5_obj_getgroupobject->none~get~3 proc~hdf5_obj_opengroup->proc~hashtbl_obj_hashstring proc~hdf5_obj_opengroup->proc~hashtbl_obj_put proc~hdf5_obj_opengroup->proc~hdf5_obj_fixgroupname proc~hdf5_obj_opengroup->proc~hdf5_obj_getgroupobject h5oopen_f h5oopen_f proc~hdf5_obj_opengroup->h5oopen_f proc~sllist_obj_put->proc~sllist_obj_put none~get~2 sllist_obj%Get proc~hashtbl_obj_get_int4->none~get~2 proc~hashtbl_obj_get_int8->none~get~2 proc~hashtbl_obj_get_real_dp->none~get~2 proc~hashtbl_obj_get_real_sp->none~get~2 proc~sllist_obj_remove sllist_obj%sllist_obj_Remove proc~hashtbl_obj_remove->proc~sllist_obj_remove proc~sllist_obj_get_int4 sllist_obj%sllist_obj_Get_int4 none~get~2->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8 sllist_obj%sllist_obj_Get_int8 none~get~2->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp sllist_obj%sllist_obj_Get_real_dp none~get~2->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp sllist_obj%sllist_obj_Get_real_sp none~get~2->proc~sllist_obj_get_real_sp proc~sllist_obj_get_int4->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp->proc~sllist_obj_get_real_sp

Source Code

    impure subroutine bc_set_Read(this,iter,time,name)
      !> Reads boundary conditions from file.
      implicit none
      class(bc_set),    intent(inout) :: this                                  !! Boundary conditions utility
      integer,          intent(out)   :: iter                                  !! Iteration read from file
      real(wp),         intent(out)   :: time                                  !! Time read from file
      character(len=*), intent(in),    &
                             optional :: name                                  !! Name of file to write
      ! Work variables
      type(hdf5_obj),   &
               pointer     :: hdf5
      type(hdf5_obj),   &
              target       :: hdf5t
      type(extent_obj)     :: extents
      character(len=str64),&
                allocatable:: gnames(:)
      character(len=str64),&
                allocatable:: vnames(:)
      character(len=:), &
                allocatable:: group_
      character(len=:), &
               allocatable :: filename
      integer              :: vtype
      integer              :: n,m

      real(wp)             :: xlo(3)
      real(wp)             :: xhi(3)
      integer              :: side
      integer              :: dir
      character(len=3)     :: normal
      real(wp), pointer    :: dataptr(:,:,:)

      ! Determine the filename
      if (present(name)) then
        filename=trim(adjustl(name))
      else
        filename='bcs'
      end if

      ! Use internal pointer, if associated.
      ! Otherwise, use a new object.
      if (associated(this%hdf5)) then
        hdf5 => this%hdf5
      else
        call hdf5t%Initialize(this%parallel)
        hdf5 => hdf5t
      end if

      ! Open file
      call hdf5%Open(filename,'R')

      ! Read iteration and time
      call hdf5%ReadAttributes('/','Iteration',    iter      )
      call hdf5%ReadAttributes('/','Time',         time      )

      ! Read groups under root
      call hdf5%ReadGroupNames('/',gnames)

      if (allocated(gnames)) then
        do n=1,size(gnames)

          ! Read extents of region
          call hdf5%ReadAttributes(trim(adjustl(gnames(n))),'xlo',  xlo )
          call hdf5%ReadAttributes(trim(adjustl(gnames(n))),'xhi',  xhi )
          call hdf5%ReadAttributes(trim(adjustl(gnames(n))),'dir',  dir )
          call hdf5%ReadAttributes(trim(adjustl(gnames(n))),'side', side)

          ! Rebuild normal
          if (side.eq.BC_RIGHT) normal(1:1) = '-'
          if (side.eq.BC_LEFT ) normal(1:1) = '+'
          write(normal(2:3),fmt='(a1,i1)')'x',dir

          ! Add region
          call this%Add(trim(adjustl(gnames(n))),xlo,xhi,normal)

          ! Get region extents
          extents = this%GetExtents(trim(adjustl(gnames(n))))

          ! Read groups under this group
          call hdf5%ReadGroupNames('/'//trim(adjustl(gnames(n)))//'/',vnames)

          do m=1,size(vnames)
            ! Path to this variable
            group_ = trim(adjustl(gnames(n)))//'/'//trim(adjustl(vnames(m)))

            ! Read type
            call hdf5%ReadAttributes(group_,'type',vtype)

            ! Read values
            select case (vtype)
              case (BC_WALL     )
                ! Nothing to do
              case (BC_INFLOW   )
                ! Nothing to do
              case (BC_OUTFLOW  )
                ! Nothing to do
              case (BC_DIRICHLET)
                ! Create boundary condition
                call this%region(n)%Add(trim(adjustl(vnames(m))), BC_DIRICHLET, extents)
                ! Get a pointer to data values
                call this%GetBCPointer(trim(adjustl(gnames(n))),trim(adjustl(vnames(m))),dataptr)
                ! Read values
                call hdf5%Read(group_,'values',dataptr,extents%lo,extents%hi)
              case (BC_NEUMANN  )
                ! Create boundary condition
                call this%region(n)%Add(trim(adjustl(vnames(m))), BC_NEUMANN,   extents)
                ! Get a pointer to data values
                call this%GetBCPointer(trim(adjustl(gnames(n))),trim(adjustl(vnames(m))),dataptr)
                ! Read values
                call hdf5%Read(group_,'values',dataptr,extents%lo,extents%hi)
              case (BC_SYMMETRY )
                ! Create boundary condition
                call this%region(n)%Add(trim(adjustl(vnames(m))), BC_SYMMETRY, extents)
            end select
          end do

        end do
      end if

      ! Close file
      call hdf5%Close()

      ! Free data
      if (associated(this%hdf5)) then
        hdf5 => null()
      else
        call hdf5t%Finalize()
        hdf5 => null()
      end if

      return
    end subroutine bc_set_Read