Writes ResPart data to file using H5HUT.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(ResPart_set), | intent(inout) | :: | this |
Set of resolved particles |
||
| integer, | intent(in) | :: | iter |
Iteration at write |
||
| real(kind=wp), | intent(in) | :: | time |
Time at write |
impure subroutine ResPart_set_WriteH5HUT(this,iter,time) !> Writes ResPart data to file using H5HUT. implicit none class(ResPart_set), intent(inout) :: this !! Set of resolved particles integer, intent(in) :: iter !! Iteration at write real(wp), intent(in) :: time !! Time at write ! Work variables type(h5hut_obj) :: h5 integer, allocatable :: buffi(:) integer(leapI8),allocatable :: buffi8(:) real(wp), allocatable :: buffr(:) ! Nothing to write, if empty call this%UpdateCount() if (this%count.eq.0) return ! Open the file if (this%overwrite) then call h5%initialize(trim(adjustl(this%write_file)),"W",this%parallel) else call h5%initialize(trim(adjustl(this%write_file)),"RW",this%parallel) end if ! Create a new time step call h5%NewTimeStep(iter,time) ! Allocate buffers allocate(buffi8(this%count_),buffi(this%count_),buffr(this%count_)) if (this%count_.eq.0) then call h5%Write("id", buffi8) call h5%Write("x", buffr ) call h5%Write("y", buffr ) call h5%Write("z", buffr ) call h5%Write("i", buffi ) call h5%Write("j", buffi ) call h5%Write("k", buffi ) call h5%Write("u", buffr ) call h5%Write("v", buffr ) call h5%Write("w", buffr ) call h5%Write("wx", buffr ) call h5%Write("wy", buffr ) call h5%Write("wz", buffr ) call h5%Write("d", buffr ) call h5%Write("s", buffi ) call h5%Write("rho", buffr ) call h5%Write("Fhx", buffr ) call h5%Write("Fhy", buffr ) call h5%Write("Fhz", buffr ) call h5%Write("Thx", buffr ) call h5%Write("Thy", buffr ) call h5%Write("Thz", buffr ) call h5%Write("Fcx", buffr ) call h5%Write("Fcy", buffr ) call h5%Write("Fcz", buffr ) call h5%Write("Tcx", buffr ) call h5%Write("Tcy", buffr ) call h5%Write("Tcz", buffr ) else ! Write the particle data select type (particles=>this%p) type is (ResPart_obj) buffi8=particles(1:this%count_)%id ; call h5%Write("id", buffi8) buffr =particles(1:this%count_)%p(1) ; call h5%Write("x", buffr ) buffr =particles(1:this%count_)%p(2) ; call h5%Write("y", buffr ) buffr =particles(1:this%count_)%p(3) ; call h5%Write("z", buffr ) buffi =particles(1:this%count_)%c(1) ; call h5%Write("i", buffi ) buffi =particles(1:this%count_)%c(2) ; call h5%Write("j", buffi ) buffi =particles(1:this%count_)%c(3) ; call h5%Write("k", buffi ) buffr =particles(1:this%count_)%v(1) ; call h5%Write("u", buffr ) buffr =particles(1:this%count_)%v(2) ; call h5%Write("v", buffr ) buffr =particles(1:this%count_)%v(3) ; call h5%Write("w", buffr ) buffr =particles(1:this%count_)%w(1) ; call h5%Write("wx", buffr ) buffr =particles(1:this%count_)%w(2) ; call h5%Write("wy", buffr ) buffr =particles(1:this%count_)%w(3) ; call h5%Write("wz", buffr ) buffr =particles(1:this%count_)%d ; call h5%Write("d", buffr ) buffi =particles(1:this%count_)%s ; call h5%Write("s", buffi ) buffr =particles(1:this%count_)%rho ; call h5%Write("rho", buffr ) buffr =particles(1:this%count_)%Fh(1); call h5%Write("Fhx", buffr ) buffr =particles(1:this%count_)%Fh(2); call h5%Write("Fhy", buffr ) buffr =particles(1:this%count_)%Fh(3); call h5%Write("Fhz", buffr ) buffr =particles(1:this%count_)%Th(1); call h5%Write("Thx", buffr ) buffr =particles(1:this%count_)%Th(2); call h5%Write("Thy", buffr ) buffr =particles(1:this%count_)%Th(3); call h5%Write("Thz", buffr ) buffr =particles(1:this%count_)%Fc(1); call h5%Write("Fcx", buffr ) buffr =particles(1:this%count_)%Fc(2); call h5%Write("Fcy", buffr ) buffr =particles(1:this%count_)%Fc(3); call h5%Write("Fcz", buffr ) buffr =particles(1:this%count_)%Tc(1); call h5%Write("Tcx", buffr ) buffr =particles(1:this%count_)%Tc(2); call h5%Write("Tcy", buffr ) buffr =particles(1:this%count_)%Tc(3); call h5%Write("Tcz", buffr ) end select end if ! Deallocate buffers deallocate(buffi8,buffi,buffr) call h5%finalize ! Write the associated IB data call this%ib%WriteH5HUT(iter,time) return end subroutine ResPart_set_WriteH5HUT