Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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!
38 use checkpoint, only : chkp_t
39 use num_types, only : rp, dp, i8
40 use field, only : field_t
41 use dofmap, only: dofmap_t
43 use space, only : space_t, gll
44 use mesh, only : mesh_t
45 use math, only : rzero
50 use comm
51 use mpi_f08
52 implicit none
53 private
54
56 type, public, extends(generic_file_t) :: chkp_file_t
57 type(space_t), pointer :: chkp_xh
58 type(space_t), pointer :: sim_xh
59 type(interpolator_t) :: space_interp
60 type(global_interpolation_t) :: global_interp
61 logical :: mesh2mesh
62 logical, private :: overwrite = .false.
63 contains
64 procedure :: read => chkp_file_read
65 procedure :: read_field => chkp_read_field
66 procedure :: write => chkp_file_write
67 procedure :: set_overwrite => chkp_file_set_overwrite
68 end type chkp_file_t
69
70contains
71
73 subroutine chkp_file_write(this, data, t)
74 class(chkp_file_t), intent(inout) :: this
75 class(*), target, intent(in) :: data
76 real(kind=rp), intent(in), optional :: t
77 real(kind=dp) :: time
78 character(len=5) :: id_str
79 character(len=1024) :: fname
80 integer :: ierr, suffix_pos, optional_fields
81 type(field_t), pointer :: u, v, w, p, s
82 type(field_t), pointer :: abx1,abx2
83 type(field_t), pointer :: aby1,aby2
84 type(field_t), pointer :: abz1,abz2
85 type(field_t), pointer :: abs1,abs2
86 type(field_series_t), pointer :: ulag => null()
87 type(field_series_t), pointer :: vlag => null()
88 type(field_series_t), pointer :: wlag => null()
89 type(field_series_t), pointer :: slag => null()
90 real(kind=rp), pointer :: dtlag(:), tlag(:)
91 type(mesh_t), pointer :: msh
92 type(mpi_status) :: status
93 type(mpi_file) :: fh
94 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
95 integer(kind=i8) :: n_glb_dofs, dof_offset
96 logical :: write_lag, write_scalar, write_dtlag, write_scalarlag, write_abvel
97 integer :: i
98
99 if (present(t)) then
100 time = real(t, kind = dp)
101 else
102 time = 0.0_dp
103 end if
104
105 select type(data)
106 type is (chkp_t)
107
108 if ( .not. associated(data%u) .or. &
109 .not. associated(data%v) .or. &
110 .not. associated(data%w) .or. &
111 .not. associated(data%p) ) then
112 call neko_error('Checkpoint not initialized')
113 end if
114
115 u => data%u
116 v => data%v
117 w => data%w
118 p => data%p
119 msh => u%msh
120
121 optional_fields = 0
122
123 if (associated(data%ulag)) then
124 ulag => data%ulag
125 vlag => data%vlag
126 wlag => data%wlag
127 write_lag = .true.
128 optional_fields = optional_fields + 1
129 else
130 write_lag = .false.
131 end if
132
133 if (associated(data%s)) then
134 s => data%s
135 write_scalar = .true.
136 optional_fields = optional_fields + 2
137 else
138 write_scalar = .false.
139 end if
140
141 if (associated(data%tlag)) then
142 tlag => data%tlag
143 dtlag => data%dtlag
144 write_dtlag = .true.
145 optional_fields = optional_fields + 4
146 else
147 write_dtlag = .false.
148 end if
149 write_abvel = .false.
150 if (associated(data%abx1)) then
151 abx1 => data%abx1
152 abx2 => data%abx2
153 aby1 => data%aby1
154 aby2 => data%aby2
155 abz1 => data%abz1
156 abz2 => data%abz2
157 optional_fields = optional_fields + 8
158 write_abvel = .true.
159 end if
160 write_scalarlag = .false.
161 if (associated(data%abs1)) then
162 slag => data%slag
163 abs1 => data%abs1
164 abs2 => data%abs2
165 optional_fields = optional_fields + 16
166 write_scalarlag = .true.
167 end if
168 class default
169 call neko_error('Invalid data')
170 end select
171
172 suffix_pos = filename_suffix_pos(this%fname)
173 if (this%overwrite) then
174 fname = trim(this%fname)
175 else
176 write(id_str, '(i5.5)') this%counter
177 fname = trim(this%fname(1:suffix_pos-1)) // id_str // '.chkp'
178 end if
179
180 dof_offset = int(msh%offset_el, i8) * int(u%Xh%lx * u%Xh%ly * u%Xh%lz, i8)
181 n_glb_dofs = int(u%Xh%lx * u%Xh%ly * u%Xh%lz, i8) * int(msh%glb_nelv, i8)
182
183 call mpi_file_open(neko_comm, trim(fname), &
184 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, ierr)
185 call mpi_file_write_all(fh, msh%glb_nelv, 1, mpi_integer, status, ierr)
186 call mpi_file_write_all(fh, msh%gdim, 1, mpi_integer, status, ierr)
187 call mpi_file_write_all(fh, u%Xh%lx, 1, mpi_integer, status, ierr)
188 call mpi_file_write_all(fh, optional_fields, 1, mpi_integer, status, ierr)
189 call mpi_file_write_all(fh, time, 1, mpi_double_precision, status, ierr)
190
191
192 !
193 ! Dump mandatory checkpoint data
194 !
195
196 byte_offset = 4_i8 * int(mpi_integer_size,i8) + int(mpi_double_precision_size,i8)
197 byte_offset = byte_offset + &
198 dof_offset * int(mpi_real_prec_size, i8)
199 call mpi_file_write_at_all(fh, byte_offset,u%x, u%dof%size(), &
200 mpi_real_precision, status, ierr)
201 mpi_offset = 4_i8 * int(mpi_integer_size,i8) + int(mpi_double_precision_size,i8)
202 mpi_offset = mpi_offset +&
203 n_glb_dofs * int(mpi_real_prec_size, i8)
204
205 byte_offset = mpi_offset + &
206 dof_offset * int(mpi_real_prec_size, i8)
207 call mpi_file_write_at_all(fh, byte_offset, v%x, v%dof%size(), &
208 mpi_real_precision, status, ierr)
209 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
210
211 byte_offset = mpi_offset + &
212 dof_offset * int(mpi_real_prec_size, i8)
213 call mpi_file_write_at_all(fh, byte_offset, w%x, w%dof%size(), &
214 mpi_real_precision, status, ierr)
215 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
216
217 byte_offset = mpi_offset + &
218 dof_offset * int(mpi_real_prec_size, i8)
219 call mpi_file_write_at_all(fh, byte_offset, p%x, p%dof%size(), &
220 mpi_real_precision, status, ierr)
221 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
222
223 !
224 ! Dump optional payload
225 !
226
227 if (write_lag) then
228
229 do i = 1, ulag%size()
230 byte_offset = mpi_offset + &
231 dof_offset * int(mpi_real_prec_size, i8)
232 ! We should not need this extra associate block, ant it works
233 ! great without it for GNU, Intel, NEC and Cray, but throws an
234 ! ICE with NAG.
235 associate(x => ulag%lf(i)%x)
236 call mpi_file_write_at_all(fh, byte_offset, x, &
237 ulag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
238 end associate
239 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
240 end do
241
242 do i = 1, vlag%size()
243 byte_offset = mpi_offset + &
244 dof_offset * int(mpi_real_prec_size, i8)
245 ! We should not need this extra associate block, ant it works
246 ! great without it for GNU, Intel, NEC and Cray, but throws an
247 ! ICE with NAG.
248 associate(x => vlag%lf(i)%x)
249 call mpi_file_write_at_all(fh, byte_offset, x, &
250 vlag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
251 end associate
252 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
253 end do
254
255 do i = 1, wlag%size()
256 byte_offset = mpi_offset + &
257 dof_offset * int(mpi_real_prec_size, i8)
258 ! We should not need this extra associate block, ant it works
259 ! great without it for GNU, Intel, NEC and Cray, but throws an
260 ! ICE with NAG.
261 associate(x => wlag%lf(i)%x)
262 call mpi_file_write_at_all(fh, byte_offset, x, &
263 wlag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
264 end associate
265 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
266 end do
267
268 end if
269
270 if (write_scalar) then
271 byte_offset = mpi_offset + &
272 dof_offset * int(mpi_real_prec_size, i8)
273 call mpi_file_write_at_all(fh, byte_offset, s%x, p%dof%size(), &
274 mpi_real_precision, status, ierr)
275 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
276 end if
277
278 if (write_dtlag) then
279 call mpi_file_write_at_all(fh, mpi_offset, tlag, 10, mpi_real_precision, status, ierr)
280 mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
281 call mpi_file_write_at_all(fh, mpi_offset, dtlag, 10, mpi_real_precision, status, ierr)
282 mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
283 end if
284
285 if (write_abvel) then
286 byte_offset = mpi_offset + &
287 dof_offset * int(mpi_real_prec_size, i8)
288 call mpi_file_write_at_all(fh, byte_offset, abx1%x, abx1%dof%size(), &
289 mpi_real_precision, status, ierr)
290 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
291 byte_offset = mpi_offset + &
292 dof_offset * int(mpi_real_prec_size, i8)
293 call mpi_file_write_at_all(fh, byte_offset, abx2%x, abx1%dof%size(), &
294 mpi_real_precision, status, ierr)
295 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
296 byte_offset = mpi_offset + &
297 dof_offset * int(mpi_real_prec_size, i8)
298 call mpi_file_write_at_all(fh, byte_offset, aby1%x, abx1%dof%size(), &
299 mpi_real_precision, status, ierr)
300 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
301 byte_offset = mpi_offset + &
302 dof_offset * int(mpi_real_prec_size, i8)
303 call mpi_file_write_at_all(fh, byte_offset, aby2%x, abx1%dof%size(), &
304 mpi_real_precision, status, ierr)
305 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
306 byte_offset = mpi_offset + &
307 dof_offset * int(mpi_real_prec_size, i8)
308 call mpi_file_write_at_all(fh, byte_offset, abz1%x, abx1%dof%size(), &
309 mpi_real_precision, status, ierr)
310 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
311 byte_offset = mpi_offset + &
312 dof_offset * int(mpi_real_prec_size, i8)
313 call mpi_file_write_at_all(fh, byte_offset, abz2%x, abx1%dof%size(), &
314 mpi_real_precision, status, ierr)
315 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
316 end if
317
318 if (write_scalarlag) then
319 do i = 1, slag%size()
320 byte_offset = mpi_offset + &
321 dof_offset * int(mpi_real_prec_size, i8)
322 ! We should not need this extra associate block, ant it works
323 ! great without it for GNU, Intel, NEC and Cray, but throws an
324 ! ICE with NAG.
325 associate(x => slag%lf(i)%x)
326 call mpi_file_write_at_all(fh, byte_offset, x, &
327 slag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
328 end associate
329 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
330 end do
331
332 byte_offset = mpi_offset + &
333 dof_offset * int(mpi_real_prec_size, i8)
334 call mpi_file_write_at_all(fh, byte_offset, abs1%x, abx1%dof%size(), &
335 mpi_real_precision, status, ierr)
336 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
337 byte_offset = mpi_offset + &
338 dof_offset * int(mpi_real_prec_size, i8)
339 call mpi_file_write_at_all(fh, byte_offset, abs2%x, abx1%dof%size(), &
340 mpi_real_precision, status, ierr)
341 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
342 end if
343
344 call mpi_file_close(fh, ierr)
345
346 if (ierr .ne. mpi_success) then
347 call neko_error('Error writing checkpoint file ' // trim(fname))
348 end if
349
350 this%counter = this%counter + 1
351
352 end subroutine chkp_file_write
353
355 subroutine chkp_file_read(this, data)
356 class(chkp_file_t) :: this
357 class(*), target, intent(inout) :: data
358 type(chkp_t), pointer :: chkp
359 character(len=5) :: id_str
360 character(len=1024) :: fname
361 integer :: ierr, suffix_pos
362 type(field_t), pointer :: u, v, w, p, s
363 type(field_series_t), pointer :: ulag => null()
364 type(field_series_t), pointer :: vlag => null()
365 type(field_series_t), pointer :: wlag => null()
366 type(field_series_t), pointer :: slag => null()
367 type(mesh_t), pointer :: msh
368 type(mpi_status) :: status
369 type(mpi_file) :: fh
370 type(field_t), pointer :: abx1,abx2
371 type(field_t), pointer :: aby1,aby2
372 type(field_t), pointer :: abz1,abz2
373 type(field_t), pointer :: abs1,abs2
374 real(kind=rp), allocatable :: x_coord(:,:,:,:)
375 real(kind=rp), allocatable :: y_coord(:,:,:,:)
376 real(kind=rp), allocatable :: z_coord(:,:,:,:)
377 real(kind=rp), pointer :: dtlag(:), tlag(:)
378 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
379 integer(kind=i8) :: n_glb_dofs, dof_offset
380 integer :: glb_nelv, gdim, lx, have_lag, have_scalar, nel, optional_fields, have_dtlag
381 integer :: have_abvel, have_scalarlag
382 logical :: read_lag, read_scalar, read_dtlag, read_abvel, read_scalarlag
383 real(kind=rp) :: tol
384 real(kind=rp) :: center_x, center_y, center_z
385 integer :: i, e
386 type(dofmap_t) :: dof
387
388 call this%check_exists()
389
390 select type(data)
391 type is (chkp_t)
392
393 if ( .not. associated(data%u) .or. &
394 .not. associated(data%v) .or. &
395 .not. associated(data%w) .or. &
396 .not. associated(data%p) ) then
397 call neko_error('Checkpoint not initialized')
398 end if
399
400 u => data%u
401 v => data%v
402 w => data%w
403 p => data%p
404 this%chkp_Xh => data%previous_Xh
406 if (allocated(data%previous_mesh%elements)) then
407 msh => data%previous_mesh
408 this%mesh2mesh = .true.
409 tol = data%mesh2mesh_tol
410 else
411 msh => u%msh
412 this%mesh2mesh = .false.
413 end if
414
415 if (associated(data%ulag)) then
416 ulag => data%ulag
417 vlag => data%vlag
418 wlag => data%wlag
419 read_lag = .true.
420 else
421 read_lag = .false.
422 end if
423
424 if (associated(data%s)) then
425 s => data%s
426 read_scalar = .true.
427 else
428 read_scalar = .false.
429 end if
430 if (associated(data%dtlag)) then
431 dtlag => data%dtlag
432 tlag => data%tlag
433 read_dtlag = .true.
434 else
435 read_dtlag = .false.
436 end if
437 read_abvel = .false.
438 if (associated(data%abx1)) then
439 abx1 => data%abx1
440 abx2 => data%abx2
441 aby1 => data%aby1
442 aby2 => data%aby2
443 abz1 => data%abz1
444 abz2 => data%abz2
445 read_abvel = .true.
446 end if
447 read_scalarlag = .false.
448 if (associated(data%abs1)) then
449 slag => data%slag
450 abs1 => data%abs1
451 abs2 => data%abs2
452 read_scalarlag = .true.
453 end if
454
455 chkp => data
456
457 class default
458 call neko_error('Invalid data')
459 end select
460
461
462
463 call mpi_file_open(neko_comm, trim(this%fname), &
464 mpi_mode_rdonly, mpi_info_null, fh, ierr)
465 call mpi_file_read_all(fh, glb_nelv, 1, mpi_integer, status, ierr)
466 call mpi_file_read_all(fh, gdim, 1, mpi_integer, status, ierr)
467 call mpi_file_read_all(fh, lx, 1, mpi_integer, status, ierr)
468 call mpi_file_read_all(fh, optional_fields, 1, mpi_integer, status, ierr)
469 call mpi_file_read_all(fh, chkp%t, 1, mpi_double_precision, status, ierr)
470
471 have_lag = mod(optional_fields,2)/1
472 have_scalar = mod(optional_fields,4)/2
473 have_dtlag = mod(optional_fields,8)/4
474 have_abvel = mod(optional_fields,16)/8
475 have_scalarlag = mod(optional_fields,32)/16
476
477 if ( ( glb_nelv .ne. msh%glb_nelv ) .or. &
478 ( gdim .ne. msh%gdim) .or. &
479 ( (have_lag .eq. 0) .and. (read_lag) ) .or. &
480 ( (have_scalar .eq. 0) .and. (read_scalar) ) ) then
481 call neko_error('Checkpoint does not match case')
482 end if
483 nel = msh%nelv
484 this%sim_Xh => u%Xh
485 if (gdim .eq. 3) then
486 call this%chkp_Xh%init(gll, lx, lx, lx)
487 else
488 call this%chkp_Xh%init(gll, lx, lx)
489 end if
490 if (this%mesh2mesh) then
491 call dof%init(msh, this%chkp_Xh)
492 allocate(x_coord(u%Xh%lx,u%Xh%ly,u%Xh%lz,u%msh%nelv))
493 allocate(y_coord(u%Xh%lx,u%Xh%ly,u%Xh%lz,u%msh%nelv))
494 allocate(z_coord(u%Xh%lx,u%Xh%ly,u%Xh%lz,u%msh%nelv))
499 do e = 1, u%dof%msh%nelv
500 center_x = 0d0
501 center_y = 0d0
502 center_z = 0d0
503 do i = 1,u%dof%Xh%lxyz
504 center_x = center_x + u%dof%x(i,1,1,e)
505 center_y = center_y + u%dof%y(i,1,1,e)
506 center_z = center_z + u%dof%z(i,1,1,e)
507 end do
508 center_x = center_x/u%Xh%lxyz
509 center_y = center_y/u%Xh%lxyz
510 center_z = center_z/u%Xh%lxyz
511 do i = 1,u%dof%Xh%lxyz
512 x_coord(i,1,1,e) = u%dof%x(i,1,1,e) - tol*(u%dof%x(i,1,1,e)-center_x)
513 y_coord(i,1,1,e) = u%dof%y(i,1,1,e) - tol*(u%dof%y(i,1,1,e)-center_y)
514 z_coord(i,1,1,e) = u%dof%z(i,1,1,e) - tol*(u%dof%z(i,1,1,e)-center_z)
515 end do
516 end do
517 call this%global_interp%init(dof,tol=tol)
518 call this%global_interp%find_points(x_coord,y_coord,z_coord,u%dof%size())
519 deallocate(x_coord)
520 deallocate(y_coord)
521 deallocate(z_coord)
522 else
523 call this%space_interp%init(this%sim_Xh, this%chkp_Xh)
524 end if
525 dof_offset = int(msh%offset_el, i8) * int(this%chkp_Xh%lxyz, i8)
526 n_glb_dofs = int(this%chkp_Xh%lxyz, i8) * int(msh%glb_nelv, i8)
527
528 !
529 ! Read mandatory checkpoint data
530 !
531
532 byte_offset = 4_i8 * int(mpi_integer_size,i8) + int(mpi_double_precision_size,i8)
533 byte_offset = byte_offset + &
534 dof_offset * int(mpi_real_prec_size, i8)
535 call this%read_field(fh, byte_offset, u%x, nel)
536 mpi_offset = 4_i8 * int(mpi_integer_size,i8) + int(mpi_double_precision_size,i8)
537 mpi_offset = mpi_offset +&
538 n_glb_dofs * int(mpi_real_prec_size, i8)
539
540 byte_offset = mpi_offset + &
541 dof_offset * int(mpi_real_prec_size, i8)
542 call this%read_field(fh, byte_offset, v%x, nel)
543 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
544
545 byte_offset = mpi_offset + &
546 dof_offset * int(mpi_real_prec_size, i8)
547 call this%read_field(fh, byte_offset, w%x, nel)
548 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
549
550 byte_offset = mpi_offset + &
551 dof_offset * int(mpi_real_prec_size, i8)
552 call this%read_field(fh, byte_offset, p%x, nel)
553 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
554
555 !
556 ! Read optional payload
557 !
558
559 if (read_lag) then
560 do i = 1, ulag%size()
561 byte_offset = mpi_offset + &
562 dof_offset * int(mpi_real_prec_size, i8)
563 call this%read_field(fh, byte_offset, ulag%lf(i)%x, nel)
564 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
565 end do
566
567 do i = 1, vlag%size()
568 byte_offset = mpi_offset + &
569 dof_offset * int(mpi_real_prec_size, i8)
570 call this%read_field(fh, byte_offset, vlag%lf(i)%x, nel)
571 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
572 end do
573
574 do i = 1, wlag%size()
575 byte_offset = mpi_offset + &
576 dof_offset * int(mpi_real_prec_size, i8)
577 call this%read_field(fh, byte_offset, wlag%lf(i)%x, nel)
578 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
579 end do
580 end if
581
582 if (read_scalar) then
583 byte_offset = mpi_offset + &
584 dof_offset * int(mpi_real_prec_size, i8)
585 call this%read_field(fh, byte_offset, s%x, nel)
586 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
587 end if
588
589 if (read_dtlag .and. have_dtlag .eq. 1) then
590 call mpi_file_read_at_all(fh, mpi_offset, tlag, 10, mpi_real_precision, status, ierr)
591 mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
592 call mpi_file_read_at_all(fh, mpi_offset, dtlag, 10, mpi_real_precision, status, ierr)
593 mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
594 end if
595
596 if (read_abvel .and. have_abvel .eq. 1) then
597 byte_offset = mpi_offset + &
598 dof_offset * int(mpi_real_prec_size, i8)
599 call this%read_field(fh, byte_offset, abx1%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, abx2%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, aby1%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, aby2%x, nel)
612 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
613 byte_offset = mpi_offset + &
614 dof_offset * int(mpi_real_prec_size, i8)
615 call this%read_field(fh, byte_offset, abz1%x, nel)
616 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
617 byte_offset = mpi_offset + &
618 dof_offset * int(mpi_real_prec_size, i8)
619 call this%read_field(fh, byte_offset, abz2%x, nel)
620 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
621 end if
622 if (read_scalarlag .and. have_scalarlag .eq. 1) then
623 do i = 1, slag%size()
624 byte_offset = mpi_offset + &
625 dof_offset * int(mpi_real_prec_size, i8)
626 call this%read_field(fh, byte_offset, slag%lf(i)%x, nel)
627 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
628 end do
629 byte_offset = mpi_offset + &
630 dof_offset * int(mpi_real_prec_size, i8)
631 call this%read_field(fh, byte_offset, abs1%x, nel)
632 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
633 byte_offset = mpi_offset + &
634 dof_offset * int(mpi_real_prec_size, i8)
635 call this%read_field(fh, byte_offset, abs2%x, nel)
636 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
637 end if
638
639 call mpi_file_close(fh, ierr)
640
641 if (ierr .ne. mpi_success) then
642 call neko_error('Error reading checkpoint file ' // trim(this%fname))
643 end if
644
645 call this%global_interp%free()
646 call this%space_interp%free()
647
648 end subroutine chkp_file_read
649
650 subroutine chkp_read_field(this, fh, byte_offset, x, nel)
651 class(chkp_file_t) :: this
652 type(mpi_file) :: fh
653 integer(kind=MPI_OFFSET_KIND) :: byte_offset
654 integer, intent(in) :: nel
655 real(kind=rp) :: x(this%sim_Xh%lxyz, nel)
656 real(kind=rp), allocatable :: read_array(:)
657 integer :: nel_stride, frac_space
658 type(mpi_status) :: status
659 integer :: ierr, i, n
660
661 n = this%chkp_xh%lxyz*nel
662 allocate(read_array(n))
663
664 call rzero(read_array,n)
665 call mpi_file_read_at_all(fh, byte_offset, read_array, &
666 n, mpi_real_precision, status, ierr)
667 if (this%mesh2mesh) then
668 x = 0.0_rp
669 call this%global_interp%evaluate(x,read_array)
670
671 else if (this%sim_Xh%lx .ne. this%chkp_Xh%lx) then
672 call this%space_interp%map_host(x, read_array, nel, this%sim_Xh)
673 else
674 do i = 1,n
675 x(i,1) = read_array(i)
676 end do
677 end if
678 deallocate(read_array)
679 end subroutine chkp_read_field
680
681 subroutine chkp_file_set_overwrite(this, overwrite)
682 class(chkp_file_t), intent(inout) :: this
683 logical, intent(in) :: overwrite
684 this%overwrite = overwrite
685 end subroutine chkp_file_set_overwrite
686end module chkp_file
double real
Defines a checkpoint.
Neko checkpoint file format.
Definition chkp_file.f90:35
subroutine chkp_file_write(this, data, t)
Write a Neko checkpoint.
Definition chkp_file.f90:74
subroutine chkp_file_read(this, data)
Load a checkpoint from file.
subroutine chkp_file_set_overwrite(this, overwrite)
subroutine chkp_read_field(this, fh, byte_offset, x, nel)
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:50
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
Implements global_interpolation given a dofmap.
Routines to interpolate between different spaces.
Definition math.f90:60
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:211
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:63
integer, public mpi_real_prec_size
Size of working precision REAL types.
Definition mpi_types.f90:67
integer, public mpi_integer_size
Size of MPI type integer.
Definition mpi_types.f90:65
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
integer, parameter, public gll
Definition space.f90:48
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:57
Interface for Neko checkpoint files.
Definition chkp_file.f90:56
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
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