Neko  0.9.0
A portable framework for high-order spectral element flow simulations
chkp_file.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021-2022, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
35 module chkp_file
36  use generic_file
37  use field_series
38  use checkpoint
39  use num_types
40  use field
41  use dofmap, only: dofmap_t
42  use device_math
43  use utils
44  use space
45  use mesh
46  use math
47  use interpolation
48  use neko_mpi_types
49  use comm
51  implicit none
52  private
53 
55  type, public, extends(generic_file_t) :: chkp_file_t
56  type(space_t), pointer :: chkp_xh
57  type(space_t), pointer :: sim_xh
58  type(interpolator_t) :: space_interp
59  type(global_interpolation_t) :: global_interp
60  logical :: mesh2mesh
61  contains
62  procedure :: read => chkp_file_read
63  procedure :: read_field => chkp_read_field
64  procedure :: write => chkp_file_write
65  end type chkp_file_t
66 
67 contains
68 
70  subroutine chkp_file_write(this, data, t)
71  class(chkp_file_t), intent(inout) :: this
72  class(*), target, intent(in) :: data
73  real(kind=rp), intent(in), optional :: t
74  real(kind=dp) :: time
75  character(len=5) :: id_str
76  character(len=1024) :: fname
77  integer :: ierr, suffix_pos, optional_fields
78  type(field_t), pointer :: u, v, w, p, s
79  type(field_t), pointer :: abx1,abx2
80  type(field_t), pointer :: aby1,aby2
81  type(field_t), pointer :: abz1,abz2
82  type(field_t), pointer :: abs1,abs2
83  type(field_series_t), pointer :: ulag => null()
84  type(field_series_t), pointer :: vlag => null()
85  type(field_series_t), pointer :: wlag => null()
86  type(field_series_t), pointer :: slag => null()
87  real(kind=rp), pointer :: dtlag(:), tlag(:)
88  type(mesh_t), pointer :: msh
89  type(mpi_status) :: status
90  type(mpi_file) :: fh
91  integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
92  integer(kind=i8) :: n_glb_dofs, dof_offset
93  logical :: write_lag, write_scalar, write_dtlag, write_scalarlag, write_abvel
94  integer :: i
95 
96  if (present(t)) then
97  time = real(t,dp)
98  else
99  time = 0d0
100  end if
101 
102  select type(data)
103  type is (chkp_t)
104 
105  if ( .not. associated(data%u) .or. &
106  .not. associated(data%v) .or. &
107  .not. associated(data%w) .or. &
108  .not. associated(data%p) ) then
109  call neko_error('Checkpoint not initialized')
110  end if
111 
112  u => data%u
113  v => data%v
114  w => data%w
115  p => data%p
116  msh => u%msh
117 
118  optional_fields = 0
119 
120  if (associated(data%ulag)) then
121  ulag => data%ulag
122  vlag => data%vlag
123  wlag => data%wlag
124  write_lag = .true.
125  optional_fields = optional_fields + 1
126  else
127  write_lag = .false.
128  end if
129 
130  if (associated(data%s)) then
131  s => data%s
132  write_scalar = .true.
133  optional_fields = optional_fields + 2
134  else
135  write_scalar = .false.
136  end if
137 
138  if (associated(data%tlag)) then
139  tlag => data%tlag
140  dtlag => data%dtlag
141  write_dtlag = .true.
142  optional_fields = optional_fields + 4
143  else
144  write_dtlag = .false.
145  end if
146  write_abvel = .false.
147  if (associated(data%abx1)) then
148  abx1 => data%abx1
149  abx2 => data%abx2
150  aby1 => data%aby1
151  aby2 => data%aby2
152  abz1 => data%abz1
153  abz2 => data%abz2
154  optional_fields = optional_fields + 8
155  write_abvel = .true.
156  end if
157  write_scalarlag = .false.
158  if (associated(data%abs1)) then
159  slag => data%slag
160  abs1 => data%abs1
161  abs2 => data%abs2
162  optional_fields = optional_fields + 16
163  write_scalarlag = .true.
164  end if
165  class default
166  call neko_error('Invalid data')
167  end select
168 
169 
170  suffix_pos = filename_suffix_pos(this%fname)
171  write(id_str, '(i5.5)') this%counter
172  fname = trim(this%fname(1:suffix_pos-1))//id_str//'.chkp'
173 
174 
175  dof_offset = int(msh%offset_el, i8) * int(u%Xh%lx * u%Xh%ly * u%Xh%lz, i8)
176  n_glb_dofs = int(u%Xh%lx * u%Xh%ly * u%Xh%lz, i8) * int(msh%glb_nelv, i8)
177 
178  call mpi_file_open(neko_comm, trim(fname), &
179  mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, ierr)
180  call mpi_file_write_all(fh, msh%glb_nelv, 1, mpi_integer, status, ierr)
181  call mpi_file_write_all(fh, msh%gdim, 1, mpi_integer, status, ierr)
182  call mpi_file_write_all(fh, u%Xh%lx, 1, mpi_integer, status, ierr)
183  call mpi_file_write_all(fh, optional_fields, 1, mpi_integer, status, ierr)
184  call mpi_file_write_all(fh, time, 1, mpi_double_precision, status, ierr)
185 
186 
187  !
188  ! Dump mandatory checkpoint data
189  !
190 
191 
192  byte_offset = 4_i8 * int(mpi_integer_size,i8) + int(mpi_double_precision_size,i8)
193  byte_offset = byte_offset + &
194  dof_offset * int(mpi_real_prec_size, i8)
195  call mpi_file_write_at_all(fh, byte_offset,u%x, u%dof%size(), &
196  mpi_real_precision, status, ierr)
197  mpi_offset = 4_i8 * int(mpi_integer_size,i8) + int(mpi_double_precision_size,i8)
198  mpi_offset = mpi_offset +&
199  n_glb_dofs * int(mpi_real_prec_size, i8)
200 
201  byte_offset = mpi_offset + &
202  dof_offset * int(mpi_real_prec_size, i8)
203  call mpi_file_write_at_all(fh, byte_offset, v%x, v%dof%size(), &
204  mpi_real_precision, status, ierr)
205  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
206 
207  byte_offset = mpi_offset + &
208  dof_offset * int(mpi_real_prec_size, i8)
209  call mpi_file_write_at_all(fh, byte_offset, w%x, w%dof%size(), &
210  mpi_real_precision, status, ierr)
211  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
212 
213  byte_offset = mpi_offset + &
214  dof_offset * int(mpi_real_prec_size, i8)
215  call mpi_file_write_at_all(fh, byte_offset, p%x, p%dof%size(), &
216  mpi_real_precision, status, ierr)
217  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
218 
219  !
220  ! Dump optional payload
221  !
222 
223  if (write_lag) then
224 
225  do i = 1, ulag%size()
226  byte_offset = mpi_offset + &
227  dof_offset * int(mpi_real_prec_size, i8)
228  ! We should not need this extra associate block, ant it works
229  ! great without it for GNU, Intel, NEC and Cray, but throws an
230  ! ICE with NAG.
231  associate(x => ulag%lf(i)%x)
232  call mpi_file_write_at_all(fh, byte_offset, x, &
233  ulag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
234  end associate
235  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
236  end do
237 
238  do i = 1, vlag%size()
239  byte_offset = mpi_offset + &
240  dof_offset * int(mpi_real_prec_size, i8)
241  ! We should not need this extra associate block, ant it works
242  ! great without it for GNU, Intel, NEC and Cray, but throws an
243  ! ICE with NAG.
244  associate(x => vlag%lf(i)%x)
245  call mpi_file_write_at_all(fh, byte_offset, x, &
246  vlag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
247  end associate
248  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
249  end do
250 
251  do i = 1, wlag%size()
252  byte_offset = mpi_offset + &
253  dof_offset * int(mpi_real_prec_size, i8)
254  ! We should not need this extra associate block, ant it works
255  ! great without it for GNU, Intel, NEC and Cray, but throws an
256  ! ICE with NAG.
257  associate(x => wlag%lf(i)%x)
258  call mpi_file_write_at_all(fh, byte_offset, x, &
259  wlag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
260  end associate
261  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
262  end do
263 
264  end if
265 
266  if (write_scalar) then
267  byte_offset = mpi_offset + &
268  dof_offset * int(mpi_real_prec_size, i8)
269  call mpi_file_write_at_all(fh, byte_offset, s%x, p%dof%size(), &
270  mpi_real_precision, status, ierr)
271  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
272  end if
273 
274  if (write_dtlag) then
275  call mpi_file_write_at_all(fh, mpi_offset, tlag, 10, mpi_real_precision, status, ierr)
276  mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
277  call mpi_file_write_at_all(fh, mpi_offset, dtlag, 10, mpi_real_precision, status, ierr)
278  mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
279  end if
280 
281  if (write_abvel) then
282  byte_offset = mpi_offset + &
283  dof_offset * int(mpi_real_prec_size, i8)
284  call mpi_file_write_at_all(fh, byte_offset, abx1%x, abx1%dof%size(), &
285  mpi_real_precision, status, ierr)
286  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
287  byte_offset = mpi_offset + &
288  dof_offset * int(mpi_real_prec_size, i8)
289  call mpi_file_write_at_all(fh, byte_offset, abx2%x, abx1%dof%size(), &
290  mpi_real_precision, status, ierr)
291  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
292  byte_offset = mpi_offset + &
293  dof_offset * int(mpi_real_prec_size, i8)
294  call mpi_file_write_at_all(fh, byte_offset, aby1%x, abx1%dof%size(), &
295  mpi_real_precision, status, ierr)
296  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
297  byte_offset = mpi_offset + &
298  dof_offset * int(mpi_real_prec_size, i8)
299  call mpi_file_write_at_all(fh, byte_offset, aby2%x, abx1%dof%size(), &
300  mpi_real_precision, status, ierr)
301  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
302  byte_offset = mpi_offset + &
303  dof_offset * int(mpi_real_prec_size, i8)
304  call mpi_file_write_at_all(fh, byte_offset, abz1%x, abx1%dof%size(), &
305  mpi_real_precision, status, ierr)
306  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
307  byte_offset = mpi_offset + &
308  dof_offset * int(mpi_real_prec_size, i8)
309  call mpi_file_write_at_all(fh, byte_offset, abz2%x, abx1%dof%size(), &
310  mpi_real_precision, status, ierr)
311  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
312  end if
313 
314  if (write_scalarlag) then
315  do i = 1, slag%size()
316  byte_offset = mpi_offset + &
317  dof_offset * int(mpi_real_prec_size, i8)
318  ! We should not need this extra associate block, ant it works
319  ! great without it for GNU, Intel, NEC and Cray, but throws an
320  ! ICE with NAG.
321  associate(x => slag%lf(i)%x)
322  call mpi_file_write_at_all(fh, byte_offset, x, &
323  slag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
324  end associate
325  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
326  end do
327 
328  byte_offset = mpi_offset + &
329  dof_offset * int(mpi_real_prec_size, i8)
330  call mpi_file_write_at_all(fh, byte_offset, abs1%x, abx1%dof%size(), &
331  mpi_real_precision, status, ierr)
332  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
333  byte_offset = mpi_offset + &
334  dof_offset * int(mpi_real_prec_size, i8)
335  call mpi_file_write_at_all(fh, byte_offset, abs2%x, abx1%dof%size(), &
336  mpi_real_precision, status, ierr)
337  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
338  end if
339 
340  call mpi_file_close(fh, ierr)
341 
342  this%counter = this%counter + 1
343 
344  end subroutine chkp_file_write
345 
347  subroutine chkp_file_read(this, data)
348  class(chkp_file_t) :: this
349  class(*), target, intent(inout) :: data
350  type(chkp_t), pointer :: chkp
351  character(len=5) :: id_str
352  character(len=1024) :: fname
353  integer :: ierr, suffix_pos
354  type(field_t), pointer :: u, v, w, p, s
355  type(field_series_t), pointer :: ulag => null()
356  type(field_series_t), pointer :: vlag => null()
357  type(field_series_t), pointer :: wlag => null()
358  type(field_series_t), pointer :: slag => null()
359  type(mesh_t), pointer :: msh
360  type(mpi_status) :: status
361  type(mpi_file) :: fh
362  type(field_t), pointer :: abx1,abx2
363  type(field_t), pointer :: aby1,aby2
364  type(field_t), pointer :: abz1,abz2
365  type(field_t), pointer :: abs1,abs2
366  real(kind=rp), allocatable :: x_coord(:,:,:,:)
367  real(kind=rp), allocatable :: y_coord(:,:,:,:)
368  real(kind=rp), allocatable :: z_coord(:,:,:,:)
369  real(kind=rp), pointer :: dtlag(:), tlag(:)
370  integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
371  integer(kind=i8) :: n_glb_dofs, dof_offset
372  integer :: glb_nelv, gdim, lx, have_lag, have_scalar, nel, optional_fields, have_dtlag
373  integer :: have_abvel, have_scalarlag
374  logical :: read_lag, read_scalar, read_dtlag, read_abvel, read_scalarlag
375  real(kind=rp) :: tol
376  real(kind=rp) :: center_x, center_y, center_z
377  integer :: i, e
378  type(dofmap_t) :: dof
379 
380  call this%check_exists()
381 
382  select type(data)
383  type is (chkp_t)
384 
385  if ( .not. associated(data%u) .or. &
386  .not. associated(data%v) .or. &
387  .not. associated(data%w) .or. &
388  .not. associated(data%p) ) then
389  call neko_error('Checkpoint not initialized')
390  end if
391 
392  u => data%u
393  v => data%v
394  w => data%w
395  p => data%p
396  this%chkp_Xh => data%previous_Xh
398  if (allocated(data%previous_mesh%elements)) then
399  msh => data%previous_mesh
400  this%mesh2mesh = .true.
401  tol = data%mesh2mesh_tol
402  else
403  msh => u%msh
404  this%mesh2mesh = .false.
405  end if
406 
407  if (associated(data%ulag)) then
408  ulag => data%ulag
409  vlag => data%vlag
410  wlag => data%wlag
411  read_lag = .true.
412  else
413  read_lag = .false.
414  end if
415 
416  if (associated(data%s)) then
417  s => data%s
418  read_scalar = .true.
419  else
420  read_scalar = .false.
421  end if
422  if (associated(data%dtlag)) then
423  dtlag => data%dtlag
424  tlag => data%tlag
425  read_dtlag = .true.
426  else
427  read_dtlag = .false.
428  end if
429  read_abvel = .false.
430  if (associated(data%abx1)) then
431  abx1 => data%abx1
432  abx2 => data%abx2
433  aby1 => data%aby1
434  aby2 => data%aby2
435  abz1 => data%abz1
436  abz2 => data%abz2
437  read_abvel = .true.
438  end if
439  read_scalarlag = .false.
440  if (associated(data%abs1)) then
441  slag => data%slag
442  abs1 => data%abs1
443  abs2 => data%abs2
444  read_scalarlag = .true.
445  end if
446 
447  chkp => data
448 
449  class default
450  call neko_error('Invalid data')
451  end select
452 
453 
454 
455  call mpi_file_open(neko_comm, trim(this%fname), &
456  mpi_mode_rdonly, mpi_info_null, fh, ierr)
457  call mpi_file_read_all(fh, glb_nelv, 1, mpi_integer, status, ierr)
458  call mpi_file_read_all(fh, gdim, 1, mpi_integer, status, ierr)
459  call mpi_file_read_all(fh, lx, 1, mpi_integer, status, ierr)
460  call mpi_file_read_all(fh, optional_fields, 1, mpi_integer, status, ierr)
461  call mpi_file_read_all(fh, chkp%t, 1, mpi_double_precision, status, ierr)
462 
463  have_lag = mod(optional_fields,2)/1
464  have_scalar = mod(optional_fields,4)/2
465  have_dtlag = mod(optional_fields,8)/4
466  have_abvel = mod(optional_fields,16)/8
467  have_scalarlag = mod(optional_fields,32)/16
468 
469  if ( ( glb_nelv .ne. msh%glb_nelv ) .or. &
470  ( gdim .ne. msh%gdim) .or. &
471  ( (have_lag .eq. 0) .and. (read_lag) ) .or. &
472  ( (have_scalar .eq. 0) .and. (read_scalar) ) ) then
473  call neko_error('Checkpoint does not match case')
474  end if
475  nel = msh%nelv
476  this%sim_Xh => u%Xh
477  if (gdim .eq. 3) then
478  call this%chkp_Xh%init(gll, lx, lx, lx)
479  else
480  call this%chkp_Xh%init(gll, lx, lx)
481  end if
482  if (this%mesh2mesh) then
483  call dof%init(msh, this%chkp_Xh)
484  allocate(x_coord(u%Xh%lx,u%Xh%ly,u%Xh%lz,u%msh%nelv))
485  allocate(y_coord(u%Xh%lx,u%Xh%ly,u%Xh%lz,u%msh%nelv))
486  allocate(z_coord(u%Xh%lx,u%Xh%ly,u%Xh%lz,u%msh%nelv))
491  do e = 1, u%dof%msh%nelv
492  center_x = 0d0
493  center_y = 0d0
494  center_z = 0d0
495  do i = 1,u%dof%Xh%lxyz
496  center_x = center_x + u%dof%x(i,1,1,e)
497  center_y = center_y + u%dof%y(i,1,1,e)
498  center_z = center_z + u%dof%z(i,1,1,e)
499  end do
500  center_x = center_x/u%Xh%lxyz
501  center_y = center_y/u%Xh%lxyz
502  center_z = center_z/u%Xh%lxyz
503  do i = 1,u%dof%Xh%lxyz
504  x_coord(i,1,1,e) = u%dof%x(i,1,1,e) - tol*(u%dof%x(i,1,1,e)-center_x)
505  y_coord(i,1,1,e) = u%dof%y(i,1,1,e) - tol*(u%dof%y(i,1,1,e)-center_y)
506  z_coord(i,1,1,e) = u%dof%z(i,1,1,e) - tol*(u%dof%z(i,1,1,e)-center_z)
507  end do
508  end do
509  call this%global_interp%init(dof,tol=tol)
510  call this%global_interp%find_points(x_coord,y_coord,z_coord,u%dof%size())
511  deallocate(x_coord)
512  deallocate(y_coord)
513  deallocate(z_coord)
514  else
515  call this%space_interp%init(this%sim_Xh, this%chkp_Xh)
516  end if
517  dof_offset = int(msh%offset_el, i8) * int(this%chkp_Xh%lxyz, i8)
518  n_glb_dofs = int(this%chkp_Xh%lxyz, i8) * int(msh%glb_nelv, i8)
519 
520  !
521  ! Read mandatory checkpoint data
522  !
523 
524  byte_offset = 4_i8 * int(mpi_integer_size,i8) + int(mpi_double_precision_size,i8)
525  byte_offset = byte_offset + &
526  dof_offset * int(mpi_real_prec_size, i8)
527  call this%read_field(fh, byte_offset, u%x, nel)
528  mpi_offset = 4_i8 * int(mpi_integer_size,i8) + int(mpi_double_precision_size,i8)
529  mpi_offset = mpi_offset +&
530  n_glb_dofs * int(mpi_real_prec_size, i8)
531 
532  byte_offset = mpi_offset + &
533  dof_offset * int(mpi_real_prec_size, i8)
534  call this%read_field(fh, byte_offset, v%x, nel)
535  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
536 
537  byte_offset = mpi_offset + &
538  dof_offset * int(mpi_real_prec_size, i8)
539  call this%read_field(fh, byte_offset, w%x, nel)
540  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
541 
542  byte_offset = mpi_offset + &
543  dof_offset * int(mpi_real_prec_size, i8)
544  call this%read_field(fh, byte_offset, p%x, nel)
545  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
546 
547  !
548  ! Read optional payload
549  !
550 
551  if (read_lag) then
552  do i = 1, ulag%size()
553  byte_offset = mpi_offset + &
554  dof_offset * int(mpi_real_prec_size, i8)
555  call this%read_field(fh, byte_offset, ulag%lf(i)%x, nel)
556  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
557  end do
558 
559  do i = 1, vlag%size()
560  byte_offset = mpi_offset + &
561  dof_offset * int(mpi_real_prec_size, i8)
562  call this%read_field(fh, byte_offset, vlag%lf(i)%x, nel)
563  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
564  end do
565 
566  do i = 1, wlag%size()
567  byte_offset = mpi_offset + &
568  dof_offset * int(mpi_real_prec_size, i8)
569  call this%read_field(fh, byte_offset, wlag%lf(i)%x, nel)
570  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
571  end do
572  end if
573 
574  if (read_scalar) then
575  byte_offset = mpi_offset + &
576  dof_offset * int(mpi_real_prec_size, i8)
577  call this%read_field(fh, byte_offset, s%x, nel)
578  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
579  end if
580 
581  if (read_dtlag .and. have_dtlag .eq. 1) then
582  call mpi_file_read_at_all(fh, mpi_offset, tlag, 10, mpi_real_precision, status, ierr)
583  mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
584  call mpi_file_read_at_all(fh, mpi_offset, dtlag, 10, mpi_real_precision, status, ierr)
585  mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
586  end if
587 
588  if (read_abvel .and. have_abvel .eq. 1) then
589  byte_offset = mpi_offset + &
590  dof_offset * int(mpi_real_prec_size, i8)
591  call this%read_field(fh, byte_offset, abx1%x, nel)
592  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
593  byte_offset = mpi_offset + &
594  dof_offset * int(mpi_real_prec_size, i8)
595  call this%read_field(fh, byte_offset, abx2%x, nel)
596  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
597  byte_offset = mpi_offset + &
598  dof_offset * int(mpi_real_prec_size, i8)
599  call this%read_field(fh, byte_offset, aby1%x, nel)
600  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
601  byte_offset = mpi_offset + &
602  dof_offset * int(mpi_real_prec_size, i8)
603  call this%read_field(fh, byte_offset, aby2%x, nel)
604  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
605  byte_offset = mpi_offset + &
606  dof_offset * int(mpi_real_prec_size, i8)
607  call this%read_field(fh, byte_offset, abz1%x, nel)
608  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
609  byte_offset = mpi_offset + &
610  dof_offset * int(mpi_real_prec_size, i8)
611  call this%read_field(fh, byte_offset, abz2%x, nel)
612  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
613  end if
614  if (read_scalarlag .and. have_scalarlag .eq. 1) then
615  do i = 1, slag%size()
616  byte_offset = mpi_offset + &
617  dof_offset * int(mpi_real_prec_size, i8)
618  call this%read_field(fh, byte_offset, slag%lf(i)%x, nel)
619  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
620  end do
621  byte_offset = mpi_offset + &
622  dof_offset * int(mpi_real_prec_size, i8)
623  call this%read_field(fh, byte_offset, abs1%x, nel)
624  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
625  byte_offset = mpi_offset + &
626  dof_offset * int(mpi_real_prec_size, i8)
627  call this%read_field(fh, byte_offset, abs2%x, nel)
628  mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
629  end if
630 
631  call mpi_file_close(fh, ierr)
632 
633  call this%global_interp%free()
634  call this%space_interp%free()
635 
636  end subroutine chkp_file_read
637 
638  subroutine chkp_read_field(this, fh, byte_offset, x, nel)
639  class(chkp_file_t) :: this
640  type(mpi_file) :: fh
641  integer(kind=MPI_OFFSET_KIND) :: byte_offset
642  integer, intent(in) :: nel
643  real(kind=rp) :: x(this%sim_Xh%lxyz, nel)
644  real(kind=rp), allocatable :: read_array(:)
645  integer :: nel_stride, frac_space
646  type(mpi_status) :: status
647  integer :: ierr, i, n
648 
649  n = this%chkp_xh%lxyz*nel
650  allocate(read_array(n))
651 
652  call rzero(read_array,n)
653  call mpi_file_read_at_all(fh, byte_offset, read_array, &
654  n, mpi_real_precision, status, ierr)
655  if (this%mesh2mesh) then
656  x = 0.0_rp
657  call this%global_interp%evaluate(x,read_array)
658 
659  else if (this%sim_Xh%lx .ne. this%chkp_Xh%lx) then
660  call this%space_interp%map_host(x, read_array, nel, this%sim_Xh)
661  else
662  do i = 1,n
663  x(i,1) = read_array(i)
664  end do
665  end if
666  deallocate(read_array)
667  end subroutine chkp_read_field
668 
669 end module chkp_file
double real
Definition: device_config.h:12
Defines a checkpoint.
Definition: checkpoint.f90:34
Neko checkpoint file format.
Definition: chkp_file.f90:35
subroutine chkp_file_write(this, data, t)
Write a Neko checkpoint.
Definition: chkp_file.f90:71
subroutine chkp_file_read(this, data)
Load a checkpoint from file.
Definition: chkp_file.f90:348
subroutine chkp_read_field(this, fh, byte_offset, x, nel)
Definition: chkp_file.f90:639
Definition: comm.F90:1
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
type(mpi_datatype) mpi_real_precision
MPI type for working precision of REAL types.
Definition: comm.F90:23
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
Stores a series fields.
Defines a field.
Definition: field.f90:34
Implements global_interpolation given a dofmap.
Routines to interpolate between different spaces.
Definition: math.f90:60
Defines a mesh.
Definition: mesh.f90:34
MPI derived types.
Definition: mpi_types.f90:34
integer, public mpi_double_precision_size
Size of MPI type double precision.
Definition: mpi_types.f90:61
integer, public mpi_real_prec_size
Size of working precision REAL types.
Definition: mpi_types.f90:65
integer, public mpi_integer_size
Size of MPI type integer.
Definition: mpi_types.f90:63
integer, parameter, public i8
Definition: num_types.f90:7
integer, parameter, public dp
Definition: num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines a function space.
Definition: space.f90:34
Utilities.
Definition: utils.f90:35
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
Definition: utils.f90:56
Interface for Neko checkpoint files.
Definition: chkp_file.f90:55
A generic file handler.
Implements global interpolation for arbitrary points in the domain.
Interpolation between two space::space_t.
The function space for the SEM solution fields.
Definition: space.f90:62