25 #include "sll_assert.h"
26 #include "sll_working_precision.h"
67 character(len=*),
intent(in) :: filename
68 sll_int32,
intent(out) :: file_id
69 sll_int32,
intent(out) :: error
73 do 100 file_id = 20, 99
74 inquire (unit=file_id, opened=lopen)
78 open (file_id, status=
'SCRATCH', err=100)
79 close (file_id, status=
'DELETE', err=100)
89 inquire (file=filename, opened=lopen)
90 sll_assert_always(.not. lopen)
92 open (file_id, file=filename, form=
'FORMATTED', iostat=error)
95 write (file_id,
"(a)")
"<?xml version='1.0' ?>"
96 write (file_id,
"(a)")
"<!DOCTYPE Xdmf SYSTEM 'Xdmf.dtd' []>"
97 write (file_id,
"(a)")
"<Xdmf xmlns:xi=""http://www.w3.org/2003/XInclude"" Version=""2.2"">"
98 write (file_id,
"(a)")
"<Domain>"
108 sll_int32,
intent(in) :: file_id
109 sll_int32,
intent(out) :: error
111 write (file_id,
"(a)")
"</Grid>"
112 write (file_id,
"(a)")
"</Domain>"
113 write (file_id,
"(a)")
"</Xdmf>"
132 sll_int32,
intent(in) :: file_id
133 character(len=*),
intent(in) :: filename
134 character(len=*),
intent(in) :: filetype
135 sll_int32,
intent(in) :: nnodes_x1
137 sll_assert_always(filetype ==
'HDF' .or. filetype ==
'Binary')
138 write (file_id,
"(a,i10,a)")
"<DataItem Dimensions='", nnodes_x1, &
139 "' NumberType='Float' Precision='8' Format='"//trim(filetype)//
"'>"
140 write (file_id,
"(a)") trim(filename)
141 write (file_id,
"(a)")
"</DataItem>"
161 sll_int32,
intent(in) :: file_id
162 character(len=*),
intent(in) :: filename
163 character(len=*),
intent(in) :: filetype
164 sll_int32,
intent(in) :: nnodes_x1
165 sll_int32,
intent(in) :: nnodes_x2
167 sll_assert_always(filetype ==
'HDF' .or. filetype ==
'Binary')
168 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nnodes_x2, nnodes_x1, &
169 "' NumberType='Float' Precision='8' Format='"//trim(filetype)//
"'>"
170 write (file_id,
"(a)") trim(filename)
171 write (file_id,
"(a)")
"</DataItem>"
182 sll_int32,
intent(in) :: file_id
183 character(len=*),
intent(in) :: filename
184 character(len=*),
intent(in) :: filetype
185 sll_int32,
intent(in) :: nnodes_x1
186 sll_int32,
intent(in) :: nnodes_x2
187 sll_int32,
intent(in) :: nnodes_x3
189 sll_assert_always(filetype ==
'HDF' .or. filetype ==
'Binary')
190 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nnodes_x3, &
191 nnodes_x2, nnodes_x1, &
192 "' NumberType='Float' Precision='8' Format='"//trim(filetype)//
"'>"
193 write (file_id,
"(a)") trim(filename)
194 write (file_id,
"(a)")
"</DataItem>"
216 sll_int32,
intent(in) :: file_id
217 character(len=*),
intent(in) :: filename
218 character(len=*),
intent(in) :: fieldname
219 character(len=*),
intent(in) :: center
220 character(len=*),
intent(in) :: filetype
221 sll_int32,
intent(in) :: npoints_1
223 write (file_id,
"(a)") &
224 "<Attribute Name='"//trim(fieldname)//
"' AttributeType='Scalar' Center='"//center//
"'>"
226 write (file_id,
"(a)")
"</Attribute>"
250 sll_int32,
intent(in) :: file_id
251 character(len=*),
intent(in) :: filename
252 character(len=*),
intent(in) :: fieldname
253 character(len=*),
intent(in) :: center
254 character(len=*),
intent(in) :: filetype
255 sll_int32,
intent(in) :: npoints_1
256 sll_int32,
intent(in) :: npoints_2
258 write (file_id,
"(a)") &
259 "<Attribute Name='"//trim(fieldname)//
"' AttributeType='Scalar' Center='"//center//
"'>"
261 write (file_id,
"(a)")
"</Attribute>"
288 sll_int32,
intent(in) :: file_id
289 character(len=*),
intent(in) :: filename
290 character(len=*),
intent(in) :: fieldname
291 character(len=*),
intent(in) :: center
292 character(len=*),
intent(in) :: filetype
293 sll_int32,
intent(in) :: npoints_1
294 sll_int32,
intent(in) :: npoints_2
295 sll_int32,
intent(in) :: npoints_3
297 write (file_id,
"(a)") &
298 "<Attribute Name='"//trim(fieldname)//
"' AttributeType='Scalar' Center='"//center//
"'>"
305 write (file_id,
"(a)")
"</Attribute>"
324 sll_int32,
intent(in) :: file_id
325 character(len=*),
intent(in) :: filename
326 sll_int32,
intent(in) :: nnodes_x1
327 sll_int32,
intent(in) :: nnodes_x2
331 trim(filename)//
"-x1.bin", nnodes_x1, &
332 trim(filename)//
"-x2.bin", nnodes_x2,
"x1",
"x2",
'Uniform')
335 trim(filename)//
"-x1.h5", nnodes_x1, &
336 trim(filename)//
"-x2.h5", nnodes_x2,
"x1",
"x2",
'Uniform')
363 sll_int32,
intent(in) :: file_id
364 character(len=*),
intent(in) :: x1filename
365 character(len=*),
intent(in) :: x2filename
366 sll_int32,
intent(in) :: nnodes_x1
367 sll_int32,
intent(in) :: nnodes_x2
368 character(len=*),
intent(in) :: x1dsetname
369 character(len=*),
intent(in) :: x2dsetname
370 character(len=*),
intent(in) :: gridtype
372 write (file_id,
"(a)")
"<Grid Name='mesh' GridType='"//gridtype//
"'>"
374 "(a,2i5,a)")
"<Topology TopologyType='2DSMesh' NumberOfElements='", &
375 nnodes_x2, nnodes_x1,
"'/>"
376 write (file_id,
"(a)")
"<Geometry GeometryType='X_Y'>"
381 nnodes_x1, nnodes_x2,
'Binary')
384 nnodes_x1, nnodes_x2,
'Binary')
387 trim(x1filename)//
":/"//x1dsetname, &
388 nnodes_x1, nnodes_x2,
'HDF')
390 trim(x2filename)//
":/"//x2dsetname, &
391 nnodes_x1, nnodes_x2,
'HDF')
394 write (file_id,
"(a)")
"</Geometry>"
401 nnodes_x1, nnodes_x2, nnodes_x3)
403 sll_int32,
intent(in) :: file_id
404 character(len=*),
intent(in) :: filename
405 sll_int32,
intent(in) :: nnodes_x1
406 sll_int32,
intent(in) :: nnodes_x2
407 sll_int32,
intent(in) :: nnodes_x3
409 write (file_id,
"(a)")
"<Grid Name='mesh' GridType='Uniform'>"
410 write (file_id,
"(a,3i5,a)")
"<Topology TopologyType='3DSMesh' NumberOfElements='", &
411 nnodes_x3, nnodes_x2, nnodes_x1,
"'/>"
412 write (file_id,
"(a)")
"<Geometry GeometryType='X_Y_Z'>"
417 nnodes_x1, nnodes_x2, nnodes_x3,
'Binary')
419 nnodes_x1, nnodes_x2, nnodes_x3,
'Binary')
421 nnodes_x1, nnodes_x2, nnodes_x3,
'Binary')
425 nnodes_x1, nnodes_x2, nnodes_x3,
'HDF')
427 nnodes_x1, nnodes_x2, nnodes_x3,
'HDF')
429 nnodes_x1, nnodes_x2, nnodes_x3,
'HDF')
433 write (file_id,
"(a)")
"</Geometry>"
441 x1filename, nnodes_x1, &
442 x2filename, nnodes_x2, &
443 x3filename, nnodes_x3, &
444 x1dsetname, x2dsetname, x3dsetname, &
447 sll_int32,
intent(in) :: file_id
448 character(len=*),
intent(in) :: x1filename
449 character(len=*),
intent(in) :: x2filename
450 character(len=*),
intent(in) :: x3filename
451 character(len=*),
intent(in) :: x1dsetname
452 character(len=*),
intent(in) :: x2dsetname
453 character(len=*),
intent(in) :: x3dsetname
454 sll_int32,
intent(in) :: nnodes_x1
455 sll_int32,
intent(in) :: nnodes_x2
456 sll_int32,
intent(in) :: nnodes_x3
457 character(len=*),
intent(in) :: gridtype
459 write (file_id,
"(a)")
"<Grid Name='mesh' GridType='"//gridtype//
"'>"
460 write (file_id,
"(a,3i5,a)") &
461 "<Topology TopologyType='3DSMesh' NumberOfElements='", &
462 nnodes_x3, nnodes_x2, nnodes_x1,
"'/>"
463 write (file_id,
"(a)")
"<Geometry GeometryType='X_Y_Z'>"
468 nnodes_x1, nnodes_x2, nnodes_x3,
'Binary')
470 nnodes_x1, nnodes_x2, nnodes_x3,
'Binary')
472 nnodes_x1, nnodes_x2, nnodes_x3,
'Binary')
476 nnodes_x1, nnodes_x2, nnodes_x3,
'HDF')
478 nnodes_x1, nnodes_x2, nnodes_x3,
'HDF')
480 nnodes_x1, nnodes_x2, nnodes_x3,
'HDF')
484 write (file_id,
"(a)")
"</Geometry>"
write a data attribute in the xml file
write grid description in the xml file
write a data item in the xml file
Implements the functions to write xml file to store light data.
subroutine sll_xml_grid_geometry_3d_high_level(file_id, filename, nnodes_x1, nnodes_x2, nnodes_x3)
Write the description of a 3D structured curvilinear grid mesh with its nodes coordinates contains in...
subroutine sll_xml_field_2d(file_id, fieldname, filename, npoints_1, npoints_2, filetype, center)
Write the description of a scalar field on a 2D mesh.
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 sll_xml_field_3d(file_id, fieldname, filename, npoints_1, npoints_2, npoints_3, filetype, center)
Write the description of a scalar field on a 3D mesh.
subroutine sll_xml_dataitem_2d(file_id, filename, nnodes_x1, nnodes_x2, filetype)
Write the description of a scalar field on a 3D mesh. Write the description of a scalar field on a 2D...
subroutine sll_xml_grid_geometry_3d_low_level(file_id, x1filename, nnodes_x1, x2filename, nnodes_x2, x3filename, nnodes_x3, x1dsetname, x2dsetname, x3dsetname, gridtype)
Write the description of a 3D structured curvilinear grid mesh with its nodes coordinates contains in...
subroutine sll_xml_grid_geometry_2d_high_level(file_id, filename, nnodes_x1, nnodes_x2)
Write the description of a 2D strutured grid mesh with its nodes coordinates contains in filename-x1 ...
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.
subroutine sll_xml_grid_geometry_2d_low_level(file_id, x1filename, nnodes_x1, x2filename, nnodes_x2, x1dsetname, x2dsetname, gridtype)
Write the description of a 2D strutured grid mesh with its nodes coordinates contains in filename-x1 ...
subroutine sll_xml_dataitem_3d(file_id, filename, nnodes_x1, nnodes_x2, nnodes_x3, filetype)
Write the description of a scalar field on a 3D mesh.
subroutine sll_xml_dataitem_1d(file_id, filename, nnodes_x1, filetype)
Write the description of a scalar field on a 1D mesh.
subroutine sll_xml_field_1d(file_id, fieldname, filename, npoints_1, filetype, center)
Write the description of a scalar field on a 1D mesh.