Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_hdf5_serial.F90
Go to the documentation of this file.
1 !------------------------------------------------------------------------------
2 ! SELALIB
3 !------------------------------------------------------------------------------
4 ! MODULE: sll_m_hdf5_serial
5 !
6 ! DESCRIPTION:
11 !------------------------------------------------------------------------------
13 
14 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
15 #ifndef NOHDF5
16 
17 #include "sll_errors.h"
18 #include "sll_working_precision.h"
19 
20  use sll_m_hdf5_io_serial, only: &
26 
27  implicit none
28 
29  public :: &
31 
32  private
33 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
34 
36  integer, parameter :: maxlen = 256
37 
38 !==============================================================================
39 
41 
42  character(len=maxlen) :: filename = ""
43  logical :: is_open = .false.
44 
45  type(sll_t_hdf5_ser_handle), private :: file_handle
46 
47  contains
48  procedure :: init => t_hdf5_serial__init
49  procedure :: create => t_hdf5_serial__create
50  procedure :: open => t_hdf5_serial__open
51  procedure :: close => t_hdf5_serial__close
52  procedure :: delete => t_hdf5_serial__delete
53 
54  procedure :: write_dble_array_1d => t_hdf5_serial__write_dble_array_1d
55  procedure :: write_dble_array_2d => t_hdf5_serial__write_dble_array_2d
56  procedure :: write_dble_array_3d => t_hdf5_serial__write_dble_array_3d
57  procedure :: write_int_array_1d => t_hdf5_serial__write_int_array_1d
58  procedure :: write_int_array_2d => t_hdf5_serial__write_int_array_2d
59  procedure :: write_int_array_3d => t_hdf5_serial__write_int_array_3d
60  generic :: write_array => write_dble_array_1d, &
61  write_dble_array_2d, &
62  write_dble_array_3d, &
63  write_int_array_1d, &
64  write_int_array_2d, &
65  write_int_array_3d
66  end type sll_t_hdf5_serial
67 
68 !==============================================================================
69 contains
70 !==============================================================================
71 
72  subroutine t_hdf5_serial__init(self, filename)
73  class(sll_t_hdf5_serial), intent(out) :: self
74  character(len=*), intent(in) :: filename
75 
76  self%filename = filename
77 
78  end subroutine t_hdf5_serial__init
79 
80  !----------------------------------------------------------------------------
81  subroutine t_hdf5_serial__create(self)
82  class(sll_t_hdf5_serial), intent(inout) :: self
83 
84  character(len=*), parameter :: this_sub_name = "t_hdf5_serial__create"
85  integer :: error
86 
87  call sll_s_hdf5_ser_file_create(self%filename, self%file_handle, error)
88 
89  if (error /= 0) then
90  sll_error(this_sub_name, "Cannot create HDF5 file "//trim(self%filename))
91  end if
92 
93  self%is_open = .true.
94 
95  end subroutine t_hdf5_serial__create
96 
97  !----------------------------------------------------------------------------
98  subroutine t_hdf5_serial__open(self)
99  class(sll_t_hdf5_serial), intent(inout) :: self
100 
101  character(len=*), parameter :: this_sub_name = "t_hdf5_serial__open"
102  integer :: error
103 
104  if (self%is_open) then
105  sll_error(this_sub_name, "File is already open")
106  end if
107 
108  call sll_s_hdf5_ser_file_open(self%filename, self%file_handle, error)
109 
110  if (error /= 0) then
111  sll_error(this_sub_name, "Cannot open HDF5 file "//trim(self%filename))
112  end if
113 
114  self%is_open = .true.
115 
116  end subroutine t_hdf5_serial__open
117 
118  !----------------------------------------------------------------------------
119  subroutine t_hdf5_serial__close(self)
120  class(sll_t_hdf5_serial), intent(inout) :: self
121 
122  character(len=*), parameter :: this_sub_name = "t_hdf5_serial__close"
123  !character(len=:), allocatable :: err_msg
124  integer :: error
125 
126  if (.not. self%is_open) then
127  sll_error(this_sub_name, "File is already closed")
128  end if
129 
130  call sll_s_hdf5_ser_file_close(self%file_handle, error)
131 
132  if (error /= 0) then
133  sll_error(this_sub_name, "Cannot close HDF5 file "//trim(self%filename))
134  end if
135 
136  self%is_open = .false.
137 
138  end subroutine t_hdf5_serial__close
139 
140  !----------------------------------------------------------------------------
141  subroutine t_hdf5_serial__delete(self)
142  class(sll_t_hdf5_serial), intent(inout) :: self
143 
144  if (self%is_open) call self%close()
145 
146  self%filename = ""
147  self%is_open = .false.
148 
149  end subroutine t_hdf5_serial__delete
150 
151  !----------------------------------------------------------------------------
152  subroutine t_hdf5_serial__write_dble_array_1d(self, array, dsetname)
153  class(sll_t_hdf5_serial), intent(inout) :: self
154  sll_real64, intent(in) :: array(:)
155  character(len=*), intent(in) :: dsetname
156 
157  character(len=*), parameter :: this_sub_name = "t_hdf5_serial__write_dble_array_1d"
158  !character(len=:), allocatable :: err_msg
159  integer :: error
160 
161  if (.not. self%is_open) then
162  sll_error(this_sub_name, "File is closed")
163  end if
164 
165  call sll_o_hdf5_ser_write_array(self%file_handle, array, dsetname, error)
166 
167  if (error /= 0) then
168  sll_error(this_sub_name, "Cannot write to file "//trim(self%filename))
169  end if
170 
172 
173  !----------------------------------------------------------------------------
174  subroutine t_hdf5_serial__write_dble_array_2d(self, array, dsetname)
175  class(sll_t_hdf5_serial), intent(inout) :: self
176  sll_real64, intent(in) :: array(:, :)
177  character(len=*), intent(in) :: dsetname
178 
179  character(len=*), parameter :: this_sub_name = "t_hdf5_serial__write_dble_array_2d"
180  integer :: error
181 
182  if (.not. self%is_open) then
183  sll_error(this_sub_name, "File is closed")
184  end if
185 
186  call sll_o_hdf5_ser_write_array(self%file_handle, array, dsetname, error)
187 
188  if (error /= 0) then
189  sll_error(this_sub_name, "Cannot write to file "//trim(self%filename))
190  end if
191 
193 
194  !----------------------------------------------------------------------------
195  subroutine t_hdf5_serial__write_dble_array_3d(self, array, dsetname)
196  class(sll_t_hdf5_serial), intent(inout) :: self
197  sll_real64, intent(in) :: array(:, :, :)
198  character(len=*), intent(in) :: dsetname
199 
200  character(len=*), parameter :: this_sub_name = "t_hdf5_serial__write_dble_array_3d"
201  !character(len=:), allocatable :: err_msg
202  integer :: error
203 
204  if (.not. self%is_open) then
205  sll_error(this_sub_name, "File is closed")
206  end if
207 
208  call sll_o_hdf5_ser_write_array(self%file_handle, array, dsetname, error)
209 
210  if (error /= 0) then
211  sll_error(this_sub_name, "Cannot write to file "//trim(self%filename))
212  end if
213 
215 
216  !----------------------------------------------------------------------------
217  subroutine t_hdf5_serial__write_int_array_1d(self, array, dsetname)
218  class(sll_t_hdf5_serial), intent(inout) :: self
219  sll_int32, intent(in) :: array(:)
220  character(len=*), intent(in) :: dsetname
221 
222  character(len=*), parameter :: this_sub_name = "t_hdf5_serial__write_int_array_1d"
223  !character(len=:), allocatable :: err_msg
224  integer :: error
225 
226  if (.not. self%is_open) then
227  sll_error(this_sub_name, "File is closed")
228  end if
229 
230  call sll_o_hdf5_ser_write_array(self%file_handle, array, dsetname, error)
231 
232  if (error /= 0) then
233  sll_error(this_sub_name, "Cannot write to file "//trim(self%filename))
234  end if
235 
236  end subroutine t_hdf5_serial__write_int_array_1d
237 
238  !----------------------------------------------------------------------------
239  subroutine t_hdf5_serial__write_int_array_2d(self, array, dsetname)
240  class(sll_t_hdf5_serial), intent(inout) :: self
241  sll_int32, intent(in) :: array(:, :)
242  character(len=*), intent(in) :: dsetname
243 
244  character(len=*), parameter :: this_sub_name = "t_hdf5_serial__write_int_array_2d"
245  !character(len=:), allocatable :: err_msg
246  integer :: error
247 
248  if (.not. self%is_open) then
249  sll_error(this_sub_name, "File is closed")
250  end if
251 
252  call sll_o_hdf5_ser_write_array(self%file_handle, array, dsetname, error)
253 
254  if (error /= 0) then
255  sll_error(this_sub_name, "Cannot write to file "//trim(self%filename))
256  end if
257 
258  end subroutine t_hdf5_serial__write_int_array_2d
259 
260  !----------------------------------------------------------------------------
261  subroutine t_hdf5_serial__write_int_array_3d(self, array, dsetname)
262  class(sll_t_hdf5_serial), intent(inout) :: self
263  sll_int32, intent(in) :: array(:, :, :)
264  character(len=*), intent(in) :: dsetname
265 
266  character(len=*), parameter :: this_sub_name = "t_hdf5_serial__write_int_array_3d"
267  !character(len=:), allocatable :: err_msg
268  integer :: error
269 
270  if (.not. self%is_open) then
271  sll_error(this_sub_name, "File is closed")
272  end if
273 
274  call sll_o_hdf5_ser_write_array(self%file_handle, array, dsetname, error)
275 
276  if (error /= 0) then
277  sll_error(this_sub_name, "Cannot write to file "//trim(self%filename))
278  end if
279 
280  end subroutine t_hdf5_serial__write_int_array_3d
281 
282 #endif /* NOHDF5 */
283 
284 end module sll_m_hdf5_serial
Write nD array of double precision floats or integers into HDF5 file.
Implements the functions to write hdf5 file to store heavy data.
subroutine, public sll_s_hdf5_ser_file_create(filename, handle, error)
Create new HDF5 file.
subroutine, public sll_s_hdf5_ser_file_close(handle, error)
Close existing HDF5 file.
subroutine, public sll_s_hdf5_ser_file_open(filename, handle, error)
Open existing HDF5 file.
Simple object-oriented wrapper to Pierre's "sll_m_hdf5_io_serial".
subroutine t_hdf5_serial__create(self)
subroutine t_hdf5_serial__close(self)
subroutine t_hdf5_serial__write_int_array_1d(self, array, dsetname)
subroutine t_hdf5_serial__write_int_array_2d(self, array, dsetname)
subroutine t_hdf5_serial__open(self)
subroutine t_hdf5_serial__delete(self)
subroutine t_hdf5_serial__init(self, filename)
subroutine t_hdf5_serial__write_dble_array_3d(self, array, dsetname)
subroutine t_hdf5_serial__write_int_array_3d(self, array, dsetname)
integer, parameter maxlen
Maximum length of non-allocatable strings in module.
subroutine t_hdf5_serial__write_dble_array_1d(self, array, dsetname)
subroutine t_hdf5_serial__write_dble_array_2d(self, array, dsetname)
    Report Typos and Errors