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