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_xdmf Module Reference

Description

Implements the functions to write xdmf file plotable by VisIt.

In XDMF (eXtensible Data Model and Format) the description of the data is separate from the values themselves. Light data is stored using XML, Heavy data is stored using HDF5 or Binary files.

Derived types and interfaces

interface  sll_o_xdmf_open
 Create a xmf file. More...
 
interface  sll_o_xdmf_write_array
 Write the field in xdmf format. More...
 

Functions/Subroutines

subroutine sll_xdmf_set_time (file_id, time)
 Add the the good value of time in VisIt plot. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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 perpendicular and spacing is constant. More...
 
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 perpendicular and spacing is constant. More...
 
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. Axis are perpendicular and spacing is define by eta1 and eta2 arrays. More...
 
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. Axis are perpendicular and spacing is define by eta1, eta2 and eta3 arrays. More...
 
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. Nodes coordinates are defined by x and y (2d arrays). More...
 
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. Nodes coordinates are defined by x,y,z (3d arrays). More...
 
subroutine, public sll_s_xdmf_close (file_id, error)
 Close the XML file and finish to write last lines. More...
 
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. More...
 
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. More...
 

Function/Subroutine Documentation

◆ sll_s_plot_f()

subroutine, public sll_m_xdmf::sll_s_plot_f ( integer(kind=i32), intent(in)  iplot,
real(kind=f64), dimension(:, :), intent(in)  f,
integer(kind=i32), intent(in)  nnodes_x1,
integer(kind=i32), intent(in)  nnodes_x2,
character(len=*), intent(in)  array_name,
character(len=*), intent(in)  mesh_name,
real(kind=f64), intent(in)  time,
real(kind=f64), dimension(:, :), intent(in), optional  x1,
real(kind=f64), dimension(:, :), intent(in), optional  x2 
)

Plot 2d distribution function for VisIt.

Parameters
[in]iplotplot counter.
[in]ffunction values .
[in]nnodes_x1nodes number on x1.
[in]nnodes_x2nodes number on x2.
[in]array_namefile name.
[in]mesh_namefile name for mesh
[in]timeAdd time value on plot title.
[in]x12d-array for x1 (create mesh if provided)
[in]x22d-array for x2 (create mesh if provided)

This routine will create a file named array_name_[iplot].xmf TODO suggestion: merge sll_s_plot_f and sll_s_plot_f_cartesian

Definition at line 998 of file sll_m_xdmf.F90.

Here is the call graph for this function:

◆ sll_s_plot_f_cartesian()

subroutine, public sll_m_xdmf::sll_s_plot_f_cartesian ( integer(kind=i32), intent(in)  iplot,
real(kind=f64), dimension(:, :), intent(in)  f,
real(kind=f64), dimension(:), intent(in)  vec_x1,
integer(kind=i32), intent(in)  nnodes_x1,
real(kind=f64), dimension(:), intent(in)  vec_x2,
integer(kind=i32), intent(in)  nnodes_x2,
character(len=*), intent(in)  array_name,
real(kind=f64)  time 
)

Plot 2d distribution function for VisIt.

Parameters
[in]iplotplot counter.
[in]ffunction values .
[in]vec_x1node positions on x1.
[in]nnodes_x1nodes number on x1.
[in]vec_x2node positions on x2.
[in]nnodes_x2nodes number on x2.
[in]array_namefile name.
[in]timeAdd time value on plot title.

This routine will create a file named array_name_[iplot].xmf

Definition at line 916 of file sll_m_xdmf.F90.

Here is the call graph for this function:

◆ sll_s_xdmf_close()

subroutine, public sll_m_xdmf::sll_s_xdmf_close ( integer(kind=i32), intent(in)  file_id,
integer(kind=i32), intent(out)  error 
)

Close the XML file and finish to write last lines.

Parameters
[in]file_idxml file unit number
[out]errorerror code

Definition at line 894 of file sll_m_xdmf.F90.

Here is the caller graph for this function:

◆ sll_s_xdmf_corect2d_nodes()

subroutine, public sll_m_xdmf::sll_s_xdmf_corect2d_nodes ( character(len=*), intent(in)  file_name,
real(kind=f64), dimension(:, :), intent(in)  array,
character(len=*), intent(in)  array_name,
real(kind=f64)  eta1_min,
real(kind=f64)  delta_eta1,
real(kind=f64)  eta2_min,
real(kind=f64)  delta_eta2,
character(len=4), optional  file_format,
integer(kind=i32), optional  iplot,
real(kind=f64), optional  time 
)

Subroutine to write a 2D array in xdmf format The field is describe on a cartesian mesh Axis are perpendicular and spacing is constant.

