bc_set_Write Subroutine

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

Writes bc_set to disk using HDF5. The file structure follows this convention:

  / (root)
   !-- Time
   !-- Iter
   !-- Region 1
      !-- xlo
      !-- xhi
      !-- dir
      !-- side
      |-- Var 1
         |-- Type
         |-- Values(:,:,:)
      !-- Var 2
          .
          .
   !-- Region 2
       .
       .

Type Bound

bc_set

Arguments

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

Boundary conditions utility

integer, intent(in) :: iter

Iteration at write

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

Time at write

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

Name of file to write


Calls

proc~~bc_set_write~~CallsGraph proc~bc_set_write bc_set%bc_set_Write none~writeattributes hdf5_obj%WriteAttributes proc~bc_set_write->none~writeattributes proc~bc_set_getextents bc_set%bc_set_GetExtents proc~bc_set_write->proc~bc_set_getextents proc~hdf5_obj_creategroup hdf5_obj%hdf5_obj_CreateGroup proc~bc_set_write->proc~hdf5_obj_creategroup proc~hdf5_obj_final hdf5_obj%hdf5_obj_Final proc~bc_set_write->proc~hdf5_obj_final proc~hdf5_obj_init hdf5_obj%hdf5_obj_Init proc~bc_set_write->proc~hdf5_obj_init proc~hdf5_obj_writeattributes0d hdf5_obj%hdf5_obj_WriteAttributes0D none~writeattributes->proc~hdf5_obj_writeattributes0d proc~hdf5_obj_writeattributes1d hdf5_obj%hdf5_obj_WriteAttributes1D none~writeattributes->proc~hdf5_obj_writeattributes1d proc~bc_set_getregionindex bc_set%bc_set_GetRegionIndex proc~bc_set_getextents->proc~bc_set_getregionindex h5gcreate_f h5gcreate_f proc~hdf5_obj_creategroup->h5gcreate_f proc~hashtbl_obj_hashstring hashtbl_obj%hashtbl_obj_HashString proc~hdf5_obj_creategroup->proc~hashtbl_obj_hashstring proc~hashtbl_obj_put hashtbl_obj%hashtbl_obj_Put proc~hdf5_obj_creategroup->proc~hashtbl_obj_put proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_creategroup->proc~hdf5_obj_closegroup proc~hdf5_obj_fixgroupname hdf5_obj%hdf5_obj_FixGroupName proc~hdf5_obj_creategroup->proc~hdf5_obj_fixgroupname 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 proc~bc_set_getregionindex->proc~hashtbl_obj_hashstring none~get~3 hashtbl_obj%Get proc~bc_set_getregionindex->none~get~3 proc~sllist_obj_put sllist_obj%sllist_obj_Put proc~hashtbl_obj_put->proc~sllist_obj_put proc~hdf5_obj_closegroup->proc~hashtbl_obj_hashstring proc~hdf5_obj_closegroup->proc~hdf5_obj_fixgroupname 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 hdf5_obj%hdf5_obj_GetGroupObject proc~hdf5_obj_closegroup->proc~hdf5_obj_getgroupobject proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_closegroup h5aclose_f h5aclose_f proc~hdf5_obj_writeattributes0d->h5aclose_f h5acreate_f h5acreate_f proc~hdf5_obj_writeattributes0d->h5acreate_f h5awrite_f h5awrite_f proc~hdf5_obj_writeattributes0d->h5awrite_f h5sclose_f h5sclose_f proc~hdf5_obj_writeattributes0d->h5sclose_f h5screate_simple_f h5screate_simple_f proc~hdf5_obj_writeattributes0d->h5screate_simple_f proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_opengroup hdf5_obj%hdf5_obj_OpenGroup proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_opengroup proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_closegroup proc~hdf5_obj_writeattributes1d->h5aclose_f proc~hdf5_obj_writeattributes1d->h5acreate_f proc~hdf5_obj_writeattributes1d->h5awrite_f proc~hdf5_obj_writeattributes1d->h5sclose_f proc~hdf5_obj_writeattributes1d->h5screate_simple_f proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_opengroup 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~sllist_obj_remove sllist_obj%sllist_obj_Remove proc~hashtbl_obj_remove->proc~sllist_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_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_Write(this,iter,time,name)
      !> Writes bc_set to disk using HDF5.
      ! The file structure follows this convention:
      !
      !       / (root)
      !        !-- Time
      !        !-- Iter
      !        !-- Region 1
      !           !-- xlo
      !           !-- xhi
      !           !-- dir
      !           !-- side
      !           |-- Var 1
      !              |-- Type
      !              |-- Values(:,:,:)
      !           !-- Var 2
      !               .
      !               .
      !        !-- Region 2
      !            .
      !            .
      !
      implicit none
      class(bc_set),    intent(inout) :: this                                  !! Boundary conditions utility
      integer,          intent(in)    :: iter                                  !! Iteration at write
      real(wp),         intent(in)    :: time                                  !! Time at write
      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=:), &
               allocatable :: groupname
      character(len=:), &
               allocatable :: varname
      character(len=:), &
               allocatable :: filename
      integer              :: n,m
      real(wp)             :: xlo(3)
      real(wp)             :: xhi(3)

      ! 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,'W')

      ! Write time and iteration
      call hdf5%WriteAttributes('/','Iteration',    iter      )
      call hdf5%WriteAttributes('/','Time',         time      )

      do n=1,this%count
        ! Create new group for this region
        groupname = this%region(n)%name
        call hdf5%CreateGroup(groupname)

        ! Write extents
        extents = this%GetExtents(groupname)
        xlo     = this%region(n)%xlo
        xhi     = this%region(n)%xhi
        call hdf5%WriteAttributes(groupname,'xlo',xlo)
        call hdf5%WriteAttributes(groupname,'xhi',xhi)
        call hdf5%WriteAttributes(groupname,'dir', this%region(n)%dir )
        call hdf5%WriteAttributes(groupname,'side',this%region(n)%side)

        ! Write boundary conditions
        do m=1,size(this%region(n)%BC)
          ! Store variable name
          varname = trim(adjustl(this%region(n)%BC(m)%name))

          ! Create group for this variable
          call hdf5%CreateGroup(trim(adjustl(groupname))//'/'//varname)

          ! Write BC type
          call hdf5%WriteAttributes(trim(adjustl(groupname))//'/'//varname,'type',this%region(n)%BC(m)%type)

          select case (this%region(n)%BC(m)%type)
            case (BC_WALL     )
              ! Nothing to do
            case (BC_INFLOW   )
              ! Nothing to do
            case (BC_OUTFLOW  )
              ! Nothing to do
            case (BC_DIRICHLET)
              call hdf5%Write((trim(adjustl(groupname))//'/'//varname),'values',this%region(n)%BC(m)%val,extents%lo,extents%hi)
            case (BC_NEUMANN  )
              call hdf5%Write((trim(adjustl(groupname))//'/'//varname),'values',this%region(n)%BC(m)%val,extents%lo,extents%hi)
            case (BC_SYMMETRY )
              ! Nothing to do
          end select
        end do
      end do

      ! 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_Write