MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_input_output_helper.t
Go to the documentation of this file.
2 
3  implicit none
4 
5  ! public methods
6  public :: count_ix
7  public :: create_output_file
8  public :: snapshot_write_header1
9  public :: block_shape_io
10  public :: get_names_from_string
11 
12  !> whether a manually inserted snapshot is saved
13  logical :: save_now
14  !> Version number of the .dat file output
15  integer, parameter :: version_number = 5
16  contains
17 
18  function get_names_from_string(aux_variable_names,nwc) result(names)
20  character(len=*),intent(in):: aux_variable_names
21  integer, intent(in) :: nwc
22  character(len=name_len) :: names(1:nwc)
23 
24  integer:: space_position,iw
25  character(len=name_len):: wname
26  character(len=std_len):: scanstring
27 
28  ! copied from subroutine getheadernames in calculate_xw
29  !TODO check for errors
30  scanstring=trim(adjustl(aux_variable_names))
31  space_position=index(scanstring,' ')
32  do iw=1,nwc
33  do while (space_position==1)
34  scanstring=scanstring(2:)
35  space_position=index(scanstring,' ')
36  enddo
37  wname=scanstring(:space_position-1)
38  scanstring=scanstring(space_position+1:)
39  space_position=index(scanstring,' ')
40 
41  names(iw)=trim(adjustl(wname))
42  enddo
43  end function get_names_from_string
44 
45  !> Compute number of elements in index range
46  pure integer function count_ix(ixO^L)
47  integer, intent(in) :: ixo^l
48 
49  count_ix = product([ ixomax^d ] - [ ixomin^d ] + 1)
50  end function count_ix
51 
52  !> Determine the shape of a block for output (whether to include ghost cells,
53  !> and on which sides)
54  subroutine block_shape_io(igrid, n_ghost, ixO^L, n_values)
56 
57  integer, intent(in) :: igrid
58  !> nghost(1:ndim) contains the number of ghost cells on the block's minimum
59  !> boundaries, and nghost(ndim+1:2*ndim) on the block's maximum boundaries
60  integer, intent(out) :: n_ghost(2*ndim)
61  !> Index range on output block
62  integer, intent(out) :: ixo^l
63  !> Number of cells/values in output
64  integer, intent(out) :: n_values
65 
66  integer :: idim
67 
68  n_ghost(:) = 0
69 
70  if(save_physical_boundary) then
71  do idim=1,ndim
72  ! Include ghost cells on lower boundary
73  if(ps(igrid)%is_physical_boundary(2*idim-1)) n_ghost(idim)=nghostcells
74  ! Include ghost cells on upper boundary
75  if(ps(igrid)%is_physical_boundary(2*idim)) n_ghost(ndim+idim)=nghostcells
76  end do
77  end if
78 
79  {ixomin^d = ixmlo^d - n_ghost(^d)\}
80  {ixomax^d = ixmhi^d + n_ghost(ndim+^d)\}
81 
82  n_values = count_ix(ixo^l) * nw
83 
84  end subroutine block_shape_io
85 
86  !> Standard method for creating a new output file
87  subroutine create_output_file(fh, ix, extension, suffix)
89  use mod_comm_lib, only: mpistop
90  integer, intent(out) :: fh !< File handle
91  integer, intent(in) :: ix !< Index of file
92  character(len=*), intent(in) :: extension !< Extension of file
93  character(len=*), intent(in), optional :: suffix !< Optional suffix
94 
95  integer :: amode
96  character(len=std_len) :: filename
97 
98  if (ix >= 10000) then
99  call mpistop("Number of output files is limited to 10000 (0...9999)")
100  end if
101 
102  if (present(suffix)) then
103  write(filename,"(a,a,i4.4,a)") trim(base_filename), &
104  trim(suffix), ix, extension
105  else
106  write(filename,"(a,i4.4,a)") trim(base_filename), ix, extension
107  end if
108 
109  ! MPI cannot easily replace existing files
110  open(unit=unitsnapshot,file=filename,status='replace')
111  close(unitsnapshot, status='delete')
112 
113  amode = ior(mpi_mode_create, mpi_mode_wronly)
114  call mpi_file_open(mpi_comm_self,filename,amode, &
115  mpi_info_null, fh, ierrmpi)
116 
117  if (ierrmpi /= 0) then
118  print *, "Error, cannot create file ", trim(filename)
119  call mpistop("Fatal error")
120  end if
121 
122  end subroutine create_output_file
123 
124  !> Write header for a snapshot, generalize cons_wnames and nw
125  !>
126  !> If you edit the header, don't forget to update: snapshot_write_header(),
127  !> snapshot_read_header(), doc/fileformat.md, tools/python/dat_reader.py
128  subroutine snapshot_write_header1(fh, offset_tree, offset_block, dataset_names, nw_vars)
129  use mod_forest
130  use mod_physics
132  use mod_slice, only: slicenext
133  integer, intent(in) :: fh !< File handle
134  integer(kind=MPI_OFFSET_KIND), intent(in) :: offset_tree !< Offset of tree info
135  integer(kind=MPI_OFFSET_KIND), intent(in) :: offset_block !< Offset of block data
136  character(len=*), intent(in) :: dataset_names(:)
137  integer, intent(in) :: nw_vars
138  integer, dimension(MPI_STATUS_SIZE) :: st
139  integer :: iw, er
140 
141  character(len=name_len) :: dname
142 
143  call mpi_file_write(fh, version_number, 1, mpi_integer, st, er)
144  call mpi_file_write(fh, int(offset_tree), 1, mpi_integer, st, er)
145  call mpi_file_write(fh, int(offset_block), 1, mpi_integer, st, er)
146  call mpi_file_write(fh, nw_vars, 1, mpi_integer, st, er)
147  call mpi_file_write(fh, ndir, 1, mpi_integer, st, er)
148  call mpi_file_write(fh, ndim, 1, mpi_integer, st, er)
149  call mpi_file_write(fh, levmax, 1, mpi_integer, st, er)
150  call mpi_file_write(fh, nleafs, 1, mpi_integer, st, er)
151  call mpi_file_write(fh, nparents, 1, mpi_integer, st, er)
152  call mpi_file_write(fh, it, 1, mpi_integer, st, er)
153  ! Note: It is nice when this double has an even number of 4 byte
154  ! integers before it (for alignment)
155  call mpi_file_write(fh, global_time, 1, mpi_double_precision, st, er)
156 
157  call mpi_file_write(fh, [ xprobmin^d ], ndim, &
158  mpi_double_precision, st, er)
159  call mpi_file_write(fh, [ xprobmax^d ], ndim, &
160  mpi_double_precision, st, er)
161  call mpi_file_write(fh, [ domain_nx^d ], ndim, &
162  mpi_integer, st, er)
163  call mpi_file_write(fh, [ block_nx^d ], ndim, &
164  mpi_integer, st, er)
165 
166  ! Periodicity (assume all variables are periodic if one is)
167  call mpi_file_write(fh, periodb, ndim, mpi_logical, st, er)
168 
169  ! Geometry
170  call mpi_file_write(fh, geometry_name(1:name_len), &
171  name_len, mpi_character, st, er)
172 
173  ! Write stagger grid mark
174  call mpi_file_write(fh, stagger_grid, 1, mpi_logical, st, er)
175 
176  do iw = 1, nw_vars
177  ! using directly trim(adjustl((dataset_names(iw)))) in MPI_FILE_WRITE call
178  ! does not work, there will be trailing characters
179  dname = trim(adjustl((dataset_names(iw))))
180  call mpi_file_write(fh, dname, name_len, mpi_character, st, er)
181  end do
182 
183  ! Physics related information
184  call mpi_file_write(fh, physics_type, name_len, mpi_character, st, er)
185 
186  ! Format:
187  ! integer :: n_par
188  ! double precision :: values(n_par)
189  ! character(n_par * name_len) :: names
190  call phys_write_info(fh)
191 
192  ! Write snapshotnext etc., which is useful for restarting.
193  ! Note we add one, since snapshotnext is updated *after* this procedure
194  if(pass_wall_time.or.save_now) then
195  call mpi_file_write(fh, snapshotnext, 1, mpi_integer, st, er)
196  else
197  call mpi_file_write(fh, snapshotnext+1, 1, mpi_integer, st, er)
198  end if
199  call mpi_file_write(fh, slicenext, 1, mpi_integer, st, er)
200  call mpi_file_write(fh, collapsenext, 1, mpi_integer, st, er)
201 
202  end subroutine snapshot_write_header1
203 
204 end module mod_input_output_helper
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
Module with basic grid data structures.
Definition: mod_forest.t:2
integer, save nleafs
Number of leaf block.
Definition: mod_forest.t:76
integer, save nparents
Number of parent blocks.
Definition: mod_forest.t:73
This module contains definitions of global parameters and variables and some generic functions/subrou...
character(len=std_len) geometry_name
integer domain_nx
number of cells for each dimension in level-one mesh
double precision global_time
The global simulation time.
integer it
Number of time steps taken.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical save_physical_boundary
True for save physical boundary cells in dat files.
logical stagger_grid
True for using stagger grid.
integer block_nx
number of cells for each dimension in grid block excluding ghostcells
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
integer ierrmpi
A global MPI error return code.
integer snapshotnext
IO: snapshot and collapsed views output numbers/labels.
logical, dimension(ndim) periodb
True for dimensions with periodic boundaries.
integer, parameter unitsnapshot
integer nghostcells
Number of ghost cells surrounding a grid.
character(len=std_len) base_filename
Base file name for simulation output, which will be followed by a number.
logical pass_wall_time
If true, wall time is up, modify snapshotnext for later overwrite.
subroutine, public create_output_file(fh, ix, extension, suffix)
Standard method for creating a new output file.
pure integer function, public count_ix(ixOL)
Compute number of elements in index range.
subroutine, public snapshot_write_header1(fh, offset_tree, offset_block, dataset_names, nw_vars)
Write header for a snapshot, generalize cons_wnames and nw.
character(len=name_len) function, dimension(1:nwc), public get_names_from_string(aux_variable_names, nwc)
subroutine, public block_shape_io(igrid, n_ghost, ixOL, n_values)
Determine the shape of a block for output (whether to include ghost cells, and on which sides)
logical save_now
whether a manually inserted snapshot is saved
integer, parameter version_number
Version number of the .dat file output.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
procedure(sub_write_info), pointer phys_write_info
Definition: mod_physics.t:81
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition: mod_physics.t:54
Writes D-1 slice, can do so in various formats, depending on slice_type.
Definition: mod_slice.t:2
integer slicenext
the file number of slices
Definition: mod_slice.t:14