Writes Eulerian/3D data to a hdf5 file with silo.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(silo_obj), | intent(inout) | :: | this |
A silo object |
||
| character(len=*), | intent(in) | :: | name |
Variable name |
||
| character(len=*), | intent(in) | :: | mesh_name |
Mesh name |
||
| class(*), | intent(in) | :: | array(:,:,:) |
3-D data array |
impure subroutine silo_obj_WriteScalar3D(this,name,mesh_name,array) !> Writes Eulerian/3D data to a hdf5 file with silo. implicit none class(silo_obj), intent(inout) :: this !! A silo object character(len=*), intent(in) :: mesh_name !! Mesh name character(len=*), intent(in) :: name !! Variable name class(*), intent(in) :: array(:,:,:) !! 3-D data array ! Work variables integer :: ierr,status integer :: n integer :: Ng(3) character(len=SILOstr), allocatable :: names(:) integer, allocatable :: lnames(:),types(:) Ng = [size(array,1),size(array,2),size(array,3)] ! Write using poor man's IO: sequential (serial) writes among a silo group do n=1,this%nproc_node if (n.eq.this%silo_rank) then ! Open file ierr=DBopen(this%siloname,len_trim(this%siloname),SILOdriver,DB_APPEND,this%fid_DAT) if (ierr.ne.0) call this%parallel%Stop('Error opening SILO file') ! Switch to appropriate directory ierr = DBsetdir(this%fid_DAT,trim(this%dirname),len_trim(this%dirname)) ! Write data in Single Precision select type(array) type is (real(leapDP)) ierr = DBputqv1(this%fid_DAT,name,len_trim(name), mesh_name, len_trim(mesh_name), real(array,leapSP), & Ng, 3, DB_F77NULL, 0, DB_FLOAT, DB_ZONECENT, DB_F77NULL, status) type is (real(leapSP)) ierr = DBputqv1(this%fid_DAT,name,len_trim(name), mesh_name, len_trim(mesh_name), real(array,leapSP), & Ng, 3, DB_F77NULL, 0, DB_FLOAT, DB_ZONECENT, DB_F77NULL, status) type is (integer) ierr = DBputqv1(this%fid_DAT,name,len_trim(name), mesh_name, len_trim(mesh_name), int(array,leapI4), & Ng, 3, DB_F77NULL, 0, DB_INT, DB_ZONECENT, DB_F77NULL, status) type is (integer(leapI8)) ierr = DBputqv1(this%fid_DAT,name,len_trim(name), mesh_name, len_trim(mesh_name), int(array,leapI4), & Ng, 3, DB_F77NULL, 0, DB_INT, DB_ZONECENT, DB_F77NULL, status) end select ! Close file so that next silo_rank can open it ierr = DBclose(this%fid_DAT) if (ierr.ne.0) call this%parallel%Stop('Error closing SILO file') end if call MPI_BARRIER(this%silo_comm) end do allocate (names(this%parallel%nproc)) allocate (lnames(this%parallel%nproc)) allocate (types(this%parallel%nproc)) if (this%parallel%RankIsRoot()) then do n=1,this%parallel%nproc write(names(n), '(2A,I5.5,2A,I5.5,2A)') trim(adjustl(this%filename)),'_g_',& this%group_ids(n),SILO_EXTENSION,':',mod(n-1,this%nproc_node)+1,'/',trim(name) lnames(n) = len_trim(names(n)) types(n) = DB_QUADVAR end do ierr = DBputmvar(this%fid_VDB, name, len_trim(name), this%parallel%nproc,names,lnames,types,DB_F77NULL,ierr) end if deallocate(names) deallocate(lnames) deallocate(types) return end subroutine silo_obj_WriteScalar3D