Neko 1.99.3
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 :: wm_x => null()
92 type(field_t), pointer :: wm_y => null()
93 type(field_t), pointer :: wm_z => null()
94 type(field_t), pointer :: abx1, abx2
95 type(field_t), pointer :: aby1, aby2
96 type(field_t), pointer :: abz1, abz2
97 type(field_t), pointer :: abs1, abs2
98 type(field_series_t), pointer :: ulag => null()
99 type(field_series_t), pointer :: vlag => null()
100 type(field_series_t), pointer :: wlag => null()
101 type(field_series_t), pointer :: wm_x_lag => null()
102 type(field_series_t), pointer :: wm_y_lag => null()
103 type(field_series_t), pointer :: wm_z_lag => null()
104 type(field_series_t), pointer :: slag => null()
105 real(kind=rp), pointer :: msh_x(:,:,:,:) => null()
106 real(kind=rp), pointer :: msh_y(:,:,:,:) => null()
107 real(kind=rp), pointer :: msh_z(:,:,:,:) => null()
108 real(kind=rp), pointer :: blag(:,:,:,:) => null()
109 real(kind=rp), pointer :: blaglag(:,:,:,:) => null()
110 real(kind=rp), pointer :: pivot_pos(:) => null()
111 real(kind=rp), pointer :: pivot_vel_lag(:,:) => null()
112 real(kind=rp), pointer :: basis_pos(:) => null()
113 real(kind=rp), pointer :: basis_vel_lag(:,:) => null()
114 real(kind=rp), pointer :: dtlag(:), tlag(:)
115 type(mesh_t), pointer :: msh
116 type(mpi_status) :: status
117 type(mpi_file) :: fh
118 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
119 integer(kind=i8) :: n_glb_dofs, dof_offset
120 logical :: write_lag, write_scalar, write_dtlag
121 logical :: write_ale
122 logical :: write_scalarlag, write_abvel
123 integer :: i
124
125 if (present(t)) then
126 time = real(t, kind = dp)
127 else
128 time = 0.0_dp
129 end if
130
131 select type (data)
132 type is (chkp_t)
133
134 if ( .not. associated(data%u) .or. &
135 .not. associated(data%v) .or. &
136 .not. associated(data%w) .or. &
137 .not. associated(data%p) ) then
138 call neko_error('Checkpoint not initialized')
139 end if
140
141 u => data%u
142 v => data%v
143 w => data%w
144 p => data%p
145 msh => u%msh
146
147 optional_fields = 0
148
149 if (associated(data%ulag)) then
150 ulag => data%ulag
151 vlag => data%vlag
152 wlag => data%wlag
153 write_lag = .true.
154 optional_fields = optional_fields + 1
155 else
156 write_lag = .false.
157 end if
158
159 if (associated(data%s)) then
160 s => data%s
161 write_scalar = .true.
162 optional_fields = optional_fields + 2
163 else
164 write_scalar = .false.
165 end if
166
167 if (associated(data%tlag)) then
168 tlag => data%tlag
169 dtlag => data%dtlag
170 write_dtlag = .true.
171 optional_fields = optional_fields + 4
172 else
173 write_dtlag = .false.
174 end if
175
176 write_abvel = .false.
177 if (associated(data%abx1)) then
178 abx1 => data%abx1
179 abx2 => data%abx2
180 aby1 => data%aby1
181 aby2 => data%aby2
182 abz1 => data%abz1
183 abz2 => data%abz2
184 optional_fields = optional_fields + 8
185 write_abvel = .true.
186 end if
187
188 write_scalarlag = .false.
189 if (associated(data%abs1)) then
190 slag => data%slag
191 abs1 => data%abs1
192 abs2 => data%abs2
193 optional_fields = optional_fields + 16
194 write_scalarlag = .true.
195 end if
196
197 write_ale = .false.
198 if (associated(data%wm_x)) then
199 msh_x => data%msh_x
200 msh_y => data%msh_y
201 msh_z => data%msh_z
202 wm_x => data%wm_x
203 wm_y => data%wm_y
204 wm_z => data%wm_z
205 wm_x_lag => data%wm_x_lag
206 wm_y_lag => data%wm_y_lag
207 wm_z_lag => data%wm_z_lag
208 blag => data%Blag
209 blaglag => data%Blaglag
210 pivot_pos => data%pivot_pos
211 pivot_vel_lag => data%pivot_vel_lag
212 basis_pos => data%basis_pos
213 basis_vel_lag => data%basis_vel_lag
214 optional_fields = optional_fields + 32
215 write_ale = .true.
216 end if
217
218 class default
219 call neko_error('Invalid data')
220 end select
221
222 dof_offset = int(msh%offset_el, i8) * int(u%Xh%lx * u%Xh%ly * u%Xh%lz, i8)
223 n_glb_dofs = int(u%Xh%lx * u%Xh%ly * u%Xh%lz, i8) * int(msh%glb_nelv, i8)
224
225 ! Retrieve the filename and increment the counter if we are not overwriting
226 if (.not. this%overwrite) call this%increment_counter()
227 fname = trim(this%get_fname())
228
229 call mpi_file_open(neko_comm, trim(fname), &
230 mpi_mode_wronly + mpi_mode_create, mpi_info_null, fh, ierr)
231 call mpi_file_write_all(fh, msh%glb_nelv, 1, mpi_integer, status, ierr)
232 call mpi_file_write_all(fh, msh%gdim, 1, mpi_integer, status, ierr)
233 call mpi_file_write_all(fh, u%Xh%lx, 1, mpi_integer, status, ierr)
234 call mpi_file_write_all(fh, optional_fields, 1, mpi_integer, status, ierr)
235 call mpi_file_write_all(fh, time, 1, mpi_double_precision, status, ierr)
236
237
238 !
239 ! Dump mandatory checkpoint data
240 !
241
242 byte_offset = 4_i8 * int(mpi_integer_size, i8) + &
244 byte_offset = byte_offset + &
245 dof_offset * int(mpi_real_prec_size, i8)
246 call mpi_file_write_at_all(fh, byte_offset,u%x, u%dof%size(), &
247 mpi_real_precision, status, ierr)
248 mpi_offset = 4_i8 * int(mpi_integer_size, i8) + &
250 mpi_offset = mpi_offset +&
251 n_glb_dofs * int(mpi_real_prec_size, i8)
252
253 byte_offset = mpi_offset + &
254 dof_offset * int(mpi_real_prec_size, i8)
255 call mpi_file_write_at_all(fh, byte_offset, v%x, v%dof%size(), &
256 mpi_real_precision, status, ierr)
257 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
258
259 byte_offset = mpi_offset + &
260 dof_offset * int(mpi_real_prec_size, i8)
261 call mpi_file_write_at_all(fh, byte_offset, w%x, w%dof%size(), &
262 mpi_real_precision, status, ierr)
263 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
264
265 byte_offset = mpi_offset + &
266 dof_offset * int(mpi_real_prec_size, i8)
267 call mpi_file_write_at_all(fh, byte_offset, p%x, p%dof%size(), &
268 mpi_real_precision, status, ierr)
269 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
270
271 !
272 ! Dump optional payload
273 !
274
275 if (write_lag) then
276
277 do i = 1, ulag%size()
278 byte_offset = mpi_offset + &
279 dof_offset * int(mpi_real_prec_size, i8)
280 ! We should not need this extra associate block, and it works
281 ! great without it for GNU, Intel, NEC and Cray, but throws an
282 ! ICE with NAG.
283 associate(x => ulag%lf(i)%x)
284 call mpi_file_write_at_all(fh, byte_offset, x, &
285 ulag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
286 end associate
287 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
288 end do
289
290 do i = 1, vlag%size()
291 byte_offset = mpi_offset + &
292 dof_offset * int(mpi_real_prec_size, i8)
293 ! We should not need this extra associate block, and it works
294 ! great without it for GNU, Intel, NEC and Cray, but throws an
295 ! ICE with NAG.
296 associate(x => vlag%lf(i)%x)
297 call mpi_file_write_at_all(fh, byte_offset, x, &
298 vlag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
299 end associate
300 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
301 end do
302
303 do i = 1, wlag%size()
304 byte_offset = mpi_offset + &
305 dof_offset * int(mpi_real_prec_size, i8)
306 ! We should not need this extra associate block, and it works
307 ! great without it for GNU, Intel, NEC and Cray, but throws an
308 ! ICE with NAG.
309 associate(x => wlag%lf(i)%x)
310 call mpi_file_write_at_all(fh, byte_offset, x, &
311 wlag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
312 end associate
313 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
314 end do
315
316 end if
317
318 if (write_scalar) then
319 byte_offset = mpi_offset + &
320 dof_offset * int(mpi_real_prec_size, i8)
321 call mpi_file_write_at_all(fh, byte_offset, s%x, p%dof%size(), &
322 mpi_real_precision, status, ierr)
323 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
324 end if
325
326 if (write_dtlag) then
327 call mpi_file_write_at_all(fh, mpi_offset, tlag, 10, &
328 mpi_real_precision, status, ierr)
329 mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
330 call mpi_file_write_at_all(fh, mpi_offset, dtlag, 10, &
331 mpi_real_precision, status, ierr)
332 mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
333 end if
334
335 if (write_abvel) then
336 byte_offset = mpi_offset + &
337 dof_offset * int(mpi_real_prec_size, i8)
338 call mpi_file_write_at_all(fh, byte_offset, abx1%x, abx1%dof%size(), &
339 mpi_real_precision, status, ierr)
340 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
341 byte_offset = mpi_offset + &
342 dof_offset * int(mpi_real_prec_size, i8)
343 call mpi_file_write_at_all(fh, byte_offset, abx2%x, abx1%dof%size(), &
344 mpi_real_precision, status, ierr)
345 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
346 byte_offset = mpi_offset + &
347 dof_offset * int(mpi_real_prec_size, i8)
348 call mpi_file_write_at_all(fh, byte_offset, aby1%x, abx1%dof%size(), &
349 mpi_real_precision, status, ierr)
350 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
351 byte_offset = mpi_offset + &
352 dof_offset * int(mpi_real_prec_size, i8)
353 call mpi_file_write_at_all(fh, byte_offset, aby2%x, abx1%dof%size(), &
354 mpi_real_precision, status, ierr)
355 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
356 byte_offset = mpi_offset + &
357 dof_offset * int(mpi_real_prec_size, i8)
358 call mpi_file_write_at_all(fh, byte_offset, abz1%x, abx1%dof%size(), &
359 mpi_real_precision, status, ierr)
360 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
361 byte_offset = mpi_offset + &
362 dof_offset * int(mpi_real_prec_size, i8)
363 call mpi_file_write_at_all(fh, byte_offset, abz2%x, abx1%dof%size(), &
364 mpi_real_precision, status, ierr)
365 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
366 end if
367
368 if (write_scalarlag) then
369 do i = 1, slag%size()
370 byte_offset = mpi_offset + &
371 dof_offset * int(mpi_real_prec_size, i8)
372 ! We should not need this extra associate block, and it works
373 ! great without it for GNU, Intel, NEC and Cray, but throws an
374 ! ICE with NAG.
375 associate(x => slag%lf(i)%x)
376 call mpi_file_write_at_all(fh, byte_offset, x, &
377 slag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
378 end associate
379 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
380 end do
381
382 byte_offset = mpi_offset + &
383 dof_offset * int(mpi_real_prec_size, i8)
384 call mpi_file_write_at_all(fh, byte_offset, abs1%x, abx1%dof%size(), &
385 mpi_real_precision, status, ierr)
386 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
387 byte_offset = mpi_offset + &
388 dof_offset * int(mpi_real_prec_size, i8)
389 call mpi_file_write_at_all(fh, byte_offset, abs2%x, abx1%dof%size(), &
390 mpi_real_precision, status, ierr)
391 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
392 end if
393
394 if (write_ale) then
395 byte_offset = mpi_offset + &
396 dof_offset * int(mpi_real_prec_size, i8)
397 call mpi_file_write_at_all(fh, byte_offset, msh_x, size(msh_x), &
398 mpi_real_precision, status, ierr)
399 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
400
401 byte_offset = mpi_offset + &
402 dof_offset * int(mpi_real_prec_size, i8)
403 call mpi_file_write_at_all(fh, byte_offset, msh_y, size(msh_y), &
404 mpi_real_precision, status, ierr)
405 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
406
407 byte_offset = mpi_offset + &
408 dof_offset * int(mpi_real_prec_size, i8)
409 call mpi_file_write_at_all(fh, byte_offset, msh_z, size(msh_z), &
410 mpi_real_precision, status, ierr)
411 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
412
413 byte_offset = mpi_offset + &
414 dof_offset * int(mpi_real_prec_size, i8)
415 call mpi_file_write_at_all(fh, byte_offset, wm_x%x, wm_x%dof%size(), &
416 mpi_real_precision, status, ierr)
417 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
418
419 byte_offset = mpi_offset + &
420 dof_offset * int(mpi_real_prec_size, i8)
421 call mpi_file_write_at_all(fh, byte_offset, wm_y%x, wm_y%dof%size(), &
422 mpi_real_precision, status, ierr)
423 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
424
425 byte_offset = mpi_offset + &
426 dof_offset * int(mpi_real_prec_size, i8)
427 call mpi_file_write_at_all(fh, byte_offset, wm_z%x, wm_z%dof%size(), &
428 mpi_real_precision, status, ierr)
429 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
430
431 do i = 1, wm_x_lag%size()
432 byte_offset = mpi_offset + &
433 dof_offset * int(mpi_real_prec_size, i8)
434 ! We should not need this extra associate block, and it works
435 ! great without it for GNU, Intel, NEC and Cray, but throws an
436 ! ICE with NAG.
437 associate(x => wm_x_lag%lf(i)%x)
438 call mpi_file_write_at_all(fh, byte_offset, x, &
439 wm_x_lag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
440 end associate
441 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
442 end do
443
444 do i = 1, wm_y_lag%size()
445 byte_offset = mpi_offset + &
446 dof_offset * int(mpi_real_prec_size, i8)
447 ! We should not need this extra associate block, and it works
448 ! great without it for GNU, Intel, NEC and Cray, but throws an
449 ! ICE with NAG.
450 associate(x => wm_y_lag%lf(i)%x)
451 call mpi_file_write_at_all(fh, byte_offset, x, &
452 wm_y_lag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
453 end associate
454 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
455 end do
456
457 do i = 1, wm_z_lag%size()
458 byte_offset = mpi_offset + &
459 dof_offset * int(mpi_real_prec_size, i8)
460 ! We should not need this extra associate block, and it works
461 ! great without it for GNU, Intel, NEC and Cray, but throws an
462 ! ICE with NAG.
463 associate(x => wm_z_lag%lf(i)%x)
464 call mpi_file_write_at_all(fh, byte_offset, x, &
465 wm_z_lag%lf(i)%dof%size(), mpi_real_precision, status, ierr)
466 end associate
467 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
468 end do
469
470 byte_offset = mpi_offset + &
471 dof_offset * int(mpi_real_prec_size, i8)
472 call mpi_file_write_at_all(fh, byte_offset, blag, size(blag), &
473 mpi_real_precision, status, ierr)
474 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
475
476 byte_offset = mpi_offset + &
477 dof_offset * int(mpi_real_prec_size, i8)
478 call mpi_file_write_at_all(fh, byte_offset, blaglag, size(blaglag), &
479 mpi_real_precision, status, ierr)
480 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
481
482 call mpi_file_write_at_all(fh, mpi_offset, pivot_pos, size(pivot_pos), &
483 mpi_real_precision, status, ierr)
484 mpi_offset = mpi_offset + &
485 int(size(pivot_pos), i8) * int(mpi_real_prec_size, i8)
486 call mpi_file_write_at_all(fh, mpi_offset, pivot_vel_lag, &
487 size(pivot_vel_lag), mpi_real_precision, status, ierr)
488 mpi_offset = mpi_offset + &
489 int(size(pivot_vel_lag), i8) * int(mpi_real_prec_size, i8)
490 call mpi_file_write_at_all(fh, mpi_offset, basis_pos, &
491 size(basis_pos), mpi_real_precision, status, ierr)
492 mpi_offset = mpi_offset + &
493 int(size(basis_pos), i8) * int(mpi_real_prec_size, i8)
494 call mpi_file_write_at_all(fh, mpi_offset, basis_vel_lag, &
495 size(basis_vel_lag), mpi_real_precision, status, ierr)
496 mpi_offset = mpi_offset + &
497 int(size(basis_vel_lag), i8) * int(mpi_real_prec_size, i8)
498 end if
499
500 call mpi_file_close(fh, ierr)
501
502 if (ierr .ne. mpi_success) then
503 call neko_error('Error writing checkpoint file ' // trim(fname))
504 end if
505
506 end subroutine chkp_file_write
507
509 subroutine chkp_file_read(this, data)
510 class(chkp_file_t) :: this
511 class(*), target, intent(inout) :: data
512 type(chkp_t), pointer :: chkp
513 character(len=5) :: id_str
514 character(len=1024) :: fname
515 integer :: ierr, suffix_pos
516 type(field_t), pointer :: u, v, w, p, s
517 type(field_t), pointer :: wm_x => null()
518 type(field_t), pointer :: wm_y => null()
519 type(field_t), pointer :: wm_z => null()
520 type(field_series_t), pointer :: ulag => null()
521 type(field_series_t), pointer :: vlag => null()
522 type(field_series_t), pointer :: wlag => null()
523 type(field_series_t), pointer :: wm_x_lag => null()
524 type(field_series_t), pointer :: wm_y_lag => null()
525 type(field_series_t), pointer :: wm_z_lag => null()
526 type(field_series_t), pointer :: slag => null()
527 type(mesh_t), pointer :: msh
528 type(mpi_status) :: status
529 type(mpi_file) :: fh
530 type(field_t), pointer :: abx1, abx2
531 type(field_t), pointer :: aby1, aby2
532 type(field_t), pointer :: abz1, abz2
533 type(field_t), pointer :: abs1, abs2
534 real(kind=rp), pointer :: blag(:,:,:,:) => null()
535 real(kind=rp), pointer :: blaglag(:,:,:,:) => null()
536 real(kind=rp), pointer :: msh_x(:,:,:,:) => null()
537 real(kind=rp), pointer :: msh_y(:,:,:,:) => null()
538 real(kind=rp), pointer :: msh_z(:,:,:,:) => null()
539 real(kind=rp), pointer :: pivot_pos(:) => null()
540 real(kind=rp), pointer :: pivot_vel_lag(:,:) => null()
541 real(kind=rp), pointer :: basis_pos(:) => null()
542 real(kind=rp), pointer :: basis_vel_lag(:,:) => null()
543 real(kind=rp), allocatable :: x_coord(:,:,:,:)
544 real(kind=rp), allocatable :: y_coord(:,:,:,:)
545 real(kind=rp), allocatable :: z_coord(:,:,:,:)
546 real(kind=rp), pointer :: dtlag(:), tlag(:)
547 integer (kind=MPI_OFFSET_KIND) :: mpi_offset, byte_offset
548 integer(kind=i8) :: n_glb_dofs, dof_offset
549 integer :: glb_nelv, gdim, lx, have_lag, have_scalar, nel
550 integer :: optional_fields, have_dtlag
551 integer :: have_abvel, have_scalarlag
552 integer :: have_ale
553 logical :: read_lag, read_scalar, read_dtlag, read_abvel, read_scalarlag
554 logical :: read_ale
555 real(kind=dp) :: tol
556 real(kind=rp) :: center_x, center_y, center_z
557 integer :: i, e
558 type(dofmap_t) :: dof
559
560 call this%check_exists()
561
562 select type (data)
563 type is (chkp_t)
564
565 if ( .not. associated(data%u) .or. &
566 .not. associated(data%v) .or. &
567 .not. associated(data%w) .or. &
568 .not. associated(data%p) ) then
569 call neko_error('Checkpoint not initialized')
570 end if
571
572 u => data%u
573 v => data%v
574 w => data%w
575 p => data%p
576 this%chkp_Xh => data%previous_Xh
578 if (allocated(data%previous_mesh%elements)) then
579 msh => data%previous_mesh
580 this%mesh2mesh = .true.
581 tol = data%mesh2mesh_tol
582 else
583 msh => u%msh
584 this%mesh2mesh = .false.
585 end if
586
587 if (associated(data%ulag)) then
588 ulag => data%ulag
589 vlag => data%vlag
590 wlag => data%wlag
591 read_lag = .true.
592 else
593 read_lag = .false.
594 end if
595
596 if (associated(data%s)) then
597 s => data%s
598 read_scalar = .true.
599 else
600 read_scalar = .false.
601 end if
602 if (associated(data%dtlag)) then
603 dtlag => data%dtlag
604 tlag => data%tlag
605 read_dtlag = .true.
606 else
607 read_dtlag = .false.
608 end if
609 read_abvel = .false.
610 if (associated(data%abx1)) then
611 abx1 => data%abx1
612 abx2 => data%abx2
613 aby1 => data%aby1
614 aby2 => data%aby2
615 abz1 => data%abz1
616 abz2 => data%abz2
617 read_abvel = .true.
618 end if
619 read_scalarlag = .false.
620 if (associated(data%abs1)) then
621 slag => data%slag
622 abs1 => data%abs1
623 abs2 => data%abs2
624 read_scalarlag = .true.
625 end if
626
627 read_ale = .false.
628 if (associated(data%wm_x)) then
629 msh_x => data%msh_x
630 msh_y => data%msh_y
631 msh_z => data%msh_z
632 wm_x => data%wm_x
633 wm_y => data%wm_y
634 wm_z => data%wm_z
635 wm_x_lag => data%wm_x_lag
636 wm_y_lag => data%wm_y_lag
637 wm_z_lag => data%wm_z_lag
638 blag => data%Blag
639 blaglag => data%Blaglag
640 pivot_pos => data%pivot_pos
641 pivot_vel_lag => data%pivot_vel_lag
642 basis_pos => data%basis_pos
643 basis_vel_lag => data%basis_vel_lag
644 read_ale = .true.
645 end if
646
647 chkp => data
648
649 class default
650 call neko_error('Invalid data')
651 end select
652
653 fname = trim(this%get_fname())
654 call neko_log%message("Reading checkpoint from file: " // trim(fname), &
655 neko_log_verbose)
656 call mpi_file_open(neko_comm, trim(fname), &
657 mpi_mode_rdonly, mpi_info_null, fh, ierr)
658 call mpi_file_read_all(fh, glb_nelv, 1, mpi_integer, status, ierr)
659 call mpi_file_read_all(fh, gdim, 1, mpi_integer, status, ierr)
660 call mpi_file_read_all(fh, lx, 1, mpi_integer, status, ierr)
661 call mpi_file_read_all(fh, optional_fields, 1, mpi_integer, status, ierr)
662 call mpi_file_read_all(fh, chkp%t, 1, mpi_double_precision, status, ierr)
663
664 have_lag = mod(optional_fields,2)/1
665 have_scalar = mod(optional_fields,4)/2
666 have_dtlag = mod(optional_fields,8)/4
667 have_abvel = mod(optional_fields,16)/8
668 have_scalarlag = mod(optional_fields,32)/16
669 have_ale = mod(optional_fields,64)/32
670
671 if ( ( glb_nelv .ne. msh%glb_nelv ) .or. &
672 ( gdim .ne. msh%gdim) .or. &
673 ( (have_lag .eq. 0) .and. (read_lag) ) .or. &
674 ( (have_scalar .eq. 0) .and. (read_scalar) ) ) then
675 call neko_error('Checkpoint does not match case')
676 end if
677 nel = msh%nelv
678 this%sim_Xh => u%Xh
679 if (gdim .eq. 3) then
680 call this%chkp_Xh%init(gll, lx, lx, lx)
681 else
682 call this%chkp_Xh%init(gll, lx, lx)
683 end if
684 if (this%mesh2mesh) then
685 if (read_ale) then
686 call neko_error('ALE does not yet support mesh2mesh ' // &
687 'interpolation for restart!')
688 end if
689 call dof%init(msh, this%chkp_Xh)
690 call this%global_interp%init(dof, neko_comm, tol = tol)
691 call this%global_interp%find_points(u%dof%x, u%dof%y, u%dof%z, &
692 u%dof%size())
693 else
694 call this%space_interp%init(this%sim_Xh, this%chkp_Xh)
695 end if
696 dof_offset = int(msh%offset_el, i8) * int(this%chkp_Xh%lxyz, i8)
697 n_glb_dofs = int(this%chkp_Xh%lxyz, i8) * int(msh%glb_nelv, i8)
698
699 !
700 ! Read mandatory checkpoint data
701 !
702
703 byte_offset = 4_i8 * int(mpi_integer_size, i8) + &
704 int(mpi_double_precision_size, i8)
705 byte_offset = byte_offset + &
706 dof_offset * int(mpi_real_prec_size, i8)
707 call this%read_field(fh, byte_offset, u%x, nel)
708 mpi_offset = 4_i8 * int(mpi_integer_size, i8) + &
709 int(mpi_double_precision_size, i8)
710 mpi_offset = mpi_offset +&
711 n_glb_dofs * int(mpi_real_prec_size, i8)
712
713 byte_offset = mpi_offset + &
714 dof_offset * int(mpi_real_prec_size, i8)
715 call this%read_field(fh, byte_offset, v%x, nel)
716 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
717
718 byte_offset = mpi_offset + &
719 dof_offset * int(mpi_real_prec_size, i8)
720 call this%read_field(fh, byte_offset, w%x, nel)
721 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
722
723 byte_offset = mpi_offset + &
724 dof_offset * int(mpi_real_prec_size, i8)
725 call this%read_field(fh, byte_offset, p%x, nel)
726 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
727
728 !
729 ! Read optional payload
730 !
731
732 if (read_lag) then
733 do i = 1, ulag%size()
734 byte_offset = mpi_offset + &
735 dof_offset * int(mpi_real_prec_size, i8)
736 call this%read_field(fh, byte_offset, ulag%lf(i)%x, nel)
737 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
738 end do
739
740 do i = 1, vlag%size()
741 byte_offset = mpi_offset + &
742 dof_offset * int(mpi_real_prec_size, i8)
743 call this%read_field(fh, byte_offset, vlag%lf(i)%x, nel)
744 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
745 end do
746
747 do i = 1, wlag%size()
748 byte_offset = mpi_offset + &
749 dof_offset * int(mpi_real_prec_size, i8)
750 call this%read_field(fh, byte_offset, wlag%lf(i)%x, nel)
751 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
752 end do
753 end if
754
755 if (read_scalar) then
756 byte_offset = mpi_offset + &
757 dof_offset * int(mpi_real_prec_size, i8)
758 call this%read_field(fh, byte_offset, s%x, nel)
759 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
760 end if
761
762 if (read_dtlag .and. have_dtlag .eq. 1) then
763 call mpi_file_read_at_all(fh, mpi_offset, tlag, 10, &
764 mpi_real_precision, status, ierr)
765 mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
766 call mpi_file_read_at_all(fh, mpi_offset, dtlag, 10, &
767 mpi_real_precision, status, ierr)
768 mpi_offset = mpi_offset + 10_i8 * int(mpi_real_prec_size, i8)
769 end if
770
771 if (read_abvel .and. have_abvel .eq. 1) then
772 byte_offset = mpi_offset + &
773 dof_offset * int(mpi_real_prec_size, i8)
774 call this%read_field(fh, byte_offset, abx1%x, nel)
775 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
776 byte_offset = mpi_offset + &
777 dof_offset * int(mpi_real_prec_size, i8)
778 call this%read_field(fh, byte_offset, abx2%x, nel)
779 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
780 byte_offset = mpi_offset + &
781 dof_offset * int(mpi_real_prec_size, i8)
782 call this%read_field(fh, byte_offset, aby1%x, nel)
783 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
784 byte_offset = mpi_offset + &
785 dof_offset * int(mpi_real_prec_size, i8)
786 call this%read_field(fh, byte_offset, aby2%x, nel)
787 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
788 byte_offset = mpi_offset + &
789 dof_offset * int(mpi_real_prec_size, i8)
790 call this%read_field(fh, byte_offset, abz1%x, nel)
791 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
792 byte_offset = mpi_offset + &
793 dof_offset * int(mpi_real_prec_size, i8)
794 call this%read_field(fh, byte_offset, abz2%x, nel)
795 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
796 end if
797 if (read_scalarlag .and. have_scalarlag .eq. 1) then
798 do i = 1, slag%size()
799 byte_offset = mpi_offset + &
800 dof_offset * int(mpi_real_prec_size, i8)
801 call this%read_field(fh, byte_offset, slag%lf(i)%x, nel)
802 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
803 end do
804 byte_offset = mpi_offset + &
805 dof_offset * int(mpi_real_prec_size, i8)
806 call this%read_field(fh, byte_offset, abs1%x, nel)
807 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
808 byte_offset = mpi_offset + &
809 dof_offset * int(mpi_real_prec_size, i8)
810 call this%read_field(fh, byte_offset, abs2%x, nel)
811 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
812 end if
813
814 if (read_ale .and. have_ale .eq. 1) then
815
816 byte_offset = mpi_offset + &
817 dof_offset * int(mpi_real_prec_size, i8)
818 call this%read_field(fh, byte_offset, msh_x, nel)
819 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
820
821 byte_offset = mpi_offset + &
822 dof_offset * int(mpi_real_prec_size, i8)
823 call this%read_field(fh, byte_offset, msh_y, nel)
824 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
825
826 byte_offset = mpi_offset + &
827 dof_offset * int(mpi_real_prec_size, i8)
828 call this%read_field(fh, byte_offset, msh_z, nel)
829 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
830
831 byte_offset = mpi_offset + &
832 dof_offset * int(mpi_real_prec_size, i8)
833 call this%read_field(fh, byte_offset, wm_x%x, nel)
834 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
835
836 byte_offset = mpi_offset + &
837 dof_offset * int(mpi_real_prec_size, i8)
838 call this%read_field(fh, byte_offset, wm_y%x, nel)
839 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
840
841 byte_offset = mpi_offset + &
842 dof_offset * int(mpi_real_prec_size, i8)
843 call this%read_field(fh, byte_offset, wm_z%x, nel)
844 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
845
846 do i = 1, wm_x_lag%size()
847 byte_offset = mpi_offset + &
848 dof_offset * int(mpi_real_prec_size, i8)
849 call this%read_field(fh, byte_offset, wm_x_lag%lf(i)%x, nel)
850 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
851 end do
852 do i = 1, wm_y_lag%size()
853 byte_offset = mpi_offset + &
854 dof_offset * int(mpi_real_prec_size, i8)
855 call this%read_field(fh, byte_offset, wm_y_lag%lf(i)%x, nel)
856 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
857 end do
858 do i = 1, wm_z_lag%size()
859 byte_offset = mpi_offset + &
860 dof_offset * int(mpi_real_prec_size, i8)
861 call this%read_field(fh, byte_offset, wm_z_lag%lf(i)%x, nel)
862 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
863 end do
864 byte_offset = mpi_offset + &
865 dof_offset * int(mpi_real_prec_size, i8)
866 call this%read_field(fh, byte_offset, blag, nel)
867 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
868
869 byte_offset = mpi_offset + &
870 dof_offset * int(mpi_real_prec_size, i8)
871 call this%read_field(fh, byte_offset, blaglag, nel)
872 mpi_offset = mpi_offset + n_glb_dofs * int(mpi_real_prec_size, i8)
873
874 call mpi_file_read_at_all(fh, mpi_offset, pivot_pos, &
875 size(pivot_pos), mpi_real_precision, status, ierr)
876 mpi_offset = mpi_offset + int(size(pivot_pos), i8) &
877 * int(mpi_real_prec_size, i8)
878
879 call mpi_file_read_at_all(fh, mpi_offset, pivot_vel_lag, &
880 size(pivot_vel_lag), mpi_real_precision, status, ierr)
881 mpi_offset = mpi_offset + int(size(pivot_vel_lag), i8) &
882 * int(mpi_real_prec_size, i8)
883
884 call mpi_file_read_at_all(fh, mpi_offset, basis_pos, &
885 size(basis_pos), mpi_real_precision, status, ierr)
886 mpi_offset = mpi_offset + &
887 int(size(basis_pos), i8) * int(mpi_real_prec_size, i8)
888
889 call mpi_file_read_at_all(fh, mpi_offset, basis_vel_lag, &
890 size(basis_vel_lag), mpi_real_precision, status, ierr)
891 mpi_offset = mpi_offset + &
892 int(size(basis_vel_lag), i8) * int(mpi_real_prec_size, i8)
893 end if
894
895
896 call mpi_file_close(fh, ierr)
897
898 if (ierr .ne. mpi_success) then
899 call neko_error('Error reading checkpoint file ' // trim(fname))
900 end if
901
902 call this%global_interp%free()
903 call this%space_interp%free()
904
905 end subroutine chkp_file_read
906
907 subroutine chkp_read_field(this, fh, byte_offset, x, nel)
908 class(chkp_file_t) :: this
909 type(mpi_file) :: fh
910 integer(kind=MPI_OFFSET_KIND) :: byte_offset
911 integer, intent(in) :: nel
912 real(kind=rp) :: x(this%sim_Xh%lxyz, nel)
913 real(kind=rp), allocatable :: read_array(:)
914 integer :: nel_stride, frac_space
915 type(mpi_status) :: status
916 integer :: ierr, i, n
917 logical :: interp_on_host = .true.
918
919 n = this%chkp_xh%lxyz*nel
920 allocate(read_array(n))
921
922 call rzero(read_array,n)
923 call mpi_file_read_at_all(fh, byte_offset, read_array, &
924 n, mpi_real_precision, status, ierr)
925 if (this%mesh2mesh) then
926 x = 0.0_rp
927 call this%global_interp%evaluate(x, read_array, interp_on_host)
928
929 else if (this%sim_Xh%lx .ne. this%chkp_Xh%lx) then
930 call this%space_interp%map_host(x, read_array, nel, this%sim_Xh)
931 else
932 do i = 1,n
933 x(i,1) = read_array(i)
934 end do
935 end if
936 deallocate(read_array)
937 end subroutine chkp_read_field
938
939 subroutine chkp_file_set_overwrite(this, overwrite)
940 class(chkp_file_t), intent(inout) :: this
941 logical, intent(in) :: overwrite
942 this%overwrite = overwrite
943 end subroutine chkp_file_set_overwrite
944end 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:51
integer, public pe_rank
MPI rank.
Definition comm.F90:56
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
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:85
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
Definition math.f90:60
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:206
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:66
integer, public mpi_real_prec_size
Size of working precision REAL types.
Definition mpi_types.f90:70
integer, public mpi_integer_size
Size of MPI type integer.
Definition mpi_types.f90:68
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