32 #include "sll_working_precision.h"
89 sll_int32,
intent(in) :: rank
90 character(len=*),
intent(in) :: file_name
91 character(len=*),
intent(in) :: mesh_name
92 sll_int32,
intent(in) :: nnodes_x1
93 sll_int32,
intent(in) :: nnodes_x2
94 sll_int32,
intent(out) :: file_id
95 sll_int32,
intent(out) :: error
100 nnodes_x1, nnodes_x2)
116 sll_int32,
intent(in) :: rank
117 character(len=*),
intent(in) :: file_name
118 character(len=*),
intent(in) :: mesh_name
119 sll_int32,
intent(in) :: nnodes_x1
120 sll_int32,
intent(in) :: nnodes_x2
121 sll_int32,
intent(in) :: nnodes_x3
122 sll_int32,
intent(out) :: file_id
123 sll_int32,
intent(out) :: error
128 nnodes_x1, nnodes_x2, nnodes_x3)
135 array, array_name, error, &
138 character(len=*),
intent(in) :: mesh_name
139 integer(i64),
intent(in) :: global_dims(2)
140 integer(i64),
intent(in) :: offset(2)
141 sll_real64,
intent(in) :: array(:, :)
142 character(len=*),
intent(in) :: array_name
143 sll_int32,
intent(out) :: error
144 sll_int32,
intent(in),
optional :: xmffile_id
145 character(len=4),
intent(in),
optional :: center
148 sll_int32 :: npoints_x1
149 sll_int32 :: npoints_x2
158 array,
"/"//trim(array_name), error)
164 if (
present(xmffile_id) .and.
present(center) .and. prank == 0)
then
165 npoints_x1 = int(global_dims(1), 4)
166 npoints_x2 = int(global_dims(2), 4)
172 trim(mesh_name)//
"-"//trim(array_name)//
".h5:/"//trim(array_name), &
179 trim(mesh_name)//
"-"//trim(array_name)//
".bin", &
180 npoints_x1, npoints_x2,
'Binary', center)
188 array, array_name, error, xmffile_id, center)
190 character(len=*),
intent(in) :: mesh_name
191 integer(i64),
intent(in) :: global_dims(3)
192 integer(i64),
intent(in) :: offset(3)
193 sll_real64,
intent(in) :: array(:, :, :)
194 character(len=*),
intent(in) :: array_name
195 sll_int32,
intent(out) :: error
196 sll_int32,
intent(in),
optional :: xmffile_id
197 character(len=4),
intent(in),
optional :: center
200 sll_int32 :: npoints_x1
201 sll_int32 :: npoints_x2
202 sll_int32 :: npoints_x3
211 "/"//trim(array_name), error)
216 if (
present(xmffile_id) .and.
present(center) .and. prank == 0)
then
217 npoints_x1 = int(global_dims(1), 4)
218 npoints_x2 = int(global_dims(2), 4)
219 npoints_x3 = int(global_dims(3), 4)
225 trim(mesh_name)//
"-"//trim(array_name)//
".h5:/"//trim(array_name), &
233 trim(mesh_name)//
"-"//trim(array_name)//
".bin", &
234 npoints_x1, npoints_x2, npoints_x3,
'Binary', center)
242 sll_int32,
intent(in) :: file_id
243 sll_int32,
intent(out) :: error
Collectively write distributed nD array into HDF5 file.
Write and array in an xmf file.
write a data attribute in the xml file
write grid description in the xml file
integer(kind=i32) function, public sll_f_get_collective_rank(col)
Determines the rank of the calling process in the communicator.
type(sll_t_collective_t), pointer, public sll_v_world_collective
Control of subset assignment. Processes with the same color are in the same new communicator.
Parallel version of sll_hdf5_io.
subroutine, public sll_s_hdf5_par_file_close(handle, error)
Close existing HDF5 file.
subroutine, public sll_s_hdf5_par_file_create(filename, comm, handle, error)
Create new HDF5 file.
Implements the functions to write xdmf file plotable by VisIt.
subroutine sll_xdmf_array_2d_parallel(mesh_name, global_dims, offset, array, array_name, error, xmffile_id, center)
Write 2d array in parallel hdf5 file and the matching line in XDMF file.
subroutine sll_xdmf_close_parallel(file_id, error)
Close the XML file and finish to write last lines.
subroutine sll_xdmf_array_3d_parallel(mesh_name, global_dims, offset, array, array_name, error, xmffile_id, center)
Write 3d array in binary or hdf5 file and the matching line in XDMF file.
subroutine sll_xdmf_open_3d_parallel(rank, file_name, mesh_name, nnodes_x1, nnodes_x2, nnodes_x3, file_id, error)
Open a XDMF format file for a 3d plot.
subroutine sll_xdmf_open_2d_parallel(rank, file_name, mesh_name, nnodes_x1, nnodes_x2, file_id, error)
Open a XDMF format file for a 2d plot.
Implements the functions to write xml file to store light data.
subroutine, public sll_s_xml_file_close(file_id, error)
Close the XML file and finish to write last lines. You give the file unit number.
subroutine, public sll_s_xml_file_create(filename, file_id, error)
Create the XML file and begin to write first lines. You get the file unit number.
Opaque object around HDF5 file id.