particle_set_WriteH5HUT Subroutine

private impure subroutine particle_set_WriteH5HUT(this, iter, time)

Writes particle data to file in parallel using H5HUT.

Type Bound

particle_set

Arguments

Type IntentOptional Attributes Name
class(particle_set), intent(inout) :: this

Lagrangian array to dump

integer, intent(in) :: iter

Iteration at write

real(kind=wp), intent(in) :: time

Time at write


Calls

proc~~particle_set_writeh5hut~~CallsGraph proc~particle_set_writeh5hut particle_set%particle_set_WriteH5HUT c c proc~particle_set_writeh5hut->c fc fc proc~particle_set_writeh5hut->fc fh fh proc~particle_set_writeh5hut->fh p p proc~particle_set_writeh5hut->p particles particles proc~particle_set_writeh5hut->particles proc~h5hut_obj_final h5hut_obj%h5hut_obj_Final proc~particle_set_writeh5hut->proc~h5hut_obj_final proc~h5hut_obj_init h5hut_obj%h5hut_obj_Init proc~particle_set_writeh5hut->proc~h5hut_obj_init proc~h5hut_obj_newtimestep h5hut_obj%h5hut_obj_NewTimeStep proc~particle_set_writeh5hut->proc~h5hut_obj_newtimestep proc~lagrangian_set_updatecount lagrangian_set%lagrangian_set_UpdateCount proc~particle_set_writeh5hut->proc~lagrangian_set_updatecount tc tc proc~particle_set_writeh5hut->tc th th proc~particle_set_writeh5hut->th v v proc~particle_set_writeh5hut->v w w proc~particle_set_writeh5hut->w proc~hdf5_obj_final hdf5_obj%hdf5_obj_Final proc~h5hut_obj_final->proc~hdf5_obj_final proc~hdf5_obj_init hdf5_obj%hdf5_obj_Init proc~h5hut_obj_init->proc~hdf5_obj_init none~writeattributes hdf5_obj%WriteAttributes proc~h5hut_obj_newtimestep->none~writeattributes proc~h5hut_obj_getnsteps h5hut_obj%h5hut_obj_GetNSteps proc~h5hut_obj_newtimestep->proc~h5hut_obj_getnsteps proc~h5hut_obj_setstep h5hut_obj%h5hut_obj_SetStep proc~h5hut_obj_newtimestep->proc~h5hut_obj_setstep mpi_allgather mpi_allgather proc~lagrangian_set_updatecount->mpi_allgather 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 proc~hdf5_obj_readgroupnames hdf5_obj%hdf5_obj_ReadGroupNames proc~h5hut_obj_getnsteps->proc~hdf5_obj_readgroupnames proc~hdf5_obj_creategroup hdf5_obj%hdf5_obj_CreateGroup proc~h5hut_obj_setstep->proc~hdf5_obj_creategroup 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 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 h5gcreate_f h5gcreate_f proc~hdf5_obj_creategroup->h5gcreate_f proc~hashtbl_obj_hashstring hashtbl_obj%hashtbl_obj_HashString proc~hdf5_obj_creategroup->proc~hashtbl_obj_hashstring proc~hashtbl_obj_put hashtbl_obj%hashtbl_obj_Put proc~hdf5_obj_creategroup->proc~hashtbl_obj_put proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_creategroup->proc~hdf5_obj_closegroup proc~hdf5_obj_fixgroupname hdf5_obj%hdf5_obj_FixGroupName proc~hdf5_obj_creategroup->proc~hdf5_obj_fixgroupname h5gn_members_f h5gn_members_f proc~hdf5_obj_readgroupnames->h5gn_members_f h5iget_name_f h5iget_name_f proc~hdf5_obj_readgroupnames->h5iget_name_f h5oclose_f h5oclose_f proc~hdf5_obj_readgroupnames->h5oclose_f h5oget_info_f h5oget_info_f proc~hdf5_obj_readgroupnames->h5oget_info_f h5oopen_by_idx_f h5oopen_by_idx_f proc~hdf5_obj_readgroupnames->h5oopen_by_idx_f proc~hdf5_obj_readgroupnames->proc~hdf5_obj_fixgroupname 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_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~sllist_obj_put sllist_obj%sllist_obj_Put proc~hashtbl_obj_put->proc~sllist_obj_put proc~hdf5_obj_closegroup->proc~hashtbl_obj_hashstring proc~hdf5_obj_closegroup->proc~hdf5_obj_fixgroupname proc~hdf5_obj_closegroup->proc~hdf5_obj_getgroupobject h5gclose_f h5gclose_f proc~hdf5_obj_closegroup->h5gclose_f proc~hashtbl_obj_remove hashtbl_obj%hashtbl_obj_Remove proc~hdf5_obj_closegroup->proc~hashtbl_obj_remove proc~hdf5_obj_getgroupobject->proc~hashtbl_obj_hashstring proc~hdf5_obj_getgroupobject->proc~hdf5_obj_fixgroupname none~get~3 hashtbl_obj%Get proc~hdf5_obj_getgroupobject->none~get~3 proc~hdf5_obj_opengroup->proc~hashtbl_obj_hashstring proc~hdf5_obj_opengroup->proc~hashtbl_obj_put proc~hdf5_obj_opengroup->proc~hdf5_obj_fixgroupname proc~hdf5_obj_opengroup->proc~hdf5_obj_getgroupobject h5oopen_f h5oopen_f proc~hdf5_obj_opengroup->h5oopen_f 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_remove sllist_obj%sllist_obj_Remove proc~hashtbl_obj_remove->proc~sllist_obj_remove proc~sllist_obj_put->proc~sllist_obj_put 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_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 particle_set_WriteH5HUT(this,iter,time)
      !> Writes particle data to file in parallel using H5HUT.
      implicit none
      class(particle_set), intent(inout) :: this                               !! Lagrangian array to dump
      integer,             intent(in)    :: iter                               !! Iteration at write
      real(wp),            intent(in)    :: time                               !! Time at write
      ! Work variables
      type(h5hut_obj) :: h5                                                    !! H5hut structure
      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_))

      ! Write the particle data
      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("vx", buffr )
        call h5%Write("vy", buffr )
        call h5%Write("vz", 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 (particle_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("vx", buffr )
          buffr =particles(1:this%count_)%v(2) ; call h5%Write("vy", buffr )
          buffr =particles(1:this%count_)%v(3) ; call h5%Write("vz", 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

      return
    end subroutine particle_set_WriteH5HUT