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 | Intent | Optional | 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 |
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