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