Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
checkpoint.f90
1! Copyright (c) 2025, 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!
33module simulation_checkpoint
34 use num_types, only: rp, sp, dp
35 use case, only: case_t
36 use json_file_module, only: json_file
37 use json_utils, only: json_get_or_default
38 use scalar_scheme, only: scalar_scheme_t
39 use time_state, only: time_state_t
40 use chkp_output, only: chkp_output_t
41 use field, only: field_t
42 use mpi_f08, only: mpi_wtime
43 use utils, only: neko_error
44 use field_math, only: field_copy, field_rzero
45 implicit none
46 private
47
49 private
50
51 ! ----------------------------------------------------------------------- !
52 ! User parameters
53
55 character(len=256) :: algorithm = "linear"
57 character(len=256) :: filename = "checkpoint"
59 integer :: n_saves_memory = 10
60
61 ! Internal parameters
62 integer :: n_saves_disc = 0
63 integer :: n_timesteps = 0
64 integer :: first_valid_timestep = 2
65 integer :: loaded_checkpoint = -1
66
67 ! Structures to hold the checkpoint data
68 type(chkp_output_t) :: chkp_output
69 type(field_t), dimension(:), allocatable :: p_list
70 type(field_t), dimension(:), allocatable :: u_list
71 type(field_t), dimension(:), allocatable :: v_list
72 type(field_t), dimension(:), allocatable :: w_list
73
74 integer :: n_scalars = 0
75 type(field_t), dimension(:), allocatable :: s_list
76
77 contains
79 generic, public :: init => init_from_json, init_from_components
81 procedure, public, pass(this) :: init_from_json => &
82 checkpoint_init_from_json
84 procedure, public, pass(this) :: init_from_components => &
85 checkpoint_init_from_components
87 procedure, public, pass(this) :: free => checkpoint_free
89 procedure, public, pass(this) :: reset => checkpoint_reset
91 procedure, public, pass(this) :: save => checkpoint_save
93 procedure, public, pass(this) :: restore => checkpoint_restore
94
96
97 ! ========================================================================== !
98 ! Module procedures for our algorithm implementations.
99
100
101 interface
102
103 module subroutine checkpoint_save_linear(this, neko_case)
104 class(simulation_checkpoint_t), intent(inout) :: this
105 class(case_t), intent(inout) :: neko_case
106 end subroutine checkpoint_save_linear
107
109 module subroutine checkpoint_restore_linear(this, neko_case, tstep)
110 class(simulation_checkpoint_t), intent(inout) :: this
111 class(case_t), target, intent(inout) :: neko_case
112 integer, intent(in) :: tstep
113 end subroutine checkpoint_restore_linear
114 end interface
115
116contains
117
118 ! ========================================================================== !
119 ! Initialization and deallocation
120
122 subroutine checkpoint_init_from_json(this, neko_case, params)
123 class(simulation_checkpoint_t), intent(inout) :: this
124 class(case_t), target, intent(inout) :: neko_case
125 type(json_file), target, intent(inout) :: params
126 integer :: n_saves_memory
127 character(len=:), allocatable :: filename, algorithm
128
129 call json_get_or_default(params, "algorithm", algorithm, "linear")
130 call json_get_or_default(params, "n_memory", n_saves_memory, 10)
131 call json_get_or_default(params, "filename", filename, "checkpoint.chkp")
132
133 call this%init_from_components(neko_case, algorithm, n_saves_memory, &
134 filename)
135 end subroutine checkpoint_init_from_json
136
138 subroutine checkpoint_init_from_components(this, neko_case, algorithm, &
139 n_saves_memory, filename)
140 class(simulation_checkpoint_t), intent(inout), target :: this
141 class(case_t), target, intent(inout) :: neko_case
142 character(len=*), optional, intent(in) :: algorithm
143 integer, optional, intent(in) :: n_saves_memory
144 character(len=*), optional, intent(in) :: filename
145
146 class(scalar_scheme_t), pointer :: scalar_i
147 integer :: i, j
148 character(len=80) :: str
149
150 call this%free()
151
152 ! Set internal parameters
153 if (present(algorithm)) this%algorithm = algorithm
154 if (present(filename)) this%filename = filename
155 if (present(n_saves_memory)) this%n_saves_memory = n_saves_memory
156 if (allocated(neko_case%scalars)) then
157 this%n_scalars = size(neko_case%scalars%scalar_fields)
158 end if
159
160 ! Initialize the Neko checkpoint output
161 call this%chkp_output%init(neko_case%chkp, this%filename)
162
163 ! Allocate the RAM Checkpoints
164 allocate(this%p_list(this%n_saves_memory))
165 allocate(this%u_list(this%n_saves_memory))
166 allocate(this%v_list(this%n_saves_memory))
167 allocate(this%w_list(this%n_saves_memory))
168 if (this%n_scalars .gt. 0) then
169 this%n_scalars = size(neko_case%scalars%scalar_fields)
170 allocate(this%s_list(this%n_saves_memory * this%n_scalars))
171 end if
172
173 do i = 1, this%n_saves_memory
174 write(str, '(A,I0)') "p_chkp_", i
175 call this%p_list(i)%init(neko_case%fluid%p%dof, str)
176 write(str, '(A,I0)') "u_chkp_", i
177 call this%u_list(i)%init(neko_case%fluid%u%dof, str)
178 write(str, '(A,I0)') "v_chkp_", i
179 call this%v_list(i)%init(neko_case%fluid%v%dof, str)
180 write(str, '(A,I0)') "w_chkp_", i
181 call this%w_list(i)%init(neko_case%fluid%w%dof, str)
182 if (this%n_scalars .gt. 0) then
183 do j = 1, this%n_scalars
184 write(str, '(A,I0,A,I0)') "s_chkp_", i, "_", j
185 scalar_i => neko_case%scalars%scalar_fields(j)
186 call this%s_list((i - 1) * this%n_scalars + j)%init(scalar_i%s%dof, str)
187 end do
188 end if
189 end do
190
191 end subroutine checkpoint_init_from_components
192
194 subroutine checkpoint_free(this)
195 class(simulation_checkpoint_t), intent(inout) :: this
196 integer :: i
197
198 ! Free the RAM Checkpoints
199 do i = 1, this%n_saves_memory
200 if (allocated(this%p_list)) call this%p_list(i)%free()
201 if (allocated(this%u_list)) call this%u_list(i)%free()
202 if (allocated(this%v_list)) call this%v_list(i)%free()
203 if (allocated(this%w_list)) call this%w_list(i)%free()
204 end do
205
206 if (allocated(this%s_list)) then
207 do i = 1, size(this%s_list)
208 call this%s_list(i)%free()
209 end do
210 end if
211
212 if (allocated(this%p_list)) deallocate(this%p_list)
213 if (allocated(this%u_list)) deallocate(this%u_list)
214 if (allocated(this%v_list)) deallocate(this%v_list)
215 if (allocated(this%w_list)) deallocate(this%w_list)
216 if (allocated(this%s_list)) deallocate(this%s_list)
217
218 end subroutine checkpoint_free
219
220 ! ========================================================================== !
221 ! Saving and Restoring
222
224 subroutine checkpoint_save(this, neko_case)
225 class(simulation_checkpoint_t), intent(inout) :: this
226 class(case_t), intent(inout) :: neko_case
227
228 ! Update the number of recorded timesteps
229 this%n_timesteps = this%n_timesteps + 1
230
231 select case (this%algorithm)
232 case ("linear")
233 call checkpoint_save_linear(this, neko_case)
234 case default
235 call neko_error("Unknown checkpoint algorithm: " // this%algorithm)
236 end select
237 end subroutine checkpoint_save
238
240 subroutine checkpoint_restore(this, neko_case, tstep)
241 class(simulation_checkpoint_t), intent(inout) :: this
242 class(case_t), target, intent(inout) :: neko_case
243 integer, intent(in) :: tstep
244 character(len=256) :: msg
245
246 if (tstep .lt. 1 .or. tstep .gt. this%n_timesteps) then
247 write(msg, '(A,I0,A,I0,A)') "Requested timestep ", tstep, &
248 " is out of range [1, ", this%n_timesteps, "]"
249 call neko_error(trim(msg))
250 end if
251
252 select case (this%algorithm)
253 case ("linear")
254 call checkpoint_restore_linear(this, neko_case, tstep)
255 case default
256 call neko_error("Unknown checkpoint algorithm: " // this%algorithm)
257 end select
258 end subroutine checkpoint_restore
259
260 ! ========================================================================== !
261 ! Meta handling
262
264 subroutine checkpoint_reset(this)
265 class(simulation_checkpoint_t), intent(inout) :: this
266 integer :: i
267
268 ! Reset our checkpoints
269 this%loaded_checkpoint = -1
270 this%n_saves_disc = 0
271 this%n_timesteps = 0
272
273 do i = 1, this%n_saves_memory
274 call field_rzero(this%p_list(i))
275 call field_rzero(this%u_list(i))
276 call field_rzero(this%v_list(i))
277 call field_rzero(this%w_list(i))
278 end do
279
280 do i = 1, size(this%s_list)
281 call field_rzero(this%s_list(i))
282 end do
283 end subroutine checkpoint_reset
284
285end module simulation_checkpoint