Parameters
[in]arraydata array
[in]file_namexmf file name
[in]array_namefield name
eta1_minx min
eta2_miny min
delta_eta1dx
delta_eta2dy
file_format"HDF5" or "Binary"
iplotplot index
timetime value

Definition at line 255 of file sll_m_xdmf.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sll_s_xdmf_corect3d_nodes()

subroutine, public sll_m_xdmf::sll_s_xdmf_corect3d_nodes ( character(len=*), intent(in)  file_name,
real(kind=f64), dimension(:, :, :), intent(in)  array,
character(len=*), intent(in)  array_name,
real(kind=f64)  eta1_min,
real(kind=f64)  delta_eta1,
real(kind=f64)  eta2_min,
real(kind=f64)  delta_eta2,
real(kind=f64)  eta3_min,
real(kind=f64)  delta_eta3,
character(len=4), optional  file_format,
integer(kind=i32), optional  iplot 
)

Subroutine to write a 3D array in xdmf format The field is describe on a cartesian mesh Axis are perpendicular and spacing is constant.

Parameters
[in]arraydata array
[in]file_namexmf file name
[in]array_namefield name
eta1_minx min
eta2_miny min
eta3_minz min
delta_eta1dx
delta_eta2dy
delta_eta3dz
file_format"HDF5" or "Binary"
iplotplot index

Definition at line 345 of file sll_m_xdmf.F90.

Here is the call graph for this function:

◆ sll_s_xdmf_curv2d_nodes()

subroutine, public sll_m_xdmf::sll_s_xdmf_curv2d_nodes ( character(len=*), intent(in)  file_name,
real(kind=f64), dimension(:, :), intent(in)  array,
character(len=*), intent(in)  array_name,
real(kind=f64), dimension(:, :), intent(in)  eta1,
real(kind=f64), dimension(:, :), intent(in)  eta2,
character(len=4), optional  file_format,
integer(kind=i32), optional  iplot 
)

Subroutine to write a 2D array in xdmf format. The field is describe on a cartesian mesh. Nodes coordinates are defined by x and y (2d arrays).

Parameters
[in]arraydata array
[in]eta1x data
[in]eta2y data
[in]file_namexmf file name
[in]array_namearray name
file_formatfile format "HDF5" or "Binary"
iplotplot index

Definition at line 628 of file sll_m_xdmf.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sll_s_xdmf_curv3d_nodes()

subroutine, public sll_m_xdmf::sll_s_xdmf_curv3d_nodes ( character(len=*), intent(in)  file_name,
real(kind=f64), dimension(:, :, :), intent(in)  array,
character(len=*), intent(in)  array_name,
real(kind=f64), dimension(:, :, :), intent(in)  eta1,
real(kind=f64), dimension(:, :, :), intent(in)  eta2,
real(kind=f64), dimension(:, :, :), intent(in)  eta3,
character(len=4), optional  file_format,
integer(kind=i32), optional  iplot 
)

Subroutine to write a 3D array in xdmf format. The field is describe on a cartesian mesh. Nodes coordinates are defined by x,y,z (3d arrays).

Example:

call sll_o_xdmf_open(file_name,mesh_name,nnodes_x1,nnodes_x2,nnodes_x3,file_id,error)
call sll_o_xdmf_write_array(mesh_name,x1,'x1',error)
call sll_o_xdmf_write_array(mesh_name,x2,'x2',error)
call sll_o_xdmf_write_array(mesh_name,x3,'x3',error)
call sll_o_xdmf_write_array("field3d",df,"NodeVal",error,file_id,"Node")
call sll_o_xdmf_write_array("field3d",df(1:ncells_x1,1:ncells_x2,1:ncells_x3), &
"CellVal",error,file_id,"Cell")
call sll_s_xdmf_close(file_id,error)
Parameters
[in]arraydata array
[in]eta1x data
[in]eta2y data
[in]eta3z data
[in]file_namexmf file name
[in]array_namearray name
file_formatfile format "HDF5" or "Binary"
iplotplot index

Definition at line 761 of file sll_m_xdmf.F90.

Here is the call graph for this function:

◆ sll_s_xdmf_rect2d_nodes()

subroutine, public sll_m_xdmf::sll_s_xdmf_rect2d_nodes ( character(len=*), intent(in)  file_name,
real(kind=f64), dimension(:, :), intent(in)  array,
character(len=*), intent(in)  array_name,
real(kind=f64), dimension(:), intent(in)  eta1,
real(kind=f64), dimension(:), intent(in)  eta2,
character(len=4), optional  file_format,
integer(kind=i32), optional  iplot,
real(kind=f64), optional  time 
)

Subroutine to write a 2D array in xdmf format. The field is describe on a cartesian mesh. Axis are perpendicular and spacing is define by eta1 and eta2 arrays.

