PB3D  [2.45]
Ideal linear high-n MHD stability in 3-D
Interfaces and Types | Functions/Subroutines | Variables
hdf5_ops Module Reference

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...
 

Detailed Description

Operations on HDF5 and XDMF variables.

Note
Notes about XDMF:
  • Collections can be spatial or temporal. If a variable is to be evolved in time or if its domain is decomposed (e.g. the same physical variable is defined in multiple HDF5 variables), then the attribute name of this variable has to be the same for all the grids in the collection.
  • If no collection is used, but just multiple grids, the attribute can be different but does not have to be, as the different grid are distinguished by their different grid names, not by the attribute of the variables they contain.
  • The selection of hyperslabs is used, as described here: http://davis.lbl.gov/Manuals/HDF5-1.8.7/UG/12_Dataspaces.html
  • Chunking is used for partial I/O, as described here: https://www.hdfgroup.org/HDF5/doc/Advanced/Chunking/ As no chunk cache is reused, the w0 setting should be 1. Furthermore, the number of slots is chosen to be equal to the number of chunks. And Finally, the chunk size is chosen to be the size of the previous dimensions, if it does not exceed 4GB
  • Individual writes or plots are done now using the standard I/O driver, as using MPI does not allow for more than 2GB. This is important for the HELENA version of PB3D where X_2 variables get larger than that. https://www.hdfgroup.org/HDF5/doc/UG/OldHtmlSource/UG_frame08TheFile.html
  • Have a look in Documentation/XDMF_HDF5.pdf for an example.

Function/Subroutine Documentation

◆ add_hdf5_item()

integer 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.

Note
This should only be used with grids (or for topologies and geometries that are used throughout)
Returns
ierr
Parameters
[in,out]file_infoinfo about HDF5 file
[in,out]xdmf_itemXDMF item to add
[in]resetif .true., XDMF_item is reset
[in]ind_plotwhether individual plot

Definition at line 340 of file HDF5_ops.f90.

+ Here is the caller graph for this function:

◆ close_hdf5_file()

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.

Returns
ierr
Parameters
[in,out]file_infoinfo about HDF5 file
[in]ind_plotwhether individual plot
[in]cont_plotcontinued plot

Definition at line 264 of file HDF5_ops.f90.

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

◆ create_output_hdf5()

integer function, public hdf5_ops::create_output_hdf5 ( character(len=*), intent(in)  HDF5_name)

Creates an HDF5 output file.

Parameters
[in]hdf5_namename of HDF5 file

Definition at line 1075 of file HDF5_ops.f90.

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

◆ merge_hdf5_3d_data_items()

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.

Parameters
[in,out]merged_idID of merged data item
[in,out]dataitem_idsID of data items
[in]var_namename of variable
[in]dim_tottotal dimensions of variable
[in]resetwhether to reset dataitems
[in]ind_plotwhether individual plot
[in]cont_plotcontinued plot

Definition at line 603 of file HDF5_ops.f90.

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

◆ open_hdf5_file()

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.

Returns
ierr
Parameters
[in,out]file_infoinfo about HDF5 file
[in]file_namename of HDF5 file
[in,out]sym_typesymmetry type
[in]descrdescription of file
[in]ind_plotwhether individual plot
[in]cont_plotcontinued plot

Definition at line 78 of file HDF5_ops.f90.

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

◆ print_hdf5_3d_data_item()

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.

Returns
ierr
Parameters
[in,out]dataitem_idID of data item
[in]file_infoinfo about HDF5 file
[in]var_namename of variable
[in]varvariable to write
[in]dim_tottotal dimensions of variable
[in]loc_dimdimensions in this group
[in]loc_offsetoffset in this group
[in]init_valinitial fill value
[in]ind_plotwhether individual plot
[in]cont_plotcontinued plot

Definition at line 396 of file HDF5_ops.f90.

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

◆ print_hdf5_arrs()

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.

Note
See https://www.hdfgroup.org/HDF5/doc/UG/12_Dataspaces.html, 7.4.2.3 for an explanation of the selection of the dataspaces.
Returns
ierr
Parameters
[in]varsvariables to write
[in]pb3d_namename of PB3D file
[in]head_namehead name of variables
[in]rich_lvlRichardson level to append to file name
[in]disp_infodisplay additional information about variable being read
[in]ind_printindividual write (possibly partial I/O)
[in]remove_previous_arrsremove previous variables if present

Definition at line 1132 of file HDF5_ops.f90.

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

◆ print_hdf5_att()

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.

Parameters
[in,out]att_idID of attribute
[in,out]att_dataitemdataitem of attribute
[in]att_namename of attribute
[in]att_centercenter of attribute
[in]att_typetype of attribute
[in]resetwhether to reset dataitems
[in]ind_plotwhether individual plot

Definition at line 696 of file HDF5_ops.f90.

+ Here is the caller graph for this function:

◆ print_hdf5_geom()

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.

Parameters
[in,out]geom_idID of geometry
[in]geom_typetype of geometry
[in,out]geom_dataitemsdata items of geometry
[in]resetwhether to reset dataitems
[in]ind_plotwhether individual plot

Definition at line 805 of file HDF5_ops.f90.

+ Here is the caller graph for this function:

◆ print_hdf5_grid()

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).

Returns
ierr
Parameters
[in,out]grid_idID of grid
[in]grid_namename
[in]grid_typetype
[in]grid_timetime of grid
grid_toptopology
grid_geomgeometry
grid_attsattributes
grid_gridsgrids
[in]resetwhether top, geom and atts and grids are reset
[in]ind_plotwhether individual plot

Definition at line 886 of file HDF5_ops.f90.

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

◆ print_hdf5_top()

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.

Note
Currently only structured grids supported.
Parameters
[in,out]top_idID of topology
[in]top_typetype
[in]top_n_elemnr. of elements
[in]ind_plotwhether individual plot

Definition at line 755 of file HDF5_ops.f90.

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

◆ read_hdf5_arr()

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.

Note
If a limit of lim_loc is a negative value, the procedure just takes the entire range. This is necessary as sometimes the calling procuderes don't have, and don't need to have, knowledge of the underlying sizes, for example in the case of having multiple parallel jobs.
Returns
ierr
Parameters
[in,out]varvariable to read
[in]pb3d_namename of PB3D file
[in]head_namehead name of variables
[in]var_namename of variable
[in]rich_lvlRichardson level to reconstruct
[in]disp_infodisplay additional information about variable being read
[in]lim_loclocal limits

Definition at line 1530 of file HDF5_ops.f90.

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

Variable Documentation

◆ debug_hdf5_ops

logical, public hdf5_ops::debug_hdf5_ops = .false.

set to true to debug general information

Note
Debug version only

Definition at line 47 of file HDF5_ops.f90.

◆ debug_print_hdf5_arrs

logical, public hdf5_ops::debug_print_hdf5_arrs = .false.

set to true to debug print_HDF5_arrs

Note
Debug version only

Definition at line 48 of file HDF5_ops.f90.