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