Reads ResPart data from file using H5HUT.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(ResPart_set), | intent(inout) | :: | this |
Set of resolved particles |
||
| integer, | intent(out) | :: | iter |
Iteration at write |
||
| real(kind=wp), | intent(out) | :: | time |
Time at write |
||
| integer, | intent(in), | optional | :: | step |
User supplied step to open |
impure subroutine ResPart_set_ReadH5HUT(this,iter,time,step) !> Reads ResPart data from file using H5HUT. implicit none class(ResPart_set), intent(inout) :: this !! Set of resolved particles integer, intent(out) :: iter !! Iteration at write real(wp), intent(out) :: time !! Time at write integer, intent(in), & optional :: step !! User supplied step to open ! Work variables type(h5hut_obj) :: h5 integer, allocatable :: buffi(:) integer(leapI8),allocatable :: buffi8(:) real(wp), allocatable :: buffr(:) integer :: n ! Open the file call h5%initialize(trim(adjustl(this%read_file)),"R",this%parallel) if (.not.present(step)) then ! Jump to last step call h5%LastTimeStep(iter,time) else ! Jump to desired step call h5%JumpToStep(step,iter,time) end if ! Get total number of point objects call h5%getnpoints(this%count) ! Distribute equally the data on each mpi rank associate(mpi=>this%parallel) this%count_=int(this%count/mpi%nproc) if (mpi%rank%mine.le.mod(this%count,mpi%nproc)) this%count_=this%count_+1 end associate ! Resize the data array to contain read Lagrangian objects call this%Resize(this%count_) ! Allocate buffers allocate(buffi8(this%count_),buffi(this%count_),buffr(this%count_)) ! Read the base Lagrangian object structure if (this%count_.eq.0) then call h5%Read("id", buffi8) call h5%Read("x", buffr ) call h5%Read("y", buffr ) call h5%Read("z", buffr ) call h5%Read("i", buffi ) call h5%Read("j", buffi ) call h5%Read("k", buffi ) call h5%Read("u", buffr ) call h5%Read("v", buffr ) call h5%Read("w", buffr ) call h5%Read("wx", buffr ) call h5%Read("wy", buffr ) call h5%Read("wz", buffr ) call h5%Read("d", buffr ) call h5%Read("s", buffi ) call h5%Read("rho", buffr ) call h5%Read("Fhx", buffr ) call h5%Read("Fhy", buffr ) call h5%Read("Fhz", buffr ) call h5%Read("Thx", buffr ) call h5%Read("Thy", buffr ) call h5%Read("Thz", buffr ) call h5%Read("Fcx", buffr ) call h5%Read("Fcy", buffr ) call h5%Read("Fcz", buffr ) call h5%Read("Tcx", buffr ) call h5%Read("Tcy", buffr ) call h5%Read("Tcz", buffr ) else select type (particles=>this%p) type is (ResPart_obj) call h5%Read("id", buffi8); particles(:)%id = buffi8 call h5%Read("x", buffr ); particles(:)%p(1) = buffr call h5%Read("y", buffr ); particles(:)%p(2) = buffr call h5%Read("z", buffr ); particles(:)%p(3) = buffr call h5%Read("i", buffi ); particles(:)%c(1) = buffi call h5%Read("j", buffi ); particles(:)%c(2) = buffi call h5%Read("k", buffi ); particles(:)%c(3) = buffi call h5%Read("u", buffr ); particles(:)%v(1) = buffr call h5%Read("v", buffr ); particles(:)%v(2) = buffr call h5%Read("w", buffr ); particles(:)%v(3) = buffr call h5%Read("wx", buffr ); particles(:)%w(1) = buffr call h5%Read("wy", buffr ); particles(:)%w(2) = buffr call h5%Read("wz", buffr ); particles(:)%w(3) = buffr call h5%Read("d", buffr ); particles(:)%d = buffr call h5%Read("s", buffi ); particles(:)%s = buffi call h5%Read("rho", buffr ); particles(:)%rho = buffr call h5%Read("Fhx", buffr ); particles(:)%Fh(1)= buffr call h5%Read("Fhy", buffr ); particles(:)%Fh(2)= buffr call h5%Read("Fhz", buffr ); particles(:)%Fh(3)= buffr call h5%Read("Thx", buffr ); particles(:)%Th(1)= buffr call h5%Read("Thy", buffr ); particles(:)%Th(2)= buffr call h5%Read("Thz", buffr ); particles(:)%Th(3)= buffr call h5%Read("Fcx", buffr ); particles(:)%Fc(1)= buffr call h5%Read("Fcy", buffr ); particles(:)%Fc(2)= buffr call h5%Read("Fcz", buffr ); particles(:)%Fc(3)= buffr call h5%Read("Tcx", buffr ); particles(:)%Tc(1)= buffr call h5%Read("Tcy", buffr ); particles(:)%Tc(2)= buffr call h5%Read("Tcz", buffr ); particles(:)%Tc(3)= buffr end select end if ! Free data deallocate(buffi8,buffi,buffr) call h5%Finalize ! Set unread variables to zero select type (particles=>this%p) class is (ResPart_obj) do n=1,this%count_ particles(n)%pold = 0.0_wp particles(n)%vold = 0.0_wp particles(n)%wold = 0.0_wp particles(n)%Fhold = 0.0_wp particles(n)%Thold = 0.0_wp particles(n)%Fcold = 0.0_wp particles(n)%Tcold = 0.0_wp end do end select ! Send Lagrangian centroids to their mpiranks call this%Communicate ! Localize markers on the grid call this%Localize ! Read the associated IB data if (.not.present(step)) then call this%ib%Read(iter,time) else call this%ib%Read(iter,time,step) end if ! Update lookup table call this%UpdateLookup return end subroutine ResPart_set_ReadH5HUT