Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
checkpoint_linear.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!
33
42submodule(simulation_checkpoint) checkpoint_linear
43 use simulation, only: simulation_step, simulation_restart
44 use file, only: file_t, file_free
45 use time_step_controller, only: time_step_controller_t
46
47contains
48
52 module subroutine checkpoint_save_linear(this, neko_case)
53 class(simulation_checkpoint_t), intent(inout) :: this
54 class(case_t), intent(inout) :: neko_case
55 logical :: save_disc
56
57 ! We save to disc only every n_saves_memory time steps
58 save_disc = modulo(neko_case%time%tstep, this%n_saves_memory) .eq. 0
59
60 ! Sample the checkpoint if needed
61 if (save_disc .or. neko_case%time%tstep .le. this%first_valid_timestep) then
62
63 call this%chkp_output%set_counter(neko_case%time%tstep -1)
64 call this%chkp_output%sample(neko_case%time%t)
65 this%n_saves_disc = this%n_saves_disc + 1
66 end if
67 end subroutine checkpoint_save_linear
68
73 module subroutine checkpoint_restore_linear(this, neko_case, tstep)
74 class(simulation_checkpoint_t), intent(inout) :: this
75 class(case_t), target, intent(inout) :: neko_case
76 integer, intent(in) :: tstep
77 type(time_step_controller_t) :: dt_controller
78 real(kind=dp) :: loop_start
79 integer :: j, k, previous_save, next_save
80 integer :: i_scalars
81 type(field_t), pointer :: u, v, w, p, s
82
83 loop_start = mpi_wtime()
84
85 u => neko_case%fluid%u
86 v => neko_case%fluid%v
87 w => neko_case%fluid%w
88 p => neko_case%fluid%p
89 s => null()
90
91 ! Determine the nearest save states on both sides
92 previous_save = tstep - modulo(tstep, this%n_saves_memory)
93 next_save = previous_save + this%n_saves_memory
94
95 ! Before the first valid state we always load from disc and we do not step
96 ! forward in time.
97 if (tstep .lt. this%first_valid_timestep) then
98 previous_save = tstep
99 next_save = previous_save + 1
100 else if (previous_save .lt. this%first_valid_timestep) then
101 previous_save = this%first_valid_timestep
102 end if
103
104 ! Load a new batch of checkpoints if needed
105 if (this%loaded_checkpoint .ne. previous_save) then
106
107 ! Restart the simulation form the checkpoint file
108 call this%chkp_output%set_counter(previous_save)
109 call this%chkp_output%file_%read(neko_case%chkp)
110 call simulation_restart(neko_case, neko_case%chkp)
111
112 ! Initialize the time step controller and set the time step
113 call dt_controller%init(neko_case%params)
114 neko_case%time%tstep = previous_save
115 this%loaded_checkpoint = neko_case%time%tstep
116
117 ! Step through the simulation and store field states in memory
118 do k = previous_save, min(next_save - 1, this%n_timesteps)
119
120 ! Do not run simulation step on the first iteration
121 if (k .ne. previous_save) then
122 if (neko_case%time%t .ge. neko_case%time%end_time) exit
123 call simulation_step(neko_case, dt_controller, loop_start)
124 end if
125
126 local_idx = modulo(k, this%n_saves_memory) + 1
127 call field_copy(this%p_list(local_idx), p)
128 call field_copy(this%u_list(local_idx), u)
129 call field_copy(this%v_list(local_idx), v)
130 call field_copy(this%w_list(local_idx), w)
131 do i_scalars = 1, this%n_scalars
132 j = (local_idx - 1) * this%n_scalars + i_scalars
133 s => neko_case%scalars%scalar_fields(i_scalars)%s
134 call field_copy(this%s_list(j), s)
135 end do
136
137 end do
138 end if
139
140 local_idx = modulo(tstep, this%n_saves_memory) + 1
141 call field_copy(p, this%p_list(local_idx))
142 call field_copy(u, this%u_list(local_idx))
143 call field_copy(v, this%v_list(local_idx))
144 call field_copy(w, this%w_list(local_idx))
145 do i_scalars = 1, this%n_scalars
146 j = (local_idx - 1) * this%n_scalars + i_scalars
147 s => neko_case%scalars%scalar_fields(i_scalars)%s
148 call field_copy(s, this%s_list(j))
149 end do
150
151
152 end subroutine checkpoint_restore_linear
153
154end submodule checkpoint_linear