Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Derived types and interfaces | Functions/Subroutines
sll_m_xml_io Module Reference

Description

Implements the functions to write xml file to store light data.

With XDMF file you can describe data to plot them with VisIt

Derived types and interfaces

interface  sll_xml_dataitem
 write a data item in the xml file More...
 
interface  sll_o_xml_field
 write a data attribute in the xml file More...
 
interface  sll_o_xml_grid_geometry
 write grid description in the xml file More...
 

Functions/Subroutines

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. More...
 
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. More...
 
subroutine sll_xml_dataitem_1d (file_id, filename, nnodes_x1, filetype)
 Write the description of a scalar field on a 1D mesh. More...
 
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 mesh. More...
 
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. More...
 
subroutine sll_xml_field_1d (file_id, fieldname, filename, npoints_1, filetype, center)
 Write the description of a scalar field on a 1D mesh. More...
 
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. More...
 
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. More...
 
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 and filename-x2. More...
 
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 and filename-x2. More...
 
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 filename-x1 and filename-x2. High level version where dataset names in hdf5 files are set automatically. More...
 
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 filename-x1 and filename-x2. Low level version where dataset names in hdf5 files must be set. More...
 

Function/Subroutine Documentation

◆ sll_s_xml_file_close()

subroutine, public sll_m_xml_io::sll_s_xml_file_close ( integer(kind=i32), intent(in)  file_id,
integer(kind=i32), intent(out)  error 
)

Close the XML file and finish to write last lines. You give the file unit number.

Parameters
[in]file_id- the unit number or your xml file
[out]error- error parameter
[in]file_idxml file unit number
[out]errorerror code

Definition at line 106 of file sll_m_xml_io.F90.

Here is the caller graph for this function:

◆ sll_s_xml_file_create()

subroutine, public sll_m_xml_io::sll_s_xml_file_create ( character(len=*), intent(in)  filename,
integer(kind=i32), intent(out)  file_id,
integer(kind=i32), intent(out)  error 
)

Create the XML file and begin to write first lines. You get the file unit number.

Parameters
[in]filenamefile name
[out]file_idfile unit number
[out]errorerror code

Definition at line 65 of file sll_m_xml_io.F90.

Here is the caller graph for this function:

◆ sll_xml_dataitem_1d()

subroutine sll_m_xml_io::sll_xml_dataitem_1d ( integer(kind=i32), intent(in)  file_id,
character(len=*), intent(in)  filename,
integer(kind=i32), intent(in)  nnodes_x1,
character(len=*), intent(in)  filetype 
)
private

Write the description of a scalar field on a 1D mesh.

Parameters
[in]file_idis the unit number or your xml file
[in]filenameis the file name where the heavy data are (bin or h5)
[in]nnodes_x1- nodes number along direction 1
[in]filetype- heavy data format 'HDF' or 'Binary'

The file named filename must exist.

Definition at line 127 of file sll_m_xml_io.F90.

Here is the caller graph for this function:

◆ sll_xml_dataitem_2d()

subroutine sll_m_xml_io::sll_xml_dataitem_2d ( integer(kind=i32), intent(in)  file_id,
character(len=*), intent(in)  filename,
integer(kind=i32), intent(in)  nnodes_x1,
integer(kind=i32), intent(in)  nnodes_x2,
character(len=*), intent(in)  filetype 
)
private

Write the description of a scalar field on a 3D mesh. Write the description of a scalar field on a 2D mesh.

Parameters
[in]file_idis the unit number or your xml file
[in]filenameis the file name where the heavy data are (bin or h5)
[in]nnodes_x1- nodes number along direction 1
[in]nnodes_x2- nodes number along direction 2
[in]filetype- heavy data format 'HDF' or 'Binary'

The file named filename must exist.

Definition at line 155 of file sll_m_xml_io.F90.

Here is the caller graph for this function:

◆ sll_xml_dataitem_3d()

subroutine sll_m_xml_io::sll_xml_dataitem_3d ( integer(kind=i32), intent(in)  file_id,
character(len=*), intent(in)  filename,
integer(kind=i32), intent(in)  nnodes_x1,
integer(kind=i32), intent(in)  nnodes_x2,
integer(kind=i32), intent(in)  nnodes_x3,
character(len=*), intent(in)  filetype 
)
private

