37 use mpi_f08,
only: mpi_wtime
38 use case,
only: case_t
39 use num_types,
only: rp, dp
40 use time_scheme_controller,
only: time_scheme_controller_t
41 use file,
only: file_t
42 use logger,
only: log_size, neko_log
43 use jobctrl,
only: jobctrl_time_limit
44 use field,
only: field_t
45 use profiler,
only: profiler_start, profiler_stop, &
46 profiler_start_region, profiler_end_region
47 use json_utils,
only: json_get_or_default
48 use time_state,
only : time_state_t
49 use time_step_controller,
only: time_step_controller_t
63 type(time_step_controller_t),
intent(inout) :: dt_controller
64 character(len=LOG_SIZE) :: log_buf
67 call neko_log%section(
'Adjoint Starting simulation')
68 write(log_buf,
'(A, E15.7,A,E15.7,A)') &
69 'T : [', c%time%t,
',', c%time%end_time,
']'
70 call neko_log%message(log_buf)
71 if (.not. dt_controller%is_variable_dt)
then
72 write(log_buf,
'(A, E15.7)')
'dt : ', c%time%dt
74 write(log_buf,
'(A, E15.7)')
'CFL : ', dt_controller%cfl_trg
76 call neko_log%message(log_buf)
79 call neko_log%section(
'Postprocessing')
80 call c%output_controller%execute(c%time)
82 call c%case%user%initialize(c%time)
83 call neko_log%end_section()
84 call neko_log%newline()
91 logical :: output_at_end
94 call json_get_or_default(c%case%params,
'case.output_at_end', &
95 output_at_end, .true.)
96 call c%output_controller%execute(c%time, output_at_end)
98 if (.not. (output_at_end) .and. c%time%t .lt. c%time%end_time)
then
99 call simulation_adjoint_joblimit_chkp(c, c%time%t)
103 call c%case%user%finalize(c%time)
105 call neko_log%end_section(
'Normal end.')
111 tstep_loop_start_time, final_time)
113 real(kind=rp),
intent(inout) :: cfl
114 type(time_step_controller_t),
intent(inout) :: dt_controller
115 real(kind=dp),
intent(in) :: tstep_loop_start_time
116 real(kind=rp),
optional,
intent(in) :: final_time
117 real(kind=rp) :: t_bkp
118 real(kind=dp) :: start_time, end_time, tstep_start_time
119 character(len=LOG_SIZE) :: log_buf
122 call profiler_start_region(
'Time-Step Adjoint')
123 c%time%tstep = c%time%tstep + 1
124 start_time = mpi_wtime()
125 tstep_start_time = start_time
130 cfl = c%fluid_adj%compute_cfl(c%time%dt)
131 call dt_controller%set_dt(c%time, cfl)
132 if (dt_controller%is_variable_dt) cfl = c%fluid_adj%compute_cfl(c%time%dt)
138 call simulation_settime(c%time, c%fluid_adj%ext_bdf)
140 if (
present(final_time))
then
142 c%time%t = final_time - t_bkp
145 call neko_log%begin()
147 write(log_buf,
'(A,E15.7,1x,A,E15.7)')
'CFL:', cfl,
'dt:', c%time%dt
148 call neko_log%message(log_buf)
152 if (
allocated(c%adjoint_scalars))
then
153 start_time = mpi_wtime()
154 call neko_log%section(
'Adjoint scalar')
155 call c%adjoint_scalars%step(c%time, &
156 c%case%fluid%ext_bdf, dt_controller)
157 end_time = mpi_wtime()
158 write(log_buf,
'(A,E15.7)') &
159 'Scalar step time: ', end_time-start_time
160 call neko_log%end_section(log_buf)
164 call neko_log%section(
'Adjoint fluid')
165 call c%fluid_adj%step(c%time, dt_controller)
166 end_time = mpi_wtime()
167 write(log_buf,
'(A,E15.7)') &
168 'Fluid step time (s): ', end_time-start_time
169 call neko_log%end_section(log_buf)
172 call neko_log%section(
'Postprocessing')
175 call c%output_controller%execute(c%time)
177 call neko_log%end_section()
180 end_time = mpi_wtime()
181 call neko_log%section(
'Step summary')
182 write(log_buf,
'(A,I8,A,E15.7)') &
183 'Total time for step ', c%time%tstep,
' (s): ', end_time-tstep_start_time
184 call neko_log%message(log_buf)
185 write(log_buf,
'(A,E15.7)') &
186 'Total elapsed time (s): ', end_time-tstep_loop_start_time
187 call neko_log%message(log_buf)
189 call neko_log%end_section()
191 call profiler_end_region
193 if (
present(final_time))
then
199 subroutine simulation_settime(time, ext_bdf)
200 type(time_state_t),
intent(inout) :: time
201 type(time_scheme_controller_t),
intent(inout),
allocatable :: ext_bdf
204 if (
allocated(ext_bdf))
then
206 time%tlag(i) = time%tlag(i-1)
207 time%dtlag(i) = time%dtlag(i-1)
210 time%dtlag(1) = time%dt
211 time%tlag(1) = time%t
212 if (ext_bdf%ndiff .eq. 0)
then
213 time%dtlag(2) = time%dt
214 time%tlag(2) = time%t
217 call ext_bdf%set_coeffs(time%dtlag)
220 time%t = time%t + time%dt
222 end subroutine simulation_settime
228 type(file_t) :: chkpf, previous_meshf
229 character(len=LOG_SIZE) :: log_buf
230 character(len=:),
allocatable :: restart_file
231 character(len=:),
allocatable :: restart_mesh_file
235 call c%case%params%get(
'case.restart_file', restart_file, found)
236 call c%case%params%get(
'case.restart_mesh_file', restart_mesh_file, found)
239 call previous_meshf%init(trim(restart_mesh_file))
240 call previous_meshf%read(c%fluid_adj%chkp%previous_mesh)
243 call c%case%params%get(
'case.mesh2mesh_tolerance', tol, found)
245 if (found) c%case%fluid%chkp%mesh2mesh_tol = tol
247 call chkpf%init(trim(restart_file))
248 call chkpf%read(c%fluid_adj%chkp)
249 c%time%dtlag = c%fluid_adj%chkp%dtlag
250 c%time%tlag = c%fluid_adj%chkp%tlag
253 do i = 1,
size(c%time%dtlag)
254 call c%case%fluid%ext_bdf%set_coeffs(c%time%dtlag)
257 call c%fluid_adj%restart(c%case%chkp)
258 call c%case%fluid%chkp%previous_mesh%free()
259 if (
allocated(c%adjoint_scalars))
then
260 call c%adjoint_scalars%restart(c%case%chkp)
263 c%time%t = real(c%case%fluid%chkp%restart_time(), kind=rp)
264 call neko_log%section(
'Restarting from checkpoint')
265 write(log_buf,
'(A,A)')
'File : ', trim(restart_file)
266 call neko_log%message(log_buf)
267 write(log_buf,
'(A,E15.7)')
'Time : ', c%time%t
268 call neko_log%message(log_buf)
269 call neko_log%end_section()
271 call c%output_controller%set_counter(c%time)
275 subroutine simulation_adjoint_joblimit_chkp(C, t)
277 real(kind=rp),
intent(inout) :: t
278 type(file_t) :: chkpf
279 character(len=:),
allocatable :: chkp_format
280 character(len=LOG_SIZE) :: log_buf
281 character(len=10) :: format_str
284 call c%case%params%get(
'case.checkpoint_format', chkp_format, found)
285 call c%case%fluid%chkp%sync_host()
288 if (chkp_format .eq.
'hdf5')
then
292 call chkpf%init(c%case%output_directory //
'joblimit' // trim(format_str))
293 call chkpf%write(c%case%fluid%chkp, t)
294 write(log_buf,
'(A)')
'! saving checkpoint >>>'
295 call neko_log%message(log_buf)
297 end subroutine simulation_adjoint_joblimit_chkp
Adjoint simulation driver.
subroutine, public simulation_adjoint_init(c, dt_controller)
Initialise a simulation_adjoint of a case.
subroutine, public simulation_adjoint_step(c, dt_controller, cfl, tstep_loop_start_time, final_time)
Compute a single time-step of an adjoint case.
subroutine, public simulation_adjoint_restart(c)
Restart a case C from a given checkpoint.
subroutine, public simulation_adjoint_finalize(c)
Finalize a simulation of a case.
Adjoint case type. Todo: This should Ideally be a subclass of case_t, however, that is not yet suppor...