WriteHDF5_MPI_Vector2D_t Subroutine

public subroutine WriteHDF5_MPI_Vector2D_t(this, fileId, group, elemoffset, nglobalelem)

Arguments

TypeIntentOptionalAttributesName
class(Vector2D_t), intent(in) :: this
integer(kind=HID_T), intent(in) :: fileId
character, intent(in) :: group
integer, intent(in) :: elemoffset
integer, intent(in) :: nglobalelem

Calls

proc~~writehdf5_mpi_vector2d_t~~CallsGraph proc~writehdf5_mpi_vector2d_t WriteHDF5_MPI_Vector2D_t proc~creategroup_hdf5 CreateGroup_HDF5 proc~writehdf5_mpi_vector2d_t->proc~creategroup_hdf5 interface~writearray_hdf5 WriteArray_HDF5 proc~writehdf5_mpi_vector2d_t->interface~writearray_hdf5 h5gcreate_f h5gcreate_f proc~creategroup_hdf5->h5gcreate_f h5lexists_f h5lexists_f proc~creategroup_hdf5->h5lexists_f h5gclose_f h5gclose_f proc~creategroup_hdf5->h5gclose_f proc~writearray_hdf5_real_r2_serial WriteArray_HDF5_real_r2_serial interface~writearray_hdf5->proc~writearray_hdf5_real_r2_serial proc~writearray_hdf5_real_r4_serial WriteArray_HDF5_real_r4_serial interface~writearray_hdf5->proc~writearray_hdf5_real_r4_serial proc~writearray_hdf5_int32_r1_serial WriteArray_HDF5_int32_r1_serial interface~writearray_hdf5->proc~writearray_hdf5_int32_r1_serial proc~writearray_hdf5_int32_r3_serial WriteArray_HDF5_int32_r3_serial interface~writearray_hdf5->proc~writearray_hdf5_int32_r3_serial proc~writearray_hdf5_real_r3_parallel WriteArray_HDF5_real_r3_parallel interface~writearray_hdf5->proc~writearray_hdf5_real_r3_parallel proc~writearray_hdf5_real_r1_serial WriteArray_HDF5_real_r1_serial interface~writearray_hdf5->proc~writearray_hdf5_real_r1_serial proc~writearray_hdf5_real_r5_serial WriteArray_HDF5_real_r5_serial interface~writearray_hdf5->proc~writearray_hdf5_real_r5_serial proc~writearray_hdf5_real_r3_serial WriteArray_HDF5_real_r3_serial interface~writearray_hdf5->proc~writearray_hdf5_real_r3_serial proc~writearray_hdf5_int32_r2_serial WriteArray_HDF5_int32_r2_serial interface~writearray_hdf5->proc~writearray_hdf5_int32_r2_serial proc~writearray_hdf5_int32_r4_serial WriteArray_HDF5_int32_r4_serial interface~writearray_hdf5->proc~writearray_hdf5_int32_r4_serial proc~writearray_hdf5_real_r4_parallel WriteArray_HDF5_real_r4_parallel interface~writearray_hdf5->proc~writearray_hdf5_real_r4_parallel h5sclose_f h5sclose_f proc~writearray_hdf5_real_r2_serial->h5sclose_f h5dwrite_f h5dwrite_f proc~writearray_hdf5_real_r2_serial->h5dwrite_f h5screate_simple_f h5screate_simple_f proc~writearray_hdf5_real_r2_serial->h5screate_simple_f h5dcreate_f h5dcreate_f proc~writearray_hdf5_real_r2_serial->h5dcreate_f h5dclose_f h5dclose_f proc~writearray_hdf5_real_r2_serial->h5dclose_f proc~writearray_hdf5_real_r4_serial->h5sclose_f proc~writearray_hdf5_real_r4_serial->h5dwrite_f proc~writearray_hdf5_real_r4_serial->h5screate_simple_f proc~writearray_hdf5_real_r4_serial->h5dcreate_f proc~writearray_hdf5_real_r4_serial->h5dclose_f proc~writearray_hdf5_int32_r1_serial->h5sclose_f proc~writearray_hdf5_int32_r1_serial->h5dwrite_f proc~writearray_hdf5_int32_r1_serial->h5screate_simple_f proc~writearray_hdf5_int32_r1_serial->h5dcreate_f proc~writearray_hdf5_int32_r1_serial->h5dclose_f proc~writearray_hdf5_int32_r3_serial->h5sclose_f proc~writearray_hdf5_int32_r3_serial->h5dwrite_f proc~writearray_hdf5_int32_r3_serial->h5screate_simple_f proc~writearray_hdf5_int32_r3_serial->h5dcreate_f proc~writearray_hdf5_int32_r3_serial->h5dclose_f proc~writearray_hdf5_real_r3_parallel->h5sclose_f proc~writearray_hdf5_real_r3_parallel->h5dwrite_f h5pset_dxpl_mpio_f h5pset_dxpl_mpio_f proc~writearray_hdf5_real_r3_parallel->h5pset_dxpl_mpio_f proc~writearray_hdf5_real_r3_parallel->h5screate_simple_f proc~writearray_hdf5_real_r3_parallel->h5dcreate_f proc~writearray_hdf5_real_r3_parallel->h5dclose_f h5pcreate_f h5pcreate_f proc~writearray_hdf5_real_r3_parallel->h5pcreate_f h5sselect_hyperslab_f h5sselect_hyperslab_f proc~writearray_hdf5_real_r3_parallel->h5sselect_hyperslab_f h5pclose_f h5pclose_f proc~writearray_hdf5_real_r3_parallel->h5pclose_f proc~writearray_hdf5_real_r1_serial->h5sclose_f proc~writearray_hdf5_real_r1_serial->h5dwrite_f proc~writearray_hdf5_real_r1_serial->h5screate_simple_f proc~writearray_hdf5_real_r1_serial->h5dcreate_f proc~writearray_hdf5_real_r1_serial->h5dclose_f proc~writearray_hdf5_real_r5_serial->h5sclose_f proc~writearray_hdf5_real_r5_serial->h5dwrite_f proc~writearray_hdf5_real_r5_serial->h5screate_simple_f proc~writearray_hdf5_real_r5_serial->h5dcreate_f proc~writearray_hdf5_real_r5_serial->h5dclose_f proc~writearray_hdf5_real_r3_serial->h5sclose_f proc~writearray_hdf5_real_r3_serial->h5dwrite_f proc~writearray_hdf5_real_r3_serial->h5screate_simple_f proc~writearray_hdf5_real_r3_serial->h5dcreate_f proc~writearray_hdf5_real_r3_serial->h5dclose_f proc~writearray_hdf5_int32_r2_serial->h5sclose_f proc~writearray_hdf5_int32_r2_serial->h5dwrite_f proc~writearray_hdf5_int32_r2_serial->h5screate_simple_f proc~writearray_hdf5_int32_r2_serial->h5dcreate_f proc~writearray_hdf5_int32_r2_serial->h5dclose_f proc~writearray_hdf5_int32_r4_serial->h5sclose_f proc~writearray_hdf5_int32_r4_serial->h5dwrite_f proc~writearray_hdf5_int32_r4_serial->h5screate_simple_f proc~writearray_hdf5_int32_r4_serial->h5dcreate_f proc~writearray_hdf5_int32_r4_serial->h5dclose_f proc~writearray_hdf5_real_r4_parallel->h5sclose_f proc~writearray_hdf5_real_r4_parallel->h5dwrite_f proc~writearray_hdf5_real_r4_parallel->h5pset_dxpl_mpio_f proc~writearray_hdf5_real_r4_parallel->h5screate_simple_f proc~writearray_hdf5_real_r4_parallel->h5dcreate_f proc~writearray_hdf5_real_r4_parallel->h5dclose_f proc~writearray_hdf5_real_r4_parallel->h5pcreate_f proc~writearray_hdf5_real_r4_parallel->h5sselect_hyperslab_f proc~writearray_hdf5_real_r4_parallel->h5pclose_f