Write the description of a scalar field on a 3D mesh.

Parameters
[in]file_idfile unit number
[in]filenamexmf file name
[in]filetype"HDF" or "Binary"
[in]nnodes_x1x nodes number
[in]nnodes_x2y nodes number
[in]nnodes_x3z nodes number

Definition at line 175 of file sll_m_xml_io.F90.

Here is the caller graph for this function:

◆ sll_xml_field_1d()

subroutine sll_m_xml_io::sll_xml_field_1d ( integer(kind=i32), intent(in)  file_id,
character(len=*), intent(in)  fieldname,
character(len=*), intent(in)  filename,
integer(kind=i32), intent(in)  npoints_1,
character(len=*), intent(in)  filetype,
character(len=*), intent(in)  center 
)
private

Write the description of a scalar field on a 1D mesh.

Parameters
[in]fieldnamethe dataset name where the heavy data are (hdf5 case)
[in]filenamethe file name where the heavy data are (bin or h5)
[in]npoints_1nodes or cells number along direction 1
[in]centervalues are centered on nodes or cells

The file named filename-fieldname.bin must exist in case of binary output. The file named filename.h5 with dataset fieldname must exist in case of hdf5 output.

Parameters
[in]file_idthe unit number or your xml file
[in]filetype"HDF" or "Binary"

Definition at line 209 of file sll_m_xml_io.F90.

◆ sll_xml_field_2d()

subroutine sll_m_xml_io::sll_xml_field_2d ( integer(kind=i32), intent(in)  file_id,
character(len=*), intent(in)  fieldname,
character(len=*), intent(in)  filename,
integer(kind=i32), intent(in)  npoints_1,
integer(kind=i32), intent(in)  npoints_2,
character(len=*), intent(in)  filetype,
character(len=*), intent(in)  center 
)
private

Write the description of a scalar field on a 2D mesh.

Parameters
[in]fieldnamethe dataset name where the heavy data are (hdf5 case)
[in]filenamethe file name where the heavy data are (bin or h5)
[in]npoints_1nodes or cells number along direction 1
[in]npoints_2nodes or cells number along direction 2
[in]centervalues are centered on nodes or cells

The file named filename-fieldname.bin must exist in case of binary output. The file named filename.h5 with dataset fieldname must exist in case of hdf5 output.

Parameters
[in]file_idthe unit number or your xml file
[in]filetype"HDF" or "Binary"

Definition at line 242 of file sll_m_xml_io.F90.

◆ sll_xml_field_3d()

subroutine sll_m_xml_io::sll_xml_field_3d ( integer(kind=i32), intent(in)  file_id,
character(len=*), intent(in)  fieldname,
character(len=*), intent(in)  filename,
integer(kind=i32), intent(in)  npoints_1,
integer(kind=i32), intent(in)  npoints_2,
integer(kind=i32), intent(in)  npoints_3,
character(len=*), intent(in)  filetype,
character(len=*), intent(in)  center 
)
private

Write the description of a scalar field on a 3D mesh.

Parameters
[in]file_idthe unit number or your xml file
[in]fieldnamethe dataset name where the heavy data are (hdf5 case)
[in]filenamethe file name where the heavy data are (bin or h5)
[in]npoints_1nodes or cells number along direction 1
[in]npoints_2nodes or cells number along direction 2
[in]npoints_3nodes or cells number along direction 3
[in]centervalues are centered on nodes or cells

The file named filename-fieldname.bin must exist in case of binary output. The file named filename.h5 with dataset fieldname must exist in case of hdf5 output.

Parameters
[in]filetype"HDF" or "Binary"

Definition at line 279 of file sll_m_xml_io.F90.

◆ sll_xml_grid_geometry_2d_high_level()

subroutine sll_m_xml_io::sll_xml_grid_geometry_2d_high_level ( integer(kind=i32), intent(in)  file_id,
character(len=*), intent(in)  filename,
integer(kind=i32), intent(in)  nnodes_x1,
integer(kind=i32), intent(in)  nnodes_x2 
)
private

Write the description of a 2D strutured grid mesh with its nodes coordinates contains in filename-x1 and filename-x2.

