PB3D
[2.45]
Ideal linear high-n MHD stability in 3-D
|
Prints variables vars
with names var_names
in an HDF5 file with name c file_name and accompanying XDMF file.
More...
Public Member Functions | |
subroutine | plot_hdf5_ind (var_name, file_name, var, tot_dim, loc_offset, X, Y, Z, cont_plot, descr) |
individual version More... | |
subroutine | plot_hdf5_arr (var_names, file_name, vars, tot_dim, loc_offset, X, Y, Z, col_id, col, sym_type, cont_plot, descr) |
array version More... | |
Prints variables vars
with names var_names
in an HDF5 file with name c file_name and accompanying XDMF file.
For XDMF collections ([18]), only the first value for var_names
is used, so it should have a size of one.
The plot is generally 3-D, but if one of the dimensions provided is equal to 1, it is checked whether there is poloidal or toroidal axisymmetry and if so, the plot becomes 2-D. This can be forced using the optional input argument sym_type
.
Optionally, the (curvilinear) grid can be provided through X
, Y
and Z
. If not, the grid is assumed to be Cartesian with discrete indices where X
corresponds to the first dimensions, Y
to the second and Z
to the third.
Additionally, the total grid size and local offset can be provided in tot_dim
and loc_offset
to run this routine in parallel automatically. Optionally, one of the dimensions (col_id
, default 4) can be associated to a collection dimension using col
different from 1:
col
= 1: no collection, just plots of different variablescol
= 2: time collectioncol
= 3: spatial collectioncol
= 4: vector fieldFurthermore, using the variable cont_plot
, a plot can be (over-)written in multiple writes. By this is meant that there should be an initial plot, with collection type 1, 2, 3 or 4, which can then be followed by an arbitrary number of additional writes. As these additional writes currently cannot modify the plot structure, nor the XDMF variables, their collection dimension should be complete from the start.
This has no implications for single plots but means that for collection types 2 to 4 all the elements in the collection have to be present, though they do not necessary need to have been completely written in the other dimensions.
Subsequent writes with cont_plot
can then, for instance, write parts of the data that had not yet been written, or overwrite ones that had. This can be useful for post-processing where the memory requirements are large so that the work has to be split.
col_id
. (This could be implemented by changing how n_plot is defined and selectively letting each processer write in the main loop at its corresponding indices.)Definition at line 137 of file output_ops.f90.
subroutine output_ops::plot_hdf5::plot_hdf5_arr | ( | character(len=*), dimension(:), intent(in) | var_names, |
character(len=*), intent(in) | file_name, | ||
real(dp), dimension(:,:,:,:), intent(in), target | vars, | ||
integer, dimension(4), intent(in), optional | tot_dim, | ||
integer, dimension(4), intent(in), optional | loc_offset, | ||
real(dp), dimension(:,:,:,:), intent(in), optional, target | X, | ||
real(dp), dimension(:,:,:,:), intent(in), optional, target | Y, | ||
real(dp), dimension(:,:,:,:), intent(in), optional, target | Z, | ||
integer, intent(in), optional | col_id, | ||
integer, intent(in), optional | col, | ||
integer, intent(in), optional | sym_type, | ||
logical, intent(in), optional | cont_plot, | ||
character(len=*), intent(in), optional | descr | ||
) |
array version
[in] | var_names | names of variable to be plot |
[in] | file_name | file name |
[in] | vars | variables to plot |
[in] | tot_dim | total dimensions of the arrays |
[in] | loc_offset | offset of local dimensions |
[in] | x | curvilinear grid X points |
[in] | y | curvilinear grid Y points |
[in] | z | curvilinear grid Z points |
[in] | col_id | index of time dimension |
[in] | col | whether a collection is made |
[in] | sym_type | type of symmetry (1: no symmetry, 2: toroidal, 3: poloidal) |
[in] | cont_plot | continued plot |
[in] | descr | description |
Definition at line 435 of file output_ops.f90.
subroutine output_ops::plot_hdf5::plot_hdf5_ind | ( | character(len=*), intent(in) | var_name, |
character(len=*), intent(in) | file_name, | ||
real(dp), dimension(:,:,:), intent(in) | var, | ||
integer, dimension(3), intent(in), optional | tot_dim, | ||
integer, dimension(3), intent(in), optional | loc_offset, | ||
real(dp), dimension(:,:,:), intent(in), optional, target | X, | ||
real(dp), dimension(:,:,:), intent(in), optional, target | Y, | ||
real(dp), dimension(:,:,:), intent(in), optional, target | Z, | ||
logical, intent(in), optional | cont_plot, | ||
character(len=*), intent(in), optional | descr | ||
) |
individual version
[in] | var_name | name of variable to be plot |
[in] | file_name | file name |
[in] | var | variable to plot |
[in] | tot_dim | total dimensions of the arrays |
[in] | loc_offset | offset of local dimensions |
[in] | x | curvilinear grid X points |
[in] | y | curvilinear grid Y points |
[in] | z | curvilinear grid Z points |
[in] | cont_plot | continued plot |
[in] | descr | description |
Definition at line 949 of file output_ops.f90.