Writes particle data to file in parallel using HDF5.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(particle_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 particle_set_WriteHDF5(this,iter,time) !> Writes particle data to file in parallel using HDF5. use leapUtils, only: stringtool_obj, sysutils_obj implicit none class(particle_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(sysutils_obj) :: sysutils type(stringtool_obj) :: stringtool character(len=:), & allocatable :: filename character(len=str64) :: tmpC type(hdf5_obj) :: hdf5 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 ! Create directory if (this%parallel%RankIsRoot()) & call sysutils%CreateDirectory('Hdf5') filename = 'Hdf5/'//stringtool%RemoveExtension(this%write_file) ! Initialize HDF5 call hdf5%Initialize(this%parallel) ! Open the file if (this%overwrite) then write(tmpC,'(ES12.4)') time filename = filename//"_"//trim(adjustl(tmpC)) call hdf5%Open(filename,"W") else call hdf5%Open(filename,"RW") end if ! Write attributes call hdf5%WriteAttributes("/","Time", time) call hdf5%WriteAttributes("/","Iter", iter) call hdf5%WriteAttributes("/","Count",this%count) ! Allocate buffers allocate(buffi8(this%count_),buffi(this%count_),buffr(this%count_)) ! Write the particle data ! Write the particle data if (this%count_.eq.0) then call hdf5%Write("/","id", buffi8) call hdf5%Write("/","x", buffr ) call hdf5%Write("/","y", buffr ) call hdf5%Write("/","z", buffr ) call hdf5%Write("/","i", buffi ) call hdf5%Write("/","j", buffi ) call hdf5%Write("/","k", buffi ) call hdf5%Write("/","vx", buffr ) call hdf5%Write("/","vy", buffr ) call hdf5%Write("/","vz", buffr ) call hdf5%Write("/","wx", buffr ) call hdf5%Write("/","wy", buffr ) call hdf5%Write("/","wz", buffr ) call hdf5%Write("/","d", buffr ) call hdf5%Write("/","s", buffi ) call hdf5%Write("/","rho",buffr ) call hdf5%Write("/","Fhx",buffr ) call hdf5%Write("/","Fhy",buffr ) call hdf5%Write("/","Fhz",buffr ) call hdf5%Write("/","Thx",buffr ) call hdf5%Write("/","Thy",buffr ) call hdf5%Write("/","Thz",buffr ) call hdf5%Write("/","Fcx",buffr ) call hdf5%Write("/","Fcy",buffr ) call hdf5%Write("/","Fcz",buffr ) call hdf5%Write("/","Tcx",buffr ) call hdf5%Write("/","Tcy",buffr ) call hdf5%Write("/","Tcz",buffr ) else ! Write the particle data select type (particles=>this%p) type is (particle_obj) buffi8=particles(1:this%count_)%id ; call hdf5%Write("/","id", buffi8) buffr =particles(1:this%count_)%p(1) ; call hdf5%Write("/","x", buffr ) buffr =particles(1:this%count_)%p(2) ; call hdf5%Write("/","y", buffr ) buffr =particles(1:this%count_)%p(3) ; call hdf5%Write("/","z", buffr ) buffi =particles(1:this%count_)%c(1) ; call hdf5%Write("/","i", buffi ) buffi =particles(1:this%count_)%c(2) ; call hdf5%Write("/","j", buffi ) buffi =particles(1:this%count_)%c(3) ; call hdf5%Write("/","k", buffi ) buffr =particles(1:this%count_)%v(1) ; call hdf5%Write("/","vx", buffr ) buffr =particles(1:this%count_)%v(2) ; call hdf5%Write("/","vy", buffr ) buffr =particles(1:this%count_)%v(3) ; call hdf5%Write("/","vz", buffr ) buffr =particles(1:this%count_)%w(1) ; call hdf5%Write("/","wx", buffr ) buffr =particles(1:this%count_)%w(2) ; call hdf5%Write("/","wy", buffr ) buffr =particles(1:this%count_)%w(3) ; call hdf5%Write("/","wz", buffr ) buffr =particles(1:this%count_)%d ; call hdf5%Write("/","d", buffr ) buffi =particles(1:this%count_)%s ; call hdf5%Write("/","s", buffi ) buffr =particles(1:this%count_)%rho ; call hdf5%Write("/","rho",buffr ) buffr =particles(1:this%count_)%Fh(1); call hdf5%Write("/","Fhx",buffr ) buffr =particles(1:this%count_)%Fh(2); call hdf5%Write("/","Fhy",buffr ) buffr =particles(1:this%count_)%Fh(3); call hdf5%Write("/","Fhz",buffr ) buffr =particles(1:this%count_)%Th(1); call hdf5%Write("/","Thx",buffr ) buffr =particles(1:this%count_)%Th(2); call hdf5%Write("/","Thy",buffr ) buffr =particles(1:this%count_)%Th(3); call hdf5%Write("/","Thz",buffr ) buffr =particles(1:this%count_)%Fc(1); call hdf5%Write("/","Fcx",buffr ) buffr =particles(1:this%count_)%Fc(2); call hdf5%Write("/","Fcy",buffr ) buffr =particles(1:this%count_)%Fc(3); call hdf5%Write("/","Fcz",buffr ) buffr =particles(1:this%count_)%Tc(1); call hdf5%Write("/","Tcx",buffr ) buffr =particles(1:this%count_)%Tc(2); call hdf5%Write("/","Tcy",buffr ) buffr =particles(1:this%count_)%Tc(3); call hdf5%Write("/","Tcz",buffr ) end select end if ! Deallocate buffers deallocate(buffi8,buffi,buffr) ! Close file call hdf5%Close() ! Finalize HDF5 call hdf5%Finalize() return end subroutine particle_set_WriteHDF5