Parameters
[in]file_idis the unit number or your xml file
[in]filenameis the file name where the coordinates data are (bin or h5)
[in]nnodes_x1- nodes number along direction 1
[in]nnodes_x2- nodes number along direction 2

The file named filename-x1.bin and filename-x2.bin must exist in case of binary output. The file named filename.h5 with dataset x1 and x2 must exist in case of hdf5 output.

Definition at line 318 of file sll_m_xml_io.F90.

◆ sll_xml_grid_geometry_2d_low_level()

subroutine sll_m_xml_io::sll_xml_grid_geometry_2d_low_level ( integer(kind=i32), intent(in)  file_id,
character(len=*), intent(in)  x1filename,
integer(kind=i32), intent(in)  nnodes_x1,
character(len=*), intent(in)  x2filename,
integer(kind=i32), intent(in)  nnodes_x2,
character(len=*), intent(in)  x1dsetname,
character(len=*), intent(in)  x2dsetname,
character(len=*), intent(in)  gridtype 
)
private

Write the description of a 2D strutured grid mesh with its nodes coordinates contains in filename-x1 and filename-x2.

Parameters
[in]file_idis the unit number or your xml file
[in]x1filenameis the file name where the coordinates x1 are (bin or h5)
[in]x2filenameis the file name where the coordinates x2 are (bin or h5)
[in]x1dsetnameis the dataset name of coordinates x1 are (bin or h5)
[in]x2dsetnameis the dataset name of coordinates x2 are (bin or h5)
[in]nnodes_x1- nodes number along direction 1
[in]nnodes_x2- nodes number along direction 2

The file named x*filename-x*dsetname.bin must exists The file named x*filename-x*dsetname.h5 with dataset x*dsetname must exists.
Low level version where you have to set dataset names in hdf5 files

Definition at line 354 of file sll_m_xml_io.F90.

Here is the caller graph for this function:

◆ sll_xml_grid_geometry_3d_high_level()

subroutine sll_m_xml_io::sll_xml_grid_geometry_3d_high_level ( integer(kind=i32), intent(in)  file_id,
character(len=*), intent(in)  filename,
integer(kind=i32), intent(in)  nnodes_x1,
integer(kind=i32), intent(in)  nnodes_x2,
integer(kind=i32), intent(in)  nnodes_x3 
)
private

Write the description of a 3D structured curvilinear grid mesh with its nodes coordinates contains in filename-x1 and filename-x2. High level version where dataset names in hdf5 files are set automatically.

Parameters
[in]file_idxmf file unit number
[in]filenamexmf file name
[in]nnodes_x1x nodes number
[in]nnodes_x2y nodes number
[in]nnodes_x3z nodes number

Definition at line 400 of file sll_m_xml_io.F90.

◆ sll_xml_grid_geometry_3d_low_level()

subroutine sll_m_xml_io::sll_xml_grid_geometry_3d_low_level ( integer(kind=i32), intent(in)  file_id,
character(len=*), intent(in)  x1filename,
integer(kind=i32), intent(in)  nnodes_x1,
character(len=*), intent(in)  x2filename,
integer(kind=i32), intent(in)  nnodes_x2,
character(len=*), intent(in)  x3filename,
integer(kind=i32), intent(in)  nnodes_x3,
character(len=*), intent(in)  x1dsetname,
character(len=*), intent(in)  x2dsetname,
character(len=*), intent(in)  x3dsetname,
character(len=*), intent(in)  gridtype 
)
private

Write the description of a 3D structured curvilinear grid mesh with its nodes coordinates contains in filename-x1 and filename-x2. Low level version where dataset names in hdf5 files must be set.

Parameters
[in]file_idxmf unif file number
[in]x1filenamex data file name
[in]x2filenamey data file name
[in]x3filenamex datz file name
[in]x1dsetnamex dataset name
[in]x2dsetnamey dataset name
[in]x3dsetnamez dataset name
[in]nnodes_x1x nodes number
[in]nnodes_x2y nodes number
[in]nnodes_x3z nodes number
[in]gridtypeuniform or collection

Definition at line 440 of file sll_m_xml_io.F90.

    Report Typos and Errors