26 #include "sll_assert.h"
27 #include "sll_errors.h"
28 #include "sll_memory.h"
29 #include "sll_working_precision.h"
91 sll_real64,
intent(in) :: time
92 sll_int32,
intent(in) :: file_id
95 inquire (file_id, opened=i_opened)
98 write (file_id,
"(a13,g15.3,a3)")
"<Time Value='", time,
"'/>"
100 sll_error(
"sll_xdmf_set_time",
"This xdmf file is not not open.")
108 character(len=*),
intent(in) :: file_name
109 character(len=*),
intent(in) :: mesh_name
110 sll_int32,
intent(in) :: nnodes_x1
111 sll_int32,
intent(in) :: nnodes_x2
112 sll_int32,
intent(out) :: file_id
113 sll_int32,
intent(out) :: error
117 nnodes_x1, nnodes_x2)
130 character(len=*),
intent(in) :: file_name
131 character(len=*),
intent(in) :: mesh_name
132 sll_int32,
intent(out) :: file_id
133 sll_int32,
intent(out) :: error
134 sll_int32,
intent(in) :: nnodes_x1
135 sll_int32,
intent(in) :: nnodes_x2
136 sll_int32,
intent(in) :: nnodes_x3
140 nnodes_x1, nnodes_x2, nnodes_x3)
147 character(len=*),
intent(in) :: mesh_name
148 sll_real64,
intent(in) :: array(:, :)
149 character(len=*),
intent(in) :: array_name
150 sll_int32,
intent(out) :: error
151 sll_int32 :: npoints_x1
152 sll_int32 :: npoints_x2
153 sll_int32,
intent(in),
optional :: xmffile_id
154 character(len=4),
optional :: center
162 npoints_x1 =
size(array, 1)
163 npoints_x2 =
size(array, 2)
177 if (
present(xmffile_id) .and.
present(center))
then
180 trim(mesh_name)//
"-"//trim(array_name)//
".bin", &
181 npoints_x1, npoints_x2,
'Binary', center)
186 trim(mesh_name)//
"-"//trim(array_name)//
".h5:/"//trim(array_name), &
198 character(len=*),
intent(in) :: mesh_name
199 sll_real64,
intent(in) :: array(:, :, :)
200 character(len=*),
intent(in) :: array_name
201 sll_int32,
intent(out) :: error
202 sll_int32,
intent(in),
optional :: xmffile_id
203 character(len=4),
optional :: center
204 sll_int32 :: npoints_x1
205 sll_int32 :: npoints_x2
206 sll_int32 :: npoints_x3
214 npoints_x1 =
size(array, 1)
215 npoints_x2 =
size(array, 2)
216 npoints_x3 =
size(array, 3)
230 if (
present(xmffile_id) .and.
present(center))
then
234 trim(mesh_name)//
"-"//trim(array_name)//
".bin", &
235 npoints_x1, npoints_x2, npoints_x3,
'Binary', center)
240 trim(mesh_name)//
"-"//trim(array_name)//
".h5:/"//trim(array_name), &
266 sll_real64,
intent(in) :: array(:, :)
267 character(len=*),
intent(in) :: file_name
268 character(len=*),
intent(in) :: array_name
270 sll_real64 :: eta1_min
271 sll_real64 :: eta2_min
272 sll_real64 :: delta_eta1
273 sll_real64 :: delta_eta2
277 character(len=4),
optional :: file_format
278 sll_int32,
optional :: iplot
279 sll_real64,
optional :: time
281 character(len=4) :: cplot
289 if (
present(iplot))
then
296 write (file_id,
"(a)")
"<Grid Name='mesh' GridType='Uniform'>"
297 if (
present(time))
then
298 write (file_id,
"(a,f12.6,a)")
"<Time Value='", time,
"'/>"
300 write (file_id,
"(a,2i5,a)")
"<Topology TopologyType='2DCoRectMesh' NumberOfElements='", &
302 write (file_id,
"(a)")
"<Geometry GeometryType='ORIGIN_DXDY'>"
303 write (file_id,
"(a)")
"<DataItem Dimensions='2' NumberType='Float' Format='XML'>"
304 write (file_id,
"(2f12.5)") eta1_min, eta2_min
305 write (file_id,
"(a)")
"</DataItem>"
306 write (file_id,
"(a)")
"<DataItem Dimensions='2' NumberType='Float' Format='XML'>"
307 write (file_id,
"(2f12.5)") delta_eta1, delta_eta2
308 write (file_id,
"(a)")
"</DataItem>"
309 write (file_id,
"(a)")
"</Geometry>"
310 write (file_id,
"(a)")
"<Attribute Name='"//array_name//
"' AttributeType='Scalar' Center='Node'>"
312 if (
present(file_format))
then
313 if (file_format ==
"HDF5")
then
314 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nx2, nx1, &
315 "' NumberType='Float' Precision='8' Format='HDF'>"
316 if (
present(iplot))
then
317 write (file_id,
"(a)") array_name//cplot//
".h5:/node_values"
319 write (file_id,
"(a)") array_name//
".h5:/node_values"
322 if (
present(iplot))
then
332 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nx2, nx1, &
333 "' NumberType='Float' Precision='8' Format='XML'>"
336 write (file_id,
"(a)")
"</DataItem>"
337 write (file_id,
"(a)")
"</Attribute>"
357 sll_real64,
intent(in) :: array(:, :, :)
358 character(len=*),
intent(in) :: file_name
359 character(len=*),
intent(in) :: array_name
361 sll_real64 :: eta1_min
362 sll_real64 :: eta2_min
363 sll_real64 :: eta3_min
364 sll_real64 :: delta_eta1
365 sll_real64 :: delta_eta2
366 sll_real64 :: delta_eta3
371 character(len=4),
optional :: file_format
372 sll_int32,
optional :: iplot
374 character(len=4) :: cplot
375 logical :: hdf5_file_format
384 if (
present(iplot))
then
391 write (file_id,
"(a)")
"<Grid Name='mesh' GridType='Uniform'>"
392 write (file_id,
"(a,3i5,a)")
"<Topology TopologyType='3DCoRectMesh' NumberOfElements='", &
394 write (file_id,
"(a)")
"<Geometry GeometryType='ORIGIN_DXDYDZ'>"
395 write (file_id,
"(a)")
"<DataItem Dimensions='3' NumberType='Float' Format='XML'>"
396 write (file_id,
"(3f12.5)") eta1_min, eta2_min, eta3_min
397 write (file_id,
"(a)")
"</DataItem>"
398 write (file_id,
"(a)")
"<DataItem Dimensions='3' NumberType='Float' Format='XML'>"
399 write (file_id,
"(3f12.5)") delta_eta1, delta_eta2, delta_eta3
400 write (file_id,
"(a)")
"</DataItem>"
401 write (file_id,
"(a)")
"</Geometry>"
402 write (file_id,
"(a)")
"<Attribute Name='"//array_name//
"' AttributeType='Scalar' Center='Node'>"
404 hdf5_file_format = .false.
405 if (
present(file_format))
then
406 if (file_format ==
"HDF5") hdf5_file_format = .true.
408 if (hdf5_file_format)
then
409 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
410 "' NumberType='Float' Precision='8' Format='HDF'>"
411 if (
present(iplot))
then
412 write (file_id,
"(a)") array_name//cplot//
".h5:/node_values"
414 write (file_id,
"(a)") array_name//
".h5:/node_values"
417 if (
present(iplot))
then
426 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
427 "' NumberType='Float' Precision='8' Format='XML'>"
430 write (file_id,
"(a)")
"</DataItem>"
431 write (file_id,
"(a)")
"</Attribute>"
448 sll_real64,
intent(in) :: array(:, :)
449 sll_real64,
intent(in) :: eta1(:)
450 sll_real64,
intent(in) :: eta2(:)
451 character(len=*),
intent(in) :: file_name
452 character(len=*),
intent(in) :: array_name
457 character(len=4),
optional :: file_format
458 sll_int32,
optional :: iplot
459 sll_real64,
optional :: time
461 character(len=4) :: cplot
462 logical :: hdf5_file_format
470 sll_assert_always(nx1 ==
size(eta1))
471 sll_assert_always(nx2 ==
size(eta2))
473 if (
present(iplot))
then
479 write (file_id,
"(a)")
"<Grid Name='mesh' GridType='Uniform'>"
480 if (
present(time))
then
481 write (file_id,
"(a,f12.6,a)")
"<Time Value='", time,
"'/>"
483 write (file_id,
"(a,2i5,a)")
"<Topology TopologyType='2DRectMesh' NumberOfElements='", &
485 write (file_id,
"(a)")
"<Geometry GeometryType='VXVY'>"
486 write (file_id,
"(a,i5,a)")
"<DataItem Dimensions='", nx1, &
487 "' NumberType='Float' Format='XML'>"
488 write (file_id, *) (eta1(i), i=1, nx1)
489 write (file_id,
"(a)")
"</DataItem>"
490 write (file_id,
"(a,i5,a)")
"<DataItem Dimensions='", nx2, &
491 "' NumberType='Float' Format='XML'>"
492 write (file_id, *) (eta2(j), j=1, nx2)
493 write (file_id,
"(a)")
"</DataItem>"
494 write (file_id,
"(a)")
"</Geometry>"
495 write (file_id,
"(a)")
"<Attribute Name='"//array_name//
"' AttributeType='Scalar' Center='Node'>"
496 hdf5_file_format = .false.
497 if (
present(file_format))
then
498 if (file_format ==
"HDF5") hdf5_file_format = .true.
500 if (hdf5_file_format)
then
501 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nx2, nx1, &
502 "' NumberType='Float' Precision='8' Format='HDF'>"
503 if (
present(iplot))
then
504 write (file_id,
"(a)") array_name//cplot//
".h5:/node_values"
506 write (file_id,
"(a)") array_name//
".h5:/node_values"
509 if (
present(iplot))
then
518 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nx2, nx1, &
519 "' NumberType='Float' Precision='4' Format='XML'>"
522 write (file_id,
"(a)")
"</DataItem>"
523 write (file_id,
"(a)")
"</Attribute>"
540 sll_real64,
intent(in) :: array(:, :, :)
541 sll_real64,
intent(in) :: eta1(:)
542 sll_real64,
intent(in) :: eta2(:)
543 sll_real64,
intent(in) :: eta3(:)
544 character(len=*),
intent(in) :: file_name
545 character(len=*),
intent(in) :: array_name
551 character(len=4),
optional :: file_format
552 sll_int32,
optional :: iplot
557 character(len=4) :: cplot
558 logical :: hdf5_file_format
564 sll_assert_always(nx1 ==
size(eta1))
565 sll_assert_always(nx2 ==
size(eta2))
566 sll_assert_always(nx3 ==
size(eta3))
568 if (
present(iplot))
then
575 write (file_id,
"(a)")
"<Grid Name='mesh' GridType='Uniform'>"
576 write (file_id,
"(a,3i5,a)")
"<Topology TopologyType='3DRectMesh' NumberOfElements='", &
578 write (file_id,
"(a)")
"<Geometry GeometryType='VXVYVZ'>"
579 write (file_id,
"(a,i5,a)")
"<DataItem Dimensions='", nx1, &
580 "' NumberType='Float' Format='XML'>"
581 write (file_id, *) (eta1(i), i=1, nx1)
582 write (file_id,
"(a)")
"</DataItem>"
583 write (file_id,
"(a,i5,a)")
"<DataItem Dimensions='", nx2, &
584 "' NumberType='Float' Format='XML'>"
585 write (file_id, *) (eta2(j), j=1, nx2)
586 write (file_id,
"(a)")
"</DataItem>"
587 write (file_id,
"(a,i5,a)")
"<DataItem Dimensions='", nx3, &
588 "' NumberType='Float' Format='XML'>"
589 write (file_id, *) (eta3(k), k=1, nx3)
590 write (file_id,
"(a)")
"</DataItem>"
591 write (file_id,
"(a)")
"</Geometry>"
592 write (file_id,
"(a)")
"<Attribute Name='"//array_name//
"' AttributeType='Scalar' Center='Node'>"
593 hdf5_file_format = .false.
594 if (
present(file_format))
then
595 if (file_format ==
"HDF5") hdf5_file_format = .true.
597 if (hdf5_file_format)
then
598 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
599 "' NumberType='Float' Precision='8' Format='HDF'>"
600 if (
present(iplot))
then
601 write (file_id,
"(a)") array_name//cplot//
".h5:/node_values"
603 write (file_id,
"(a)") array_name//
".h5:/node_values"
606 if (
present(iplot))
then
615 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
616 "' NumberType='Float' Precision='4' Format='XML'>"
619 write (file_id,
"(a)")
"</DataItem>"
620 write (file_id,
"(a)")
"</Attribute>"
636 sll_real64,
intent(in) :: array(:, :)
637 sll_real64,
intent(in) :: eta1(:, :)
638 sll_real64,
intent(in) :: eta2(:, :)
639 character(len=*),
intent(in) :: file_name
640 character(len=*),
intent(in) :: array_name
645 character(len=4),
optional :: file_format
646 sll_int32,
optional :: iplot
647 character(len=4) :: cplot
648 logical :: hdf5_file_format
656 sll_assert_always(nx1 ==
size(eta1, 1))
657 sll_assert_always(nx2 ==
size(eta2, 2))
659 if (
present(iplot))
then
666 write (file_id,
"(a)")
"<Grid Name='mesh' GridType='Uniform'>"
667 write (file_id,
"(a,2i5,a)")
"<Topology TopologyType='2DSMesh' NumberOfElements='", &
669 write (file_id,
"(a)")
"<Geometry GeometryType='X_Y'>"
671 hdf5_file_format = .false.
672 if (
present(file_format))
then
673 if (file_format ==
"HDF5") hdf5_file_format = .true.
675 if (hdf5_file_format)
then
677 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nx2, nx1, &
678 "' NumberType='Float' Precision='8' Format='HDF'>"
679 write (file_id,
"(a)") array_name//
".h5:/x1_values"
680 write (file_id,
"(a)")
"</DataItem>"
681 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nx2, nx1, &
682 "' NumberType='Float' Precision='8' Format='HDF'>"
683 if (
present(iplot))
then
684 write (file_id,
"(a)") array_name//cplot//
".h5:/x2_values"
686 write (file_id,
"(a)") array_name//
".h5:/x2_values"
688 write (file_id,
"(a)")
"</DataItem>"
698 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nx2, nx1, &
699 "' NumberType='Float' Precision='4' Format='XML'>"
701 write (file_id,
"(a)")
"</DataItem>"
702 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nx2, nx1, &
703 "' NumberType='Float' Precision='4' Format='XML'>"
705 write (file_id,
"(a)")
"</DataItem>"
709 write (file_id,
"(a)")
"</Geometry>"
710 write (file_id,
"(a)") &
711 "<Attribute Name='"//array_name//
"' AttributeType='Scalar' Center='Node'>"
714 if (hdf5_file_format)
then
716 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nx2, nx1, &
717 "' NumberType='Float' Precision='8' Format='HDF'>"
718 if (
present(iplot))
then
719 write (file_id,
"(a)") array_name//cplot//
".h5:/node_values"
721 write (file_id,
"(a)") array_name//
".h5:/node_values"
730 write (file_id,
"(a,2i5,a)")
"<DataItem Dimensions='", nx2, nx1, &
731 "' NumberType='Float' Precision='4' Format='XML'>"
738 write (file_id,
"(a)")
"</DataItem>"
739 write (file_id,
"(a)")
"</Attribute>"
770 sll_real64,
intent(in) :: array(:, :, :)
771 sll_real64,
intent(in) :: eta1(:, :, :)
772 sll_real64,
intent(in) :: eta2(:, :, :)
773 sll_real64,
intent(in) :: eta3(:, :, :)
774 character(len=*),
intent(in) :: file_name
775 character(len=*),
intent(in) :: array_name
781 character(len=4),
optional :: file_format
782 sll_int32,
optional :: iplot
784 character(len=4) :: cplot
785 logical :: hdf5_file_format
794 sll_assert_always(nx1 ==
size(eta1, 1))
795 sll_assert_always(nx2 ==
size(eta2, 2))
796 sll_assert_always(nx3 ==
size(eta3, 3))
798 if (
present(iplot))
then
805 write (file_id,
"(a)")
"<Grid Name='mesh' GridType='Uniform'>"
806 write (file_id,
"(a,3i5,a)")
"<Topology TopologyType='3DSMesh' NumberOfElements='", &
808 write (file_id,
"(a)")
"<Geometry GeometryType='X_Y_Z'>"
810 hdf5_file_format = .false.
811 if (
present(file_format))
then
812 if (file_format ==
"HDF5") hdf5_file_format = .true.
814 if (hdf5_file_format)
then
816 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
817 "' NumberType='Float' Precision='8' Format='HDF'>"
818 write (file_id,
"(a)") array_name//
".h5:/x1_values"
819 write (file_id,
"(a)")
"</DataItem>"
820 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
821 "' NumberType='Float' Precision='8' Format='HDF'>"
822 write (file_id,
"(a)") array_name//
".h5:/x2_values"
823 write (file_id,
"(a)")
"</DataItem>"
824 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
825 "' NumberType='Float' Precision='8' Format='HDF'>"
826 write (file_id,
"(a)") array_name//
".h5:/x3_values"
827 write (file_id,
"(a)")
"</DataItem>"
838 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
839 "' NumberType='Float' Precision='4' Format='XML'>"
841 write (file_id,
"(a)")
"</DataItem>"
842 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
843 "' NumberType='Float' Precision='4' Format='XML'>"
845 write (file_id,
"(a)")
"</DataItem>"
846 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
847 "' NumberType='Float' Precision='4' Format='XML'>"
849 write (file_id,
"(a)")
"</DataItem>"
853 write (file_id,
"(a)")
"</Geometry>"
854 if (
present(iplot))
then
855 write (file_id,
"(a)") &
856 "<Attribute Name='"//array_name//cplot//
"' AttributeType='Scalar' Center='Node'>"
858 write (file_id,
"(a)") &
859 "<Attribute Name='"//array_name//
"' AttributeType='Scalar' Center='Node'>"
863 if (hdf5_file_format)
then
865 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
866 "' NumberType='Float' Precision='8' Format='HDF'>"
867 if (
present(iplot))
then
868 write (file_id,
"(a)") array_name//cplot//
".h5:/node_values"
870 write (file_id,
"(a)") array_name//
".h5:/node_values"
879 write (file_id,
"(a,3i5,a)")
"<DataItem Dimensions='", nx3, nx2, nx1, &
880 "' NumberType='Float' Precision='4' Format='XML'>"
887 write (file_id,
"(a)")
"</DataItem>"
888 write (file_id,
"(a)")
"</Attribute>"
895 sll_int32,
intent(in) :: file_id
896 sll_int32,
intent(out) :: error
898 write (file_id,
"(a)")
"</Grid>"
899 write (file_id,
"(a)")
"</Domain>"
900 write (file_id,
"(a)")
"</Xdmf>"
925 sll_int32,
intent(in) :: iplot
926 sll_real64,
dimension(:, :),
intent(in) :: f
927 sll_real64,
dimension(:),
intent(in) :: vec_x1
928 sll_int32,
intent(in) :: nnodes_x1
929 sll_real64,
dimension(:),
intent(in) :: vec_x2
930 sll_int32,
intent(in) :: nnodes_x2
931 character(len=*),
intent(in) :: array_name
936 sll_real64,
dimension(:, :),
allocatable :: x1
937 sll_real64,
dimension(:, :),
allocatable :: x2
940 character(len=4) :: cplot
946 sll_assert_always(iplot > 0)
949 sll_allocate(x1(nnodes_x1, nnodes_x2), error)
950 sll_allocate(x2(nnodes_x1, nnodes_x2), error)
973 call sll_o_xdmf_open(trim(array_name)//cplot//
".xmf",
"cartesian_mesh", &
974 nnodes_x1, nnodes_x2, file_id, error)
976 write (file_id,
"(a,f8.3,a)")
"<Time Value='", time,
"'/>"
979 error, file_id,
"Node")
1009 sll_int32,
intent(in) :: iplot
1010 sll_real64,
intent(in) :: f(:, :)
1011 sll_int32,
intent(in) :: nnodes_x1
1012 sll_int32,
intent(in) :: nnodes_x2
1013 character(len=*),
intent(in) :: array_name
1014 character(len=*),
intent(in) :: mesh_name
1015 sll_real64,
intent(in) :: time
1016 sll_real64,
intent(in),
optional :: x1(:, :)
1017 sll_real64,
intent(in),
optional :: x2(:, :)
1019 sll_int32 :: file_id
1021 character(len=4) :: cplot
1026 if (
present(x1) .and.
present(x2))
then
1041 trim(array_name)//cplot//
".xmf", &
1048 write (file_id,
"(a,f8.3,a)")
"<Time Value='", time,
"'/>"
1051 trim(array_name)//cplot, &
Write nD array in ascii file.
Write a nD array in a binary file.
Write nD array of double precision floats or integers into HDF5 file.
Write the field in xdmf format.
write a data attribute in the xml file
write grid description in the xml file
Module that contains routines to write data in ASCII format file.
Implements the functions to write binary file to store heavy data.
subroutine, public sll_s_binary_file_close(file_id, error)
Open binary file.
subroutine, public sll_s_binary_file_create(filename, file_id, error)
Create binary 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 sll_xdmf_open_2d(file_name, mesh_name, nnodes_x1, nnodes_x2, file_id, error)
Open a XDMF format file for a 2d plot.
subroutine, public sll_s_plot_f(iplot, f, nnodes_x1, nnodes_x2, array_name, mesh_name, time, x1, x2)
Plot 2d distribution function for VisIt.
subroutine, public sll_s_plot_f_cartesian(iplot, f, vec_x1, nnodes_x1, vec_x2, nnodes_x2, array_name, time)
Plot 2d distribution function for VisIt.
subroutine sll_xdmf_open_3d(file_name, mesh_name, nnodes_x1, nnodes_x2, nnodes_x3, file_id, error)
Open a XDMF format file for a 3d plot.
subroutine, public sll_s_xdmf_rect3d_nodes(file_name, array, array_name, eta1, eta2, eta3, file_format, iplot)
Subroutine to write a 3D array in xdmf format. The field is describe on a cartesian mesh....
subroutine, public sll_s_xdmf_close(file_id, error)
Close the XML file and finish to write last lines.
subroutine sll_xdmf_array_2d(mesh_name, array, array_name, error, xmffile_id, center)
Write 2d array in binary or hdf5 file and the matching line in XDMF file.
subroutine sll_xdmf_array_3d(mesh_name, array, array_name, error, xmffile_id, center)
Write 3d array in binary or hdf5 file and the matching line in XDMF file.
subroutine, public sll_s_xdmf_curv2d_nodes(file_name, array, array_name, eta1, eta2, file_format, iplot)
Subroutine to write a 2D array in xdmf format. The field is describe on a cartesian mesh....
subroutine, public sll_s_xdmf_corect2d_nodes(file_name, array, array_name, eta1_min, delta_eta1, eta2_min, delta_eta2, file_format, iplot, time)
Subroutine to write a 2D array in xdmf format The field is describe on a cartesian mesh Axis are perp...
subroutine sll_xdmf_set_time(file_id, time)
Add the the good value of time in VisIt plot.
subroutine, public sll_s_xdmf_rect2d_nodes(file_name, array, array_name, eta1, eta2, file_format, iplot, time)
Subroutine to write a 2D array in xdmf format. The field is describe on a cartesian mesh....
subroutine, public sll_s_xdmf_corect3d_nodes(file_name, array, array_name, eta1_min, delta_eta1, eta2_min, delta_eta2, eta3_min, delta_eta3, file_format, iplot)
Subroutine to write a 3D array in xdmf format The field is describe on a cartesian mesh Axis are perp...
subroutine, public sll_s_xdmf_curv3d_nodes(file_name, array, array_name, eta1, eta2, eta3, file_format, iplot)
Subroutine to write a 3D array in xdmf format. The field is describe on a cartesian mesh....
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.