Contents


Source Code

  subroutine WriteHDF5_MPI_Vector2D_t(this,fileId,group,elemoffset,nglobalelem)
    implicit none
    class(Vector2D_t),intent(in) :: this
    character(*),intent(in) :: group
    integer(HID_T),intent(in) :: fileId
    integer,intent(in) :: elemoffset
    integer,intent(in) :: nglobalelem
    ! Local
    integer(HID_T) :: offset(1:3)
    integer(HID_T) :: globalDims(1:3)
    integer :: ivar,idir
    character(4) :: dimvar

    offset(1:3) = (/0,0,elemoffset/)
    globalDims(1:3) = (/this%interp%N+1, &
                        this%interp%N+1, &
                        nglobalelem/)

    call CreateGroup_HDF5(fileId,trim(group))

    do idir = 1,2
      write(dimvar,'(I1)') idir
      dimvar = "dim"//trim(dimvar)
      do ivar = 1,this%nVar
        !call this%meta(ivar)%WriteHDF5(group,ivar,fileId)
        call WriteArray_HDF5(fileId, &
                             trim(group)//"/"//trim(this%meta(ivar)%name)//"_"//dimvar, &
                             this%interior(:,:,:,ivar,idir),offset,globalDims)
      enddo
    enddo

  endsubroutine WriteHDF5_MPI_Vector2D_t