Reads boundary conditions from file.
| Type | Intent | Optional | 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 |
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