ResPart_set_WriteHDF5 Subroutine

private impure subroutine ResPart_set_WriteHDF5(this, iter, time)

Uses

  • proc~~respart_set_writehdf5~~UsesGraph proc~respart_set_writehdf5 ResPart_set%ResPart_set_WriteHDF5 module~leaputils leapUtils proc~respart_set_writehdf5->module~leaputils module~leapkinds leapKinds module~leaputils->module~leapkinds iso_fortran_env iso_fortran_env module~leapkinds->iso_fortran_env

Writes ResPart data to file using HDF5.

Type Bound

ResPart_set

Arguments

Type IntentOptional 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


Calls

proc~~respart_set_writehdf5~~CallsGraph proc~respart_set_writehdf5 ResPart_set%ResPart_set_WriteHDF5 c c proc~respart_set_writehdf5->c fc fc proc~respart_set_writehdf5->fc fh fh proc~respart_set_writehdf5->fh none~writeattributes hdf5_obj%WriteAttributes proc~respart_set_writehdf5->none~writeattributes p p proc~respart_set_writehdf5->p particles particles proc~respart_set_writehdf5->particles proc~hdf5_obj_final hdf5_obj%hdf5_obj_Final proc~respart_set_writehdf5->proc~hdf5_obj_final proc~lagrangian_set_updatecount lagrangian_set%lagrangian_set_UpdateCount proc~respart_set_writehdf5->proc~lagrangian_set_updatecount proc~marker_set_writehdf5 marker_set%marker_set_WriteHDF5 proc~respart_set_writehdf5->proc~marker_set_writehdf5 proc~parallel_obj_rankisroot parallel_obj%parallel_obj_RankIsRoot proc~respart_set_writehdf5->proc~parallel_obj_rankisroot proc~stringtool_obj_removeextension stringtool_obj%stringtool_obj_RemoveExtension proc~respart_set_writehdf5->proc~stringtool_obj_removeextension proc~sysutils_obj_createdirectory sysutils_obj%sysutils_obj_CreateDirectory proc~respart_set_writehdf5->proc~sysutils_obj_createdirectory tc tc proc~respart_set_writehdf5->tc th th proc~respart_set_writehdf5->th v v proc~respart_set_writehdf5->v w w proc~respart_set_writehdf5->w proc~hdf5_obj_writeattributes0d hdf5_obj%hdf5_obj_WriteAttributes0D none~writeattributes->proc~hdf5_obj_writeattributes0d proc~hdf5_obj_writeattributes1d hdf5_obj%hdf5_obj_WriteAttributes1D none~writeattributes->proc~hdf5_obj_writeattributes1d h5garbage_collect_f h5garbage_collect_f proc~hdf5_obj_final->h5garbage_collect_f proc~hashtbl_obj_final hashtbl_obj%hashtbl_obj_Final proc~hdf5_obj_final->proc~hashtbl_obj_final mpi_allgather mpi_allgather proc~lagrangian_set_updatecount->mpi_allgather proc~marker_set_writehdf5->c proc~marker_set_writehdf5->none~writeattributes proc~marker_set_writehdf5->p proc~marker_set_writehdf5->proc~hdf5_obj_final proc~marker_set_writehdf5->proc~lagrangian_set_updatecount proc~marker_set_writehdf5->proc~parallel_obj_rankisroot proc~marker_set_writehdf5->proc~stringtool_obj_removeextension proc~marker_set_writehdf5->proc~sysutils_obj_createdirectory proc~marker_set_writehdf5->v f f proc~marker_set_writehdf5->f markers markers proc~marker_set_writehdf5->markers n n proc~marker_set_writehdf5->n proc~hdf5_obj_init hdf5_obj%hdf5_obj_Init proc~marker_set_writehdf5->proc~hdf5_obj_init h5open_f h5open_f proc~hdf5_obj_init->h5open_f proc~hashtbl_obj_init hashtbl_obj%hashtbl_obj_Init proc~hdf5_obj_init->proc~hashtbl_obj_init h5aclose_f h5aclose_f proc~hdf5_obj_writeattributes0d->h5aclose_f h5acreate_f h5acreate_f proc~hdf5_obj_writeattributes0d->h5acreate_f h5awrite_f h5awrite_f proc~hdf5_obj_writeattributes0d->h5awrite_f h5sclose_f h5sclose_f proc~hdf5_obj_writeattributes0d->h5sclose_f h5screate_simple_f h5screate_simple_f proc~hdf5_obj_writeattributes0d->h5screate_simple_f proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_closegroup proc~hdf5_obj_getgroupobject hdf5_obj%hdf5_obj_GetGroupObject proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_opengroup hdf5_obj%hdf5_obj_OpenGroup proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_opengroup proc~hdf5_obj_writeattributes1d->h5aclose_f proc~hdf5_obj_writeattributes1d->h5acreate_f proc~hdf5_obj_writeattributes1d->h5awrite_f proc~hdf5_obj_writeattributes1d->h5sclose_f proc~hdf5_obj_writeattributes1d->h5screate_simple_f proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_closegroup proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_opengroup proc~hdf5_obj_closegroup->proc~hdf5_obj_getgroupobject h5gclose_f h5gclose_f proc~hdf5_obj_closegroup->h5gclose_f proc~hashtbl_obj_hashstring hashtbl_obj%hashtbl_obj_HashString proc~hdf5_obj_closegroup->proc~hashtbl_obj_hashstring proc~hashtbl_obj_remove hashtbl_obj%hashtbl_obj_Remove proc~hdf5_obj_closegroup->proc~hashtbl_obj_remove proc~hdf5_obj_fixgroupname hdf5_obj%hdf5_obj_FixGroupName proc~hdf5_obj_closegroup->proc~hdf5_obj_fixgroupname none~get~3 hashtbl_obj%Get proc~hdf5_obj_getgroupobject->none~get~3 proc~hdf5_obj_getgroupobject->proc~hashtbl_obj_hashstring proc~hdf5_obj_getgroupobject->proc~hdf5_obj_fixgroupname proc~hdf5_obj_opengroup->proc~hdf5_obj_getgroupobject h5oopen_f h5oopen_f proc~hdf5_obj_opengroup->h5oopen_f proc~hdf5_obj_opengroup->proc~hashtbl_obj_hashstring proc~hashtbl_obj_put hashtbl_obj%hashtbl_obj_Put proc~hdf5_obj_opengroup->proc~hashtbl_obj_put proc~hdf5_obj_opengroup->proc~hdf5_obj_fixgroupname proc~hashtbl_obj_get_int4 hashtbl_obj%hashtbl_obj_Get_int4 none~get~3->proc~hashtbl_obj_get_int4 proc~hashtbl_obj_get_int8 hashtbl_obj%hashtbl_obj_Get_int8 none~get~3->proc~hashtbl_obj_get_int8 proc~hashtbl_obj_get_real_dp hashtbl_obj%hashtbl_obj_Get_real_dp none~get~3->proc~hashtbl_obj_get_real_dp proc~hashtbl_obj_get_real_sp hashtbl_obj%hashtbl_obj_Get_real_sp none~get~3->proc~hashtbl_obj_get_real_sp proc~sllist_obj_put sllist_obj%sllist_obj_Put proc~hashtbl_obj_put->proc~sllist_obj_put proc~sllist_obj_remove sllist_obj%sllist_obj_Remove proc~hashtbl_obj_remove->proc~sllist_obj_remove none~get~2 sllist_obj%Get proc~hashtbl_obj_get_int4->none~get~2 proc~hashtbl_obj_get_int8->none~get~2 proc~hashtbl_obj_get_real_dp->none~get~2 proc~hashtbl_obj_get_real_sp->none~get~2 proc~sllist_obj_put->proc~sllist_obj_put proc~sllist_obj_get_int4 sllist_obj%sllist_obj_Get_int4 none~get~2->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8 sllist_obj%sllist_obj_Get_int8 none~get~2->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp sllist_obj%sllist_obj_Get_real_dp none~get~2->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp sllist_obj%sllist_obj_Get_real_sp none~get~2->proc~sllist_obj_get_real_sp proc~sllist_obj_get_int4->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp->proc~sllist_obj_get_real_sp

Source Code

    impure subroutine ResPart_set_WriteHDF5(this,iter,time)
      !> Writes ResPart data to file using HDF5.
      use leapUtils, only: stringtool_obj, sysutils_obj
      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(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)

      ! 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 centroid 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("/","u",   buffr )
        call hdf5%Write("/","v",   buffr )
        call hdf5%Write("/","w",   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 (ResPart_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("/","u",   buffr )
          buffr =particles(1:this%count_)%v(2) ; call hdf5%Write("/","v",   buffr )
          buffr =particles(1:this%count_)%v(3) ; call hdf5%Write("/","w",   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()

      ! Write the associated IB data
      call this%ib%WriteHDF5(iter,time)

      return
    end subroutine ResPart_set_WriteHDF5