Parameters
[in]arraydata array
[in]eta1x data
[in]eta2y data
[in]file_namexmf file name
[in]array_namearray name
file_formatfile format "HDF5" or "Binary"
iplotplot index

Definition at line 439 of file sll_m_xdmf.F90.

Here is the call graph for this function:

◆ sll_s_xdmf_rect3d_nodes()

subroutine, public sll_m_xdmf::sll_s_xdmf_rect3d_nodes ( character(len=*), intent(in)  file_name,
real(kind=f64), dimension(:, :, :), intent(in)  array,
character(len=*), intent(in)  array_name,
real(kind=f64), dimension(:), intent(in)  eta1,
real(kind=f64), dimension(:), intent(in)  eta2,
real(kind=f64), dimension(:), intent(in)  eta3,
character(len=4), optional  file_format,
integer(kind=i32), optional  iplot 
)

Subroutine to write a 3D array in xdmf format. The field is describe on a cartesian mesh. Axis are perpendicular and spacing is define by eta1, eta2 and eta3 arrays.

Parameters
[in]arraydata array
[in]eta1x data
[in]eta2y data
[in]eta3z data
[in]file_namexmf file name
[in]array_namearray name
file_formatfile format "HDF5" or "Binary"
iplotplot index

Definition at line 531 of file sll_m_xdmf.F90.

Here is the call graph for this function:

◆ sll_xdmf_array_2d()

subroutine sll_m_xdmf::sll_xdmf_array_2d ( character(len=*), intent(in)  mesh_name,
real(kind=f64), dimension(:, :), intent(in)  array,
character(len=*), intent(in)  array_name,
integer(kind=i32), intent(out)  error,
integer(kind=i32), intent(in), optional  xmffile_id,
character(len=4), optional  center 
)
private

Write 2d array in binary or hdf5 file and the matching line in XDMF file.

Parameters
[in]mesh_namemesh name
[in]arraydata array
[in]array_namearray name (hdf5 dataset)
[out]errorerror code
[in]xmffile_idxmf file unit number
center"Node" or "Cell"

Definition at line 145 of file sll_m_xdmf.F90.

◆ sll_xdmf_array_3d()

subroutine sll_m_xdmf::sll_xdmf_array_3d ( character(len=*), intent(in)  mesh_name,
real(kind=f64), dimension(:, :, :), intent(in)  array,
character(len=*), intent(in)  array_name,
integer(kind=i32), intent(out)  error,
integer(kind=i32), intent(in), optional  xmffile_id,
character(len=4), optional  center 
)
private

Write 3d array in binary or hdf5 file and the matching line in XDMF file.

Parameters
[in]mesh_namemesh name
[in]arraydata array
[in]array_namehdf5 dataset name
[out]errorerror code
[in]xmffile_idxmf file unit number
center"Node" or "Cell"

Definition at line 196 of file sll_m_xdmf.F90.

◆ sll_xdmf_open_2d()

subroutine sll_m_xdmf::sll_xdmf_open_2d ( character(len=*), intent(in)  file_name,
character(len=*), intent(in)  mesh_name,
integer(kind=i32), intent(in)  nnodes_x1,
integer(kind=i32), intent(in)  nnodes_x2,
integer(kind=i32), intent(out)  file_id,
integer(kind=i32), intent(out)  error 
)
private

Open a XDMF format file for a 2d plot.

Parameters
[in]file_namexmf file name
[in]mesh_namemesh file name
[in]nnodes_x1x node number
[in]nnodes_x2y node number
[out]file_idfile unit number
[out]errorerror code

Definition at line 106 of file sll_m_xdmf.F90.

◆ sll_xdmf_open_3d()

subroutine sll_m_xdmf::sll_xdmf_open_3d ( character(len=*), intent(in)  file_name,
character(len=*), intent(in)  mesh_name,
integer(kind=i32), intent(in)  nnodes_x1,
integer(kind=i32), intent(in)  nnodes_x2,
integer(kind=i32), intent(in)  nnodes_x3,
integer(kind=i32), intent(out)  file_id,
integer(kind=i32), intent(out)  error 
)
private

Open a XDMF format file for a 3d plot.

Parameters
[in]file_namexmf file name
[in]mesh_namemesh name
[out]file_idfile unit number
[out]errorerror code
[in]nnodes_x1x nodes number
[in]nnodes_x2y nodes number
[in]nnodes_x3z nodes number

Definition at line 122 of file sll_m_xdmf.F90.

◆ sll_xdmf_set_time()

subroutine sll_m_xdmf::sll_xdmf_set_time ( integer(kind=i32), intent(in)  file_id,
real(kind=f64), intent(in)  time 
)
private

Add the the good value of time in VisIt plot.

Parameters
[in]timeinput time
[in]file_idfile unit number

Definition at line 90 of file sll_m_xdmf.F90.

    Report Typos and Errors