Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
checkpoint.f90
Go to the documentation of this file.
1! Copyright (c) 2021, 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!
36 use num_types, only : rp, dp
39 use space, only : space_t, operator(.ne.)
42 use field, only : field_t, field_ptr_t
43 use utils, only : neko_error
44 use mesh, only : mesh_t
45 use math, only : neko_eps
46 implicit none
47 private
48
49 type, public :: chkp_t
50 type(field_t), pointer :: u => null()
51 type(field_t), pointer :: v => null()
52 type(field_t), pointer :: w => null()
53 type(field_t), pointer :: p => null()
54
55
56 !
57 ! Optional payload
58 !
59 type(field_series_t), pointer :: ulag => null()
60 type(field_series_t), pointer :: vlag => null()
61 type(field_series_t), pointer :: wlag => null()
62
63 real(kind=rp), pointer :: tlag(:) => null()
64 real(kind=rp), pointer :: dtlag(:) => null()
65
67 type(field_t), pointer :: abx1 => null()
68 type(field_t), pointer :: abx2 => null()
69 type(field_t), pointer :: aby1 => null()
70 type(field_t), pointer :: aby2 => null()
71 type(field_t), pointer :: abz1 => null()
72 type(field_t), pointer :: abz2 => null()
73
74 type(field_t), pointer :: s => null()
75 type(field_series_t), pointer :: slag => null()
76 type(field_t), pointer :: abs1 => null()
77 type(field_t), pointer :: abs2 => null()
78
79 type(field_series_list_t) :: scalar_lags
80
82 type(field_ptr_t), allocatable :: scalar_abx1(:)
83 type(field_ptr_t), allocatable :: scalar_abx2(:)
84
85 real(kind=dp) :: t = 0d0
86 type(mesh_t) :: previous_mesh
87 type(space_t) :: previous_xh
88 real(kind=rp) :: mesh2mesh_tol = neko_eps*1e3_rp
89
90 contains
91 procedure, pass(this) :: init => chkp_init
92 procedure, pass(this) :: sync_host => chkp_sync_host
93 procedure, pass(this) :: sync_device => chkp_sync_device
94 procedure, pass(this) :: add_fluid => chkp_add_fluid
95 procedure, pass(this) :: add_lag => chkp_add_lag
96 procedure, pass(this) :: add_scalar => chkp_add_scalar
97 procedure, pass(this) :: restart_time => chkp_restart_time
98 procedure, pass(this) :: free => chkp_free
99 end type chkp_t
100
101contains
102
104 subroutine chkp_init(this)
105 class(chkp_t), intent(inout) :: this
106
107 ! Make sure the object is clean
108 call this%free()
109
110 end subroutine chkp_init
111
113 subroutine chkp_free(this)
114 class(chkp_t), intent(inout) :: this
115
116 this%t = 0d0
117
118 if (associated(this%u)) nullify(this%u)
119 if (associated(this%v)) nullify(this%v)
120 if (associated(this%w)) nullify(this%w)
121 if (associated(this%p)) nullify(this%p)
122
123 if (associated(this%ulag)) nullify(this%ulag)
124 if (associated(this%vlag)) nullify(this%vlag)
125 if (associated(this%wlag)) nullify(this%wlag)
126
127 if (associated(this%tlag)) nullify(this%tlag)
128 if (associated(this%dtlag)) nullify(this%dtlag)
129
130 if (associated(this%abx1)) nullify(this%abx1)
131 if (associated(this%abx2)) nullify(this%abx2)
132 if (associated(this%aby1)) nullify(this%aby1)
133 if (associated(this%aby2)) nullify(this%aby2)
134 if (associated(this%abz1)) nullify(this%abz1)
135 if (associated(this%abz2)) nullify(this%abz2)
136
137 if (associated(this%s)) nullify(this%s)
138 if (associated(this%slag)) nullify(this%slag)
139 if (associated(this%abs1)) nullify(this%abs1)
140 if (associated(this%abs2)) nullify(this%abs2)
141
142 call this%scalar_lags%free()
143
144 if (allocated(this%scalar_abx1)) then
145 deallocate(this%scalar_abx1)
146 end if
147
148 if (allocated(this%scalar_abx2)) then
149 deallocate(this%scalar_abx2)
150 end if
151
152 call this%previous_mesh%free()
153 call this%previous_Xh%free()
154
155 end subroutine chkp_free
156
158 subroutine chkp_sync_host(this)
159 class(chkp_t), intent(inout) :: this
160 integer :: i, j
161
162 if (neko_bcknd_device .eq. 1) then
163 associate(u => this%u, v => this%v, w => this%w, &
164 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag, &
165 p => this%p)
166
167 if (associated(this%u) .and. associated(this%v) .and. &
168 associated(this%w) .and. associated(this%p)) then
169 call u%copy_from(device_to_host, sync = .false.)
170 call v%copy_from(device_to_host, sync = .false.)
171 call w%copy_from(device_to_host, sync = .false.)
172 call p%copy_from(device_to_host, sync = .false.)
173 end if
174
175 if (associated(this%ulag) .and. associated(this%vlag) .and. &
176 associated(this%wlag)) then
177 call ulag%lf(1)%copy_from(device_to_host, sync = .false.)
178 call ulag%lf(2)%copy_from(device_to_host, sync = .false.)
179
180 call vlag%lf(1)%copy_from(device_to_host, sync = .false.)
181 call vlag%lf(2)%copy_from(device_to_host, sync = .false.)
182
183 call wlag%lf(1)%copy_from(device_to_host, sync = .false.)
184 call wlag%lf(2)%copy_from(device_to_host, sync = .false.)
185 call this%abx1%copy_from(device_to_host, sync = .false.)
186 call this%abx2%copy_from(device_to_host, sync = .false.)
187 call this%aby1%copy_from(device_to_host, sync = .false.)
188 call this%aby2%copy_from(device_to_host, sync = .false.)
189 call this%abz1%copy_from(device_to_host, sync = .false.)
190 call this%abz2%copy_from(device_to_host, sync = .false.)
191 end if
192
193 if (associated(this%s)) then
194 call this%s%copy_from(device_to_host, sync = .false.)
195 call this%slag%lf(1)%copy_from(device_to_host, sync = .false.)
196 call this%slag%lf(2)%copy_from(device_to_host, sync = .false.)
197 call this%abs1%copy_from(device_to_host, sync = .false.)
198 call this%abs2%copy_from(device_to_host, sync = .false.)
199 end if
200
201 ! Multi-scalar lag field synchronization
202 do i = 1, this%scalar_lags%size()
203 block
204 type(field_series_t), pointer :: slag
205 integer :: slag_size
206 slag => this%scalar_lags%get(i)
207 slag_size = slag%size()
208 do j = 1, slag_size
209 call slag%lf(j)%copy_from(device_to_host, sync = .false.)
210 end do
211 end block
212 end do
213
214 ! Multi-scalar ABX field synchronization
215 if (allocated(this%scalar_abx1) .and. allocated(this%scalar_abx2)) then
216 do i = 1, size(this%scalar_abx1)
217 call this%scalar_abx1(i)%ptr%copy_from(device_to_host, &
218 sync = .false.)
219 call this%scalar_abx2(i)%ptr%copy_from(device_to_host, &
220 sync = .false.)
221 end do
222 end if
223 end associate
224 call device_sync(glb_cmd_queue)
225 end if
226
227 end subroutine chkp_sync_host
228
230 subroutine chkp_sync_device(this)
231 class(chkp_t), intent(inout) :: this
232 integer :: i, j
233
234 if (neko_bcknd_device .eq. 1) then
235 associate(u => this%u, v => this%v, w => this%w, &
236 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag, &
237 p => this%p)
238
239 if (associated(this%u) .and. associated(this%v) .and. &
240 associated(this%w)) then
241 call u%copy_from(host_to_device, sync = .false.)
242 call v%copy_from(host_to_device, sync = .false.)
243 call w%copy_from(host_to_device, sync = .false.)
244 call p%copy_from(host_to_device, sync = .false.)
245 end if
246
247 if (associated(this%ulag) .and. associated(this%vlag) .and. &
248 associated(this%wlag)) then
249 call ulag%lf(1)%copy_from(host_to_device, sync = .false.)
250 call ulag%lf(2)%copy_from(host_to_device, sync = .false.)
251
252 call vlag%lf(1)%copy_from(host_to_device, sync = .false.)
253 call vlag%lf(2)%copy_from(host_to_device, sync = .false.)
254
255 call wlag%lf(1)%copy_from(host_to_device, sync = .false.)
256 call wlag%lf(2)%copy_from(host_to_device, sync = .false.)
257 end if
258
259 if (associated(this%s)) then
260 call this%s%copy_from(host_to_device, sync = .false.)
261
262 call this%slag%lf(1)%copy_from(host_to_device, sync = .false.)
263 call this%slag%lf(2)%copy_from(host_to_device, sync = .false.)
264 call this%abs1%copy_from(host_to_device, sync = .false.)
265 call this%abs2%copy_from(host_to_device, sync = .false.)
266 end if
267
268 ! Multi-scalar lag field synchronization
269 if (allocated(this%scalar_lags%items) .and. &
270 this%scalar_lags%size() > 0) then
271 do i = 1, this%scalar_lags%size()
272 block
273 type(field_series_t), pointer :: slag
274 integer :: slag_size, dof_size
275 slag => this%scalar_lags%get(i)
276 slag_size = slag%size()
277 dof_size = slag%f%dof%size()
278 do j = 1, slag_size
279 call slag%lf(j)%copy_from(host_to_device, sync = .false.)
280 end do
281 end block
282 end do
283 end if
284
285 ! Multi-scalar ABX field synchronization
286 if (allocated(this%scalar_abx1) .and. allocated(this%scalar_abx2)) then
287 do i = 1, size(this%scalar_abx1)
288 call this%scalar_abx1(i)%ptr%copy_from(host_to_device, &
289 sync = .false.)
290 call this%scalar_abx2(i)%ptr%copy_from(host_to_device, &
291 sync = .false.)
292 end do
293 end if
294 end associate
295 end if
296
297 end subroutine chkp_sync_device
298
300 subroutine chkp_add_fluid(this, u, v, w, p)
301 class(chkp_t), intent(inout) :: this
302 type(field_t), target :: u
303 type(field_t), target :: v
304 type(field_t), target :: w
305 type(field_t), target :: p
306
307 ! Check that all velocity components are defined on the same
308 ! function space
309 if ( u%Xh .ne. v%Xh .or. &
310 u%Xh .ne. w%Xh ) then
311 call neko_error('Different function spaces for each velocity component')
312 end if
313
314 ! Check that both velocity and pressure is defined on the same mesh
315 if ( u%msh%nelv .ne. p%msh%nelv ) then
316 call neko_error('Velocity and pressure defined on different meshes')
317 end if
318
319 this%u => u
320 this%v => v
321 this%w => w
322 this%p => p
323
324 end subroutine chkp_add_fluid
325
327 subroutine chkp_add_lag(this, ulag, vlag, wlag)
328 class(chkp_t), intent(inout) :: this
329 type(field_series_t), target :: ulag
330 type(field_series_t), target :: vlag
331 type(field_series_t), target :: wlag
332
333 this%ulag => ulag
334 this%vlag => vlag
335 this%wlag => wlag
336
337 end subroutine chkp_add_lag
338
339
340
342 subroutine chkp_add_scalar(this, s, slag, abs1, abs2)
343 class(chkp_t), intent(inout) :: this
344 type(field_t), target, intent(in) :: s
345 type(field_series_t), target, intent(in) :: slag
346 type(field_t), target, intent(in), optional :: abs1, abs2
347
348 this%s => s
349 this%slag => slag
350
351 if (present(abs1)) this%abs1 => abs1
352 if (present(abs2)) this%abs2 => abs2
353
354 end subroutine chkp_add_scalar
355
356
358 pure function chkp_restart_time(this) result(rtime)
359 class(chkp_t), intent(in) :: this
360 real(kind=dp) :: rtime
361
362 rtime = this%t
363 end function chkp_restart_time
364
365end module checkpoint
Copy data between host and device (or device and device)
Definition device.F90:71
Synchronize a device or stream.
Definition device.F90:107
Defines a checkpoint.
subroutine chkp_sync_host(this)
Synchronize checkpoint with device.
subroutine chkp_sync_device(this)
Synchronize device with checkpoint.
pure real(kind=dp) function chkp_restart_time(this)
Return restart time from a loaded checkpoint.
subroutine chkp_free(this)
Reset checkpoint.
subroutine chkp_add_lag(this, ulag, vlag, wlag)
Add lagged velocity terms.
subroutine chkp_add_fluid(this, u, v, w, p)
Add a fluid to the checkpoint.
subroutine chkp_init(this)
Initialize checkpoint structure with mandatory data.
subroutine chkp_add_scalar(this, s, slag, abs1, abs2)
Add a scalar to checkpointing.
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
integer, parameter, public device_to_host
Definition device.F90:47
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
Definition device.F90:51
Contains the field_series_list_t type for managing multiple field series.
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
Definition math.f90:60
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:69
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
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
field_ptr_t, To easily obtain a pointer to a field
Definition field.f90:82
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
A list of field series pointers, used for managing multiple scalar lag fields.
The function space for the SEM solution fields.
Definition space.f90:63