28 #include "sll_assert.h"
29 #include "sll_working_precision.h"
89 character(len=*),
intent(in) :: prefix
90 sll_real64,
intent(in) :: x1(:, :)
91 sll_real64,
intent(in) :: x2(:, :)
92 sll_int32,
intent(out) :: xmf
95 character(len=4) :: crank
100 call mpi_comm_rank(mpi_comm_world, prank, error)
116 array, array_name, error, &
119 sll_int32,
intent(in) :: xmf
120 character(len=*),
intent(in) :: prefix
121 sll_real64,
intent(in) :: array(:, :)
122 character(len=*),
intent(in) :: array_name
128 sll_int32,
intent(out) :: error
129 character(len=4),
optional :: center
132 character(len=4) :: crank
134 character(len=4) :: cproc
135 sll_int32,
parameter :: tag = 1111
136 sll_int32 :: statut(mpi_status_size)
138 call mpi_comm_rank(mpi_comm_world, prank, error)
139 call mpi_comm_size(mpi_comm_world, psize, error)
148 npts_x1 =
size(array, 1)
149 npts_x2 =
size(array, 2)
153 "<Grid Name=""AllDomain"" GridType=""Collection"" CollectionType=""Spatial"">"
154 do iproc = 0, psize - 1
157 call mpi_recv(npts_x1, 1, mpi_integer, iproc, tag, mpi_comm_world, statut, error)
158 call mpi_recv(npts_x2, 1, mpi_integer, iproc, tag, mpi_comm_world, statut, error)
161 write (xmf,
'(a)')
"<Grid Name=""SubDomain"" GridType=""Uniform"">"
162 write (xmf,
'(a,2i6,a)') &
163 "<Topology TopologyType=""2DSMesh"" NumberOfElements='", npts_x2, npts_x1,
"'/>"
164 write (xmf,
'(a)')
"<Geometry GeometryType=""X_Y"">"
165 write (xmf,
'(a,2i6,a)')
"<DataItem Dimensions='", npts_x2, npts_x1, &
166 "' NumberType=""Float"" Precision=""4"" Format=""HDF"">"
167 write (xmf,
'(a)') trim(prefix)//
"-mesh-"//cproc//
".h5:/x1"
168 write (xmf,
'(a)')
"</DataItem>"
169 write (xmf,
'(a,2i6,a)')
"<DataItem Dimensions='", npts_x2, npts_x1, &
170 "' NumberType=""Float"" Precision=""4"" Format=""HDF"">"
171 write (xmf,
'(a)') trim(prefix)//
"-mesh-"//cproc//
".h5:/x2"
172 write (xmf,
'(a)')
"</DataItem>"
173 write (xmf,
'(a)')
"</Geometry>"
175 "<Attribute Name='"//trim(array_name)//
"' AttributeType=""Scalar"" Center=""Node"">"
176 write (xmf,
'(a,2i6,a)')
"<DataItem Dimensions='", npts_x2, npts_x1, &
177 "' NumberType=""Float"" Precision=""8"" Format=""HDF"">"
178 write (xmf,
'(a)') trim(prefix)//
"-"//cproc//
".h5:/"//trim(array_name)
179 write (xmf,
'(a)')
"</DataItem>"
180 write (xmf,
'(a)')
"</Attribute>"
181 write (xmf,
'(a)')
"</Grid>"
184 call mpi_send(npts_x1, 1, mpi_integer, 0, tag, mpi_comm_world, error)
185 call mpi_send(npts_x2, 1, mpi_integer, 0, tag, mpi_comm_world, error)
189 sll_assert_always(
present(center))
195 sll_int32,
intent(in) :: file_id
196 sll_int32,
intent(out) :: error
198 logical :: i_opened, i_exist
199 character(len=255) :: i_name
201 call mpi_comm_rank(mpi_comm_world, prank, error)
203 INQUIRE (file_id, opened=i_opened, name=i_name, exist=i_exist)
Write nD array of double precision floats or integers into HDF5 file.
Create a new xmdf file to plot parallel array using hdf5 serial library.
Write and array in an xmf file.
Implements the functions to write hdf5 file to store heavy data.
subroutine, public sll_s_hdf5_ser_file_create(filename, handle, error)
Create new HDF5 file.
subroutine, public sll_s_hdf5_ser_file_close(handle, error)
Close existing HDF5 file.
Some common numerical utilities.
subroutine, public sll_s_int2string(istep, cstep)
Convert an integer < 9999 to a 4 characters string.
Implements the functions to write xdmf file plotable by VisIt.
subroutine, public sll_s_xdmf_array_2d_serial_blocks(xmf, prefix, array, array_name, error, center)
Write 2d array in parallel hdf5 file and the matching line in XDMF file.
subroutine, public sll_s_xdmf_open_serial_blocks(prefix, x1, x2, xmf, error)
Open a XDMF format file for a 2d plot.
subroutine, public sll_s_xdmf_close_serial_blocks(file_id, error)
Close the XML file and finish to write last lines.
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.