PB3D
[2.45]
Ideal linear high-n MHD stability in 3-D
|
Operations on HDF5 and XDMF variables. More...
Interfaces and Types | |
interface | reset_hdf5_item |
Resets an HDF5 XDMF item. More... | |
Functions/Subroutines | |
integer function, public | open_hdf5_file (file_info, file_name, sym_type, descr, ind_plot, cont_plot) |
Opens an HDF5 file and accompanying xmf file and returns the handles. More... | |
integer function, public | close_hdf5_file (file_info, ind_plot, cont_plot) |
Closes an HDF5 file and writes the accompanying xmf file. More... | |
integer function, public | add_hdf5_item (file_info, XDMF_item, reset, ind_plot) |
Add an XDMF item to a HDF5 file. More... | |
integer function, public | print_hdf5_3d_data_item (dataitem_id, file_info, var_name, var, dim_tot, loc_dim, loc_offset, init_val, ind_plot, cont_plot) |
Prints an HDF5 data item. More... | |
subroutine, public | merge_hdf5_3d_data_items (merged_id, dataitem_ids, var_name, dim_tot, reset, ind_plot, cont_plot) |
Joins dataitems in a vector. More... | |
subroutine, public | print_hdf5_att (att_id, att_dataitem, att_name, att_center, att_type, reset, ind_plot) |
Prints an HDF5 attribute. More... | |
subroutine, public | print_hdf5_top (top_id, top_type, top_n_elem, ind_plot) |
Prints an HDF5 topology. More... | |
subroutine, public | print_hdf5_geom (geom_id, geom_type, geom_dataitems, reset, ind_plot) |
Prints an HDF5 geometry. More... | |
integer function, public | print_hdf5_grid (grid_id, grid_name, grid_type, grid_time, grid_top, grid_geom, grid_atts, grid_grids, reset, ind_plot) |
Prints an HDF5 grid. More... | |
integer function, public | create_output_hdf5 (HDF5_name) |
Creates an HDF5 output file. More... | |
integer function, public | print_hdf5_arrs (vars, PB3D_name, head_name, rich_lvl, disp_info, ind_print, remove_previous_arrs) |
Prints a series of arrays, in the form of an array of pointers, to an HDF5 file. More... | |
integer function, public | read_hdf5_arr (var, PB3D_name, head_name, var_name, rich_lvl, disp_info, lim_loc) |
Reads a PB3D output file in HDF5 format. More... | |
Variables | |
logical, public | debug_hdf5_ops = .false. |
set to true to debug general information More... | |
logical, public | debug_print_hdf5_arrs = .false. |
set to true to debug print_HDF5_arrs More... | |
Operations on HDF5 and XDMF variables.
X_2
variables get larger than that. https://www.hdfgroup.org/HDF5/doc/UG/OldHtmlSource/UG_frame08TheFile.htmlinteger function, public hdf5_ops::add_hdf5_item | ( | type(hdf5_file_type), intent(inout) | file_info, |
type(xml_str_type), intent(inout) | XDMF_item, | ||
logical, intent(in), optional | reset, | ||
logical, intent(in), optional | ind_plot | ||
) |
Add an XDMF item to a HDF5 file.
[in,out] | file_info | info about HDF5 file |
[in,out] | xdmf_item | XDMF item to add |
[in] | reset | if .true., XDMF_item is reset |
[in] | ind_plot | whether individual plot |
Definition at line 340 of file HDF5_ops.f90.
integer function, public hdf5_ops::close_hdf5_file | ( | type(hdf5_file_type), intent(inout) | file_info, |
logical, intent(in), optional | ind_plot, | ||
logical, intent(in), optional | cont_plot | ||
) |
Closes an HDF5 file and writes the accompanying xmf file.
[in,out] | file_info | info about HDF5 file |
[in] | ind_plot | whether individual plot |
[in] | cont_plot | continued plot |
Definition at line 264 of file HDF5_ops.f90.
integer function, public hdf5_ops::create_output_hdf5 | ( | character(len=*), intent(in) | HDF5_name | ) |
Creates an HDF5 output file.
[in] | hdf5_name | name of HDF5 file |
Definition at line 1075 of file HDF5_ops.f90.
subroutine, public hdf5_ops::merge_hdf5_3d_data_items | ( | type(xml_str_type), intent(inout) | merged_id, |
type(xml_str_type), dimension(:), intent(inout) | dataitem_ids, | ||
character(len=*), intent(in) | var_name, | ||
integer, dimension(3), intent(in) | dim_tot, | ||
logical, intent(in), optional | reset, | ||
logical, intent(in), optional | ind_plot, | ||
logical, intent(in), optional | cont_plot | ||
) |
Joins dataitems in a vector.
[in,out] | merged_id | ID of merged data item |
[in,out] | dataitem_ids | ID of data items |
[in] | var_name | name of variable |
[in] | dim_tot | total dimensions of variable |
[in] | reset | whether to reset dataitems |
[in] | ind_plot | whether individual plot |
[in] | cont_plot | continued plot |
Definition at line 603 of file HDF5_ops.f90.
integer function, public hdf5_ops::open_hdf5_file | ( | type(hdf5_file_type), intent(inout) | file_info, |
character(len=*), intent(in) | file_name, | ||
integer, intent(inout), optional | sym_type, | ||
character(len=*), intent(in), optional | descr, | ||
logical, intent(in), optional | ind_plot, | ||
logical, intent(in), optional | cont_plot | ||
) |
Opens an HDF5 file and accompanying xmf file and returns the handles.
Optionally, a description of the file can be provided.
Also, the plot can be done for only one process, setting the variable ind_plot
.
Furthermore, if the plot is a continuation, using cont_plot
, the previous plot is opened and sym_type
is returned.
[in,out] | file_info | info about HDF5 file |
[in] | file_name | name of HDF5 file |
[in,out] | sym_type | symmetry type |
[in] | descr | description of file |
[in] | ind_plot | whether individual plot |
[in] | cont_plot | continued plot |
Definition at line 78 of file HDF5_ops.f90.
integer function, public hdf5_ops::print_hdf5_3d_data_item | ( | type(xml_str_type), intent(inout) | dataitem_id, |
type(hdf5_file_type), intent(in) | file_info, | ||
character(len=*), intent(in) | var_name, | ||
real(dp), dimension(:,:,:), intent(in) | var, | ||
integer, dimension(3), intent(in) | dim_tot, | ||
integer, dimension(3), intent(in), optional | loc_dim, | ||
integer, dimension(3), intent(in), optional | loc_offset, | ||
real(dp), intent(in), optional | init_val, | ||
logical, intent(in), optional | ind_plot, | ||
logical, intent(in), optional | cont_plot | ||
) |
Prints an HDF5 data item.
If this is a parallel data item, the group dimension and offset have to be specified as well.
[in,out] | dataitem_id | ID of data item |
[in] | file_info | info about HDF5 file |
[in] | var_name | name of variable |
[in] | var | variable to write |
[in] | dim_tot | total dimensions of variable |
[in] | loc_dim | dimensions in this group |
[in] | loc_offset | offset in this group |
[in] | init_val | initial fill value |
[in] | ind_plot | whether individual plot |
[in] | cont_plot | continued plot |
Definition at line 396 of file HDF5_ops.f90.
integer function, public hdf5_ops::print_hdf5_arrs | ( | type(var_1d_type), dimension(:), intent(in) | vars, |
character(len=*), intent(in) | PB3D_name, | ||
character(len=*), intent(in) | head_name, | ||
integer, intent(in), optional | rich_lvl, | ||
logical, intent(in), optional | disp_info, | ||
logical, intent(in), optional | ind_print, | ||
logical, intent(in), optional | remove_previous_arrs | ||
) |
Prints a series of arrays, in the form of an array of pointers, to an HDF5 file.
Optionally, output can be given about the variable being written. Also, if rich_lvl
is provided, _R_[rich_lvl]
is appended to the head name if it is > 0.
Also optionally, previously encountered arrays can be removed. This is to be used with care, as it may disturb the internal workings of PB3D. Currently (v2.15) it is used only for the specific case of jumping to solutions for X_grid_style
1 or 3, when writing solutions and solution grids.
[in] | vars | variables to write |
[in] | pb3d_name | name of PB3D file |
[in] | head_name | head name of variables |
[in] | rich_lvl | Richardson level to append to file name |
[in] | disp_info | display additional information about variable being read |
[in] | ind_print | individual write (possibly partial I/O) |
[in] | remove_previous_arrs | remove previous variables if present |
Definition at line 1132 of file HDF5_ops.f90.
subroutine, public hdf5_ops::print_hdf5_att | ( | type(xml_str_type), intent(inout) | att_id, |
type(xml_str_type), intent(inout) | att_dataitem, | ||
character(len=*), intent(in) | att_name, | ||
integer, intent(in) | att_center, | ||
integer, intent(in) | att_type, | ||
logical, intent(in), optional | reset, | ||
logical, intent(in), optional | ind_plot | ||
) |
Prints an HDF5 attribute.
[in,out] | att_id | ID of attribute |
[in,out] | att_dataitem | dataitem of attribute |
[in] | att_name | name of attribute |
[in] | att_center | center of attribute |
[in] | att_type | type of attribute |
[in] | reset | whether to reset dataitems |
[in] | ind_plot | whether individual plot |
Definition at line 696 of file HDF5_ops.f90.
subroutine, public hdf5_ops::print_hdf5_geom | ( | type(xml_str_type), intent(inout) | geom_id, |
integer, intent(in) | geom_type, | ||
type(xml_str_type), dimension(:), intent(inout) | geom_dataitems, | ||
logical, intent(in), optional | reset, | ||
logical, intent(in), optional | ind_plot | ||
) |
Prints an HDF5 geometry.
[in,out] | geom_id | ID of geometry |
[in] | geom_type | type of geometry |
[in,out] | geom_dataitems | data items of geometry |
[in] | reset | whether to reset dataitems |
[in] | ind_plot | whether individual plot |
Definition at line 805 of file HDF5_ops.f90.
integer function, public hdf5_ops::print_hdf5_grid | ( | type(xml_str_type), intent(inout) | grid_id, |
character(len=*), intent(in) | grid_name, | ||
integer, intent(in) | grid_type, | ||
real(dp), intent(in), optional | grid_time, | ||
type(xml_str_type), optional | grid_top, | ||
type(xml_str_type), optional | grid_geom, | ||
type(xml_str_type), dimension(:), optional | grid_atts, | ||
type(xml_str_type), dimension(:), optional | grid_grids, | ||
logical, intent(in), optional | reset, | ||
logical, intent(in), optional | ind_plot | ||
) |
Prints an HDF5 grid.
If this is is a uniform grid, the geometry and topology have to be specified, or else it will be assumed that it is already present in the XDMF domain, and reused.
If the grid is a collection grid, the grids in the collection have to be specified.
Optionally, also a time can be specified (for the grids in a collection grid).
[in,out] | grid_id | ID of grid |
[in] | grid_name | name |
[in] | grid_type | type |
[in] | grid_time | time of grid |
grid_top | topology | |
grid_geom | geometry | |
grid_atts | attributes | |
grid_grids | grids | |
[in] | reset | whether top, geom and atts and grids are reset |
[in] | ind_plot | whether individual plot |
Definition at line 886 of file HDF5_ops.f90.
subroutine, public hdf5_ops::print_hdf5_top | ( | type(xml_str_type), intent(inout) | top_id, |
integer, intent(in) | top_type, | ||
integer, dimension(:), intent(in) | top_n_elem, | ||
logical, intent(in), optional | ind_plot | ||
) |
Prints an HDF5 topology.
[in,out] | top_id | ID of topology |
[in] | top_type | type |
[in] | top_n_elem | nr. of elements |
[in] | ind_plot | whether individual plot |
Definition at line 755 of file HDF5_ops.f90.
integer function, public hdf5_ops::read_hdf5_arr | ( | type(var_1d_type), intent(inout) | var, |
character(len=*), intent(in) | PB3D_name, | ||
character(len=*), intent(in) | head_name, | ||
character(len=*), intent(in) | var_name, | ||
integer, intent(in), optional | rich_lvl, | ||
logical, intent(in), optional | disp_info, | ||
integer, dimension(:,:), intent(in), optional | lim_loc | ||
) |
Reads a PB3D output file in HDF5 format.
This happens in a non-parallel way. By default, all variables are read, but an array of strings with acceptable variable names can be passed.
Optionally, output can be given about the variable being read. Also, if rich_lvl
is provided, _R_rich_lvl
is appended to the head name if it is > 0.
Furthermore, using lim_loc
, a hyperslab of the variable can be read.
[in,out] | var | variable to read |
[in] | pb3d_name | name of PB3D file |
[in] | head_name | head name of variables |
[in] | var_name | name of variable |
[in] | rich_lvl | Richardson level to reconstruct |
[in] | disp_info | display additional information about variable being read |
[in] | lim_loc | local limits |
Definition at line 1530 of file HDF5_ops.f90.
logical, public hdf5_ops::debug_hdf5_ops = .false. |
set to true to debug general information
Definition at line 47 of file HDF5_ops.f90.
logical, public hdf5_ops::debug_print_hdf5_arrs = .false. |
set to true to debug print_HDF5_arrs
Definition at line 48 of file HDF5_ops.f90.