Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_xml_io.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! Pierre Navaro
5 !
6 ! This code SeLaLib (for Semi-Lagrangian-Library)
7 ! is a parallel library for simulating the plasma turbulence
8 ! in a tokamak.
9 !
10 ! This software is governed by the CeCILL-B license
11 ! under French law and abiding by the rules of distribution
12 ! of free software. You can use, modify and redistribute
13 ! the software under the terms of the CeCILL-B license as
14 ! circulated by CEA, CNRS and INRIA at the following URL
15 ! "http://www.cecill.info".
16 !**************************************************************
17 
24 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
25 #include "sll_assert.h"
26 #include "sll_working_precision.h"
27 
28  implicit none
29 
30  public :: &
35 
36  private
37 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
38 
40  interface sll_xml_dataitem
41  module procedure sll_xml_dataitem_1d
42  module procedure sll_xml_dataitem_2d
43  module procedure sll_xml_dataitem_3d
44  end interface
45 
47  interface sll_o_xml_field
48  module procedure sll_xml_field_1d
49  module procedure sll_xml_field_2d
50  module procedure sll_xml_field_3d
51  end interface
52 
56  module procedure sll_xml_grid_geometry_2d_low_level
58  module procedure sll_xml_grid_geometry_3d_low_level
59  end interface
60 
61 contains
62 
65  subroutine sll_s_xml_file_create(filename, file_id, error)
66 
67  character(len=*), intent(in) :: filename
68  sll_int32, intent(out) :: file_id
69  sll_int32, intent(out) :: error
70  logical :: lopen
71 
72  error = 0
73  do 100 file_id = 20, 99
74  inquire (unit=file_id, opened=lopen)
75  if (lopen) then
76  cycle
77  else
78  open (file_id, status='SCRATCH', err=100)
79  close (file_id, status='DELETE', err=100)
80  goto 200
81  end if
82 
83 100 continue
84  error = 1
85 200 continue
86  error = 0
87 
88  ! Create a new file using default properties
89  inquire (file=filename, opened=lopen)
90  sll_assert_always(.not. lopen)
91 
92  open (file_id, file=filename, form='FORMATTED', iostat=error)
93  rewind(file_id)
94 
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>"
99 
100  end subroutine sll_s_xml_file_create
101 
106  subroutine sll_s_xml_file_close(file_id, error)
107 
108  sll_int32, intent(in) :: file_id
109  sll_int32, intent(out) :: error
110 
111  write (file_id, "(a)") "</Grid>"
112  write (file_id, "(a)") "</Domain>"
113  write (file_id, "(a)") "</Xdmf>"
114  close (file_id)
115  error = 0
116  end subroutine sll_s_xml_file_close
117 
127  subroutine sll_xml_dataitem_1d(file_id, &
128  filename, &
129  nnodes_x1, &
130  filetype)
131 
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
136 
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>"
142  end subroutine sll_xml_dataitem_1d
143 
155  subroutine sll_xml_dataitem_2d(file_id, &
156  filename, &
157  nnodes_x1, &
158  nnodes_x2, &
159  filetype)
160 
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
166 
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>"
172  end subroutine sll_xml_dataitem_2d
173 
175  subroutine sll_xml_dataitem_3d(file_id, &
176  filename, &
177  nnodes_x1, &
178  nnodes_x2, &
179  nnodes_x3, &
180  filetype)
181 
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
188 
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>"
195  end subroutine sll_xml_dataitem_3d
196 
209  subroutine sll_xml_field_1d(file_id, &
210  fieldname, &
211  filename, &
212  npoints_1, &
213  filetype, &
214  center)
215 
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
222 
223  write (file_id, "(a)") &
224  "<Attribute Name='"//trim(fieldname)//"' AttributeType='Scalar' Center='"//center//"'>"
225  call sll_xml_dataitem_1d(file_id, filename, npoints_1, filetype)
226  write (file_id, "(a)") "</Attribute>"
227  end subroutine sll_xml_field_1d
228 
242  subroutine sll_xml_field_2d(file_id, &
243  fieldname, &
244  filename, &
245  npoints_1, &
246  npoints_2, &
247  filetype, &
248  center)
249 
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
257 
258  write (file_id, "(a)") &
259  "<Attribute Name='"//trim(fieldname)//"' AttributeType='Scalar' Center='"//center//"'>"
260  call sll_xml_dataitem_2d(file_id, filename, npoints_1, npoints_2, filetype)
261  write (file_id, "(a)") "</Attribute>"
262  end subroutine sll_xml_field_2d
263 
279  subroutine sll_xml_field_3d(file_id, &
280  fieldname, &
281  filename, &
282  npoints_1, &
283  npoints_2, &
284  npoints_3, &
285  filetype, &
286  center)
287 
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
296 
297  write (file_id, "(a)") &
298  "<Attribute Name='"//trim(fieldname)//"' AttributeType='Scalar' Center='"//center//"'>"
299  call sll_xml_dataitem_3d(file_id, &
300  filename, &
301  npoints_1, &
302  npoints_2, &
303  npoints_3, &
304  filetype)
305  write (file_id, "(a)") "</Attribute>"
306  end subroutine sll_xml_field_3d
307 
319  file_id, &
320  filename, &
321  nnodes_x1, &
322  nnodes_x2)
323 
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
328 
329 #ifdef NOHDF5
330  call sll_xml_grid_geometry_2d_low_level(file_id, &
331  trim(filename)//"-x1.bin", nnodes_x1, &
332  trim(filename)//"-x2.bin", nnodes_x2, "x1", "x2", 'Uniform')
333 #else
334  call sll_xml_grid_geometry_2d_low_level(file_id, &
335  trim(filename)//"-x1.h5", nnodes_x1, &
336  trim(filename)//"-x2.h5", nnodes_x2, "x1", "x2", 'Uniform')
337 #endif
338 
340 
355  x1filename, &
356  nnodes_x1, &
357  x2filename, &
358  nnodes_x2, &
359  x1dsetname, &
360  x2dsetname, &
361  gridtype)
362 
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
371 
372  write (file_id, "(a)") "<Grid Name='mesh' GridType='"//gridtype//"'>"
373  write (file_id, &
374  "(a,2i5,a)") "<Topology TopologyType='2DSMesh' NumberOfElements='", &
375  nnodes_x2, nnodes_x1, "'/>"
376  write (file_id, "(a)") "<Geometry GeometryType='X_Y'>"
377 
378 #ifdef NOHDF5
379  call sll_xml_dataitem_2d(file_id, &
380  trim(x1filename), &
381  nnodes_x1, nnodes_x2, 'Binary')
382  call sll_xml_dataitem_2d(file_id, &
383  trim(x2filename), &
384  nnodes_x1, nnodes_x2, 'Binary')
385 #else
386  call sll_xml_dataitem_2d(file_id, &
387  trim(x1filename)//":/"//x1dsetname, &
388  nnodes_x1, nnodes_x2, 'HDF')
389  call sll_xml_dataitem_2d(file_id, &
390  trim(x2filename)//":/"//x2dsetname, &
391  nnodes_x1, nnodes_x2, 'HDF')
392 #endif
393 
394  write (file_id, "(a)") "</Geometry>"
396 
400  subroutine sll_xml_grid_geometry_3d_high_level(file_id, filename, &
401  nnodes_x1, nnodes_x2, nnodes_x3)
402 
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
408 
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'>"
413 
414 #ifdef NOHDF5
415 
416  call sll_xml_dataitem_3d(file_id, trim(filename)//"-x1.bin", &
417  nnodes_x1, nnodes_x2, nnodes_x3, 'Binary')
418  call sll_xml_dataitem_3d(file_id, trim(filename)//"-x2.bin", &
419  nnodes_x1, nnodes_x2, nnodes_x3, 'Binary')
420  call sll_xml_dataitem_3d(file_id, trim(filename)//"-x3.bin", &
421  nnodes_x1, nnodes_x2, nnodes_x3, 'Binary')
422 #else
423 
424  call sll_xml_dataitem_3d(file_id, trim(filename)//"-x1.h5:/x1", &
425  nnodes_x1, nnodes_x2, nnodes_x3, 'HDF')
426  call sll_xml_dataitem_3d(file_id, trim(filename)//"-x2.h5:/x2", &
427  nnodes_x1, nnodes_x2, nnodes_x3, 'HDF')
428  call sll_xml_dataitem_3d(file_id, trim(filename)//"-x3.h5:/x3", &
429  nnodes_x1, nnodes_x2, nnodes_x3, 'HDF')
430 
431 #endif
432 
433  write (file_id, "(a)") "</Geometry>"
434 
436 
441  x1filename, nnodes_x1, &
442  x2filename, nnodes_x2, &
443  x3filename, nnodes_x3, &
444  x1dsetname, x2dsetname, x3dsetname, &
445  gridtype)
446 
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
458 
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'>"
464 
465 #ifdef NOHDF5
466 
467  call sll_xml_dataitem_3d(file_id, x1filename, &
468  nnodes_x1, nnodes_x2, nnodes_x3, 'Binary')
469  call sll_xml_dataitem_3d(file_id, x2filename, &
470  nnodes_x1, nnodes_x2, nnodes_x3, 'Binary')
471  call sll_xml_dataitem_3d(file_id, x3filename, &
472  nnodes_x1, nnodes_x2, nnodes_x3, 'Binary')
473 #else
474 
475  call sll_xml_dataitem_3d(file_id, x1filename//":/"//x1dsetname, &
476  nnodes_x1, nnodes_x2, nnodes_x3, 'HDF')
477  call sll_xml_dataitem_3d(file_id, x2filename//":/"//x2dsetname, &
478  nnodes_x1, nnodes_x2, nnodes_x3, 'HDF')
479  call sll_xml_dataitem_3d(file_id, x3filename//":/"//x3dsetname, &
480  nnodes_x1, nnodes_x2, nnodes_x3, 'HDF')
481 
482 #endif
483 
484  write (file_id, "(a)") "</Geometry>"
485 
487 
488  end module sll_m_xml_io
489 
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.
    Report Typos and Errors