Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
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... | |
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.
[in] | iplot | plot counter. |
[in] | f | function values . |
[in] | nnodes_x1 | nodes number on x1. |
[in] | nnodes_x2 | nodes number on x2. |
[in] | array_name | file name. |
[in] | mesh_name | file name for mesh |
[in] | time | Add time value on plot title. |
[in] | x1 | 2d-array for x1 (create mesh if provided) |
[in] | x2 | 2d-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.
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.
[in] | iplot | plot counter. |
[in] | f | function values . |
[in] | vec_x1 | node positions on x1. |
[in] | nnodes_x1 | nodes number on x1. |
[in] | vec_x2 | node positions on x2. |
[in] | nnodes_x2 | nodes number on x2. |
[in] | array_name | file name. |
[in] | time | Add 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.
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.
[in] | file_id | xml file unit number |
[out] | error | error code |
Definition at line 894 of file sll_m_xdmf.F90.
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.
[in] | array | data array |
[in] | file_name | xmf file name |
[in] | array_name | field name |
eta1_min | x min | |
eta2_min | y min | |
delta_eta1 | dx | |
delta_eta2 | dy | |
file_format | "HDF5" or "Binary" | |
iplot | plot index | |
time | time value |
Definition at line 255 of file sll_m_xdmf.F90.
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.
[in] | array | data array |
[in] | file_name | xmf file name |
[in] | array_name | field name |
eta1_min | x min | |
eta2_min | y min | |
eta3_min | z min | |
delta_eta1 | dx | |
delta_eta2 | dy | |
delta_eta3 | dz | |
file_format | "HDF5" or "Binary" | |
iplot | plot index |
Definition at line 345 of file sll_m_xdmf.F90.
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).
[in] | array | data array |
[in] | eta1 | x data |
[in] | eta2 | y data |
[in] | file_name | xmf file name |
[in] | array_name | array name |
file_format | file format "HDF5" or "Binary" | |
iplot | plot index |
Definition at line 628 of file sll_m_xdmf.F90.
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:
[in] | array | data array |
[in] | eta1 | x data |
[in] | eta2 | y data |
[in] | eta3 | z data |
[in] | file_name | xmf file name |
[in] | array_name | array name |
file_format | file format "HDF5" or "Binary" | |
iplot | plot index |
Definition at line 761 of file sll_m_xdmf.F90.
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.
[in] | array | data array |
[in] | eta1 | x data |
[in] | eta2 | y data |
[in] | file_name | xmf file name |
[in] | array_name | array name |
file_format | file format "HDF5" or "Binary" | |
iplot | plot index |
Definition at line 439 of file sll_m_xdmf.F90.
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.
[in] | array | data array |
[in] | eta1 | x data |
[in] | eta2 | y data |
[in] | eta3 | z data |
[in] | file_name | xmf file name |
[in] | array_name | array name |
file_format | file format "HDF5" or "Binary" | |
iplot | plot index |
Definition at line 531 of file sll_m_xdmf.F90.
|
private |
Write 2d array in binary or hdf5 file and the matching line in XDMF file.
[in] | mesh_name | mesh name |
[in] | array | data array |
[in] | array_name | array name (hdf5 dataset) |
[out] | error | error code |
[in] | xmffile_id | xmf file unit number |
center | "Node" or "Cell" |
Definition at line 145 of file sll_m_xdmf.F90.
|
private |
Write 3d array in binary or hdf5 file and the matching line in XDMF file.
[in] | mesh_name | mesh name |
[in] | array | data array |
[in] | array_name | hdf5 dataset name |
[out] | error | error code |
[in] | xmffile_id | xmf file unit number |
center | "Node" or "Cell" |
Definition at line 196 of file sll_m_xdmf.F90.
|
private |
Open a XDMF format file for a 2d plot.
[in] | file_name | xmf file name |
[in] | mesh_name | mesh file name |
[in] | nnodes_x1 | x node number |
[in] | nnodes_x2 | y node number |
[out] | file_id | file unit number |
[out] | error | error code |
Definition at line 106 of file sll_m_xdmf.F90.
|
private |
Open a XDMF format file for a 3d plot.
[in] | file_name | xmf file name |
[in] | mesh_name | mesh name |
[out] | file_id | file unit number |
[out] | error | error code |
[in] | nnodes_x1 | x nodes number |
[in] | nnodes_x2 | y nodes number |
[in] | nnodes_x3 | z nodes number |
Definition at line 122 of file sll_m_xdmf.F90.
|
private |
Add the the good value of time in VisIt plot.
[in] | time | input time |
[in] | file_id | file unit number |
Definition at line 90 of file sll_m_xdmf.F90.