Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
simulation_adjoint.f90
1! Copyright (c) 2020-2021, 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!
35 use mpi_f08, only: mpi_wtime
36 use case, only: case_t
37 use num_types, only: rp, dp
38 use time_scheme_controller, only: time_scheme_controller_t
39 use file, only: file_t
40 use logger, only: log_size, neko_log
41 use jobctrl, only: jobctrl_time_limit
42 use field, only: field_t
43 use profiler, only: profiler_start, profiler_stop, &
44 profiler_start_region, profiler_end_region
45 use json_utils, only: json_get_or_default
46 use time_state, only : time_state_t
47 use time_step_controller, only: time_step_controller_t
48 use adjoint_case, only: adjoint_case_t
49 implicit none
50 private
51
52 public :: simulation_adjoint_init, simulation_adjoint_step, &
53 simulation_adjoint_finalize, simulation_adjoint_restart
54
55contains
56
57
59 subroutine simulation_adjoint_init(C, dt_controller)
60 type(adjoint_case_t), intent(inout) :: c
61 type(time_step_controller_t), intent(inout) :: dt_controller
62 character(len=LOG_SIZE) :: log_buf
63
64 ! Write the initial logging message
65 call neko_log%section('Adjoint Starting simulation')
66 write(log_buf, '(A, E15.7,A,E15.7,A)') &
67 'T : [', c%time%t, ',', c%time%end_time, ']'
68 call neko_log%message(log_buf)
69 if (.not. dt_controller%is_variable_dt) then
70 write(log_buf, '(A, E15.7)') 'dt : ', c%time%dt
71 else
72 write(log_buf, '(A, E15.7)') 'CFL : ', dt_controller%cfl_trg
73 end if
74 call neko_log%message(log_buf)
75
76 ! Call stats, samplers and user-init before time loop
77 call neko_log%section('Postprocessing')
78 call c%output_controller%execute(c%time)
79
80 call c%case%user%user_init_modules(c%time%t, c%fluid_adj%u_adj, &
81 c%fluid_adj%v_adj, c%fluid_adj%w_adj, &
82 c%fluid_adj%p_adj, c%fluid_adj%c_Xh, c%case%params)
83 call neko_log%end_section()
84 call neko_log%newline()
85
86 end subroutine simulation_adjoint_init
87
89 subroutine simulation_adjoint_finalize(C)
90 type(adjoint_case_t), intent(inout) :: c
91 logical :: output_at_end
92
93 ! Run a final output if specified in the json
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)
97
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)
100 end if
101
102 ! Finalize the user modules
103 call c%case%user%user_finalize_modules(c%time%t, c%case%params)
104
105 call neko_log%end_section('Normal end.')
106
107 end subroutine simulation_adjoint_finalize
108
110 subroutine simulation_adjoint_step(C, dt_controller, cfl, tstep_loop_start_time)
111 type(adjoint_case_t), intent(inout) :: c
112 real(kind=rp), intent(inout) :: cfl
113 type(time_step_controller_t), intent(inout) :: dt_controller
114 real(kind=dp), intent(in) :: tstep_loop_start_time
115 real(kind=dp) :: start_time, end_time, tstep_start_time
116 character(len=LOG_SIZE) :: log_buf
117
118 ! Setup the time step, and start time
119 call profiler_start_region('Time-Step Adjoint')
120 c%time%tstep = c%time%tstep + 1
121 start_time = mpi_wtime()
122 tstep_start_time = start_time
123
124 ! Compute the next time step
125 ! NOTE. we should be wary here since CFL is based on the convective velocity
126 ! not the adjoint velocity
127 cfl = c%fluid_adj%compute_cfl(c%time%dt)
128 call dt_controller%set_dt(c%time, cfl)
129 if (dt_controller%is_variable_dt) cfl = c%fluid_adj%compute_cfl(c%time%dt)
130
131 ! Calculate the cfl after the possibly varied dt
132 ! cfl = C%fluid_adj%compute_cfl(C%time%dt)
133
134 ! Advance time step from t to t+dt and print the status
135 call simulation_settime(c%time, c%fluid_adj%ext_bdf)
136 call c%time%status()
137 call neko_log%begin()
138
139 write(log_buf, '(A,E15.7,1x,A,E15.7)') 'CFL:', cfl, 'dt:', c%time%dt
140 call neko_log%message(log_buf)
141
142 ! Scalar step
143 ! (Note that for the adjoint we should the adjoint_scalars first)
144 if (allocated(c%adjoint_scalars)) then
145 start_time = mpi_wtime()
146 call neko_log%section('Adjoint scalar')
147 call c%adjoint_scalars%step(c%time, &
148 c%case%fluid%ext_bdf, dt_controller)
149 end_time = mpi_wtime()
150 write(log_buf, '(A,E15.7)') &
151 'Scalar step time: ', end_time-start_time
152 call neko_log%end_section(log_buf)
153 end if
154
155 ! Fluid step
156 call neko_log%section('Adjoint fluid')
157 call c%fluid_adj%step(c%time, dt_controller)
158 end_time = mpi_wtime()
159 write(log_buf, '(A,E15.7)') &
160 'Fluid step time (s): ', end_time-start_time
161 call neko_log%end_section(log_buf)
162
163 ! Postprocessing
164 call neko_log%section('Postprocessing')
165
166 ! Run any IO needed.
167 call c%output_controller%execute(c%time)
168
169 call neko_log%end_section()
170
171 ! End the step and print summary
172 end_time = mpi_wtime()
173 call neko_log%section('Step summary')
174 write(log_buf, '(A,I8,A,E15.7)') &
175 'Total time for step ', c%time%tstep, ' (s): ', end_time-tstep_start_time
176 call neko_log%message(log_buf)
177 write(log_buf, '(A,E15.7)') &
178 'Total elapsed time (s): ', end_time-tstep_loop_start_time
179 call neko_log%message(log_buf)
180
181 call neko_log%end_section()
182 call neko_log%end()
183 call profiler_end_region
184
185 end subroutine simulation_adjoint_step
186
187 subroutine simulation_settime(time, ext_bdf)
188 type(time_state_t), intent(inout) :: time
189 type(time_scheme_controller_t), intent(inout), allocatable :: ext_bdf
190 integer :: i
191
192 if (allocated(ext_bdf)) then
193 do i = 10, 2, -1
194 time%tlag(i) = time%tlag(i-1)
195 time%dtlag(i) = time%dtlag(i-1)
196 end do
197
198 time%dtlag(1) = time%dt
199 time%tlag(1) = time%t
200 if (ext_bdf%ndiff .eq. 0) then
201 time%dtlag(2) = time%dt
202 time%tlag(2) = time%t
203 end if
204
205 call ext_bdf%set_coeffs(time%dtlag)
206 end if
207
208 time%t = time%t + time%dt
209
210 end subroutine simulation_settime
211
213 subroutine simulation_adjoint_restart(C)
214 type(adjoint_case_t), intent(inout) :: c
215 integer :: i
216 type(file_t) :: chkpf, previous_meshf
217 character(len=LOG_SIZE) :: log_buf
218 character(len=:), allocatable :: restart_file
219 character(len=:), allocatable :: restart_mesh_file
220 real(kind=rp) :: tol
221 logical :: found
222
223 call c%case%params%get('case.restart_file', restart_file, found)
224 call c%case%params%get('case.restart_mesh_file', restart_mesh_file, found)
225
226 if (found) then
227 call previous_meshf%init(trim(restart_mesh_file))
228 call previous_meshf%read(c%fluid_adj%chkp%previous_mesh)
229 end if
230
231 call c%case%params%get('case.mesh2mesh_tolerance', tol, found)
232
233 if (found) c%case%fluid%chkp%mesh2mesh_tol = tol
234
235 call chkpf%init(trim(restart_file))
236 call chkpf%read(c%fluid_adj%chkp)
237 c%time%dtlag = c%fluid_adj%chkp%dtlag
238 c%time%tlag = c%fluid_adj%chkp%tlag
239
240 ! Free the previous mesh, dont need it anymore
241 do i = 1, size(c%time%dtlag)
242 call c%case%fluid%ext_bdf%set_coeffs(c%time%dtlag)
243 end do
244
245 call c%fluid_adj%restart(c%case%chkp)
246 call c%case%fluid%chkp%previous_mesh%free()
247 if (allocated(c%adjoint_scalars)) then
248 call c%adjoint_scalars%restart(c%case%chkp)
249 end if
250
251 c%time%t = c%case%fluid%chkp%restart_time()
252 call neko_log%section('Restarting from checkpoint')
253 write(log_buf, '(A,A)') 'File : ', trim(restart_file)
254 call neko_log%message(log_buf)
255 write(log_buf, '(A,E15.7)') 'Time : ', c%time%t
256 call neko_log%message(log_buf)
257 call neko_log%end_section()
258
259 call c%output_controller%set_counter(c%time)
260 end subroutine simulation_adjoint_restart
261
263 subroutine simulation_adjoint_joblimit_chkp(C, t)
264 type(adjoint_case_t), intent(inout) :: c
265 real(kind=rp), intent(inout) :: t
266 type(file_t) :: chkpf
267 character(len=:), allocatable :: chkp_format
268 character(len=LOG_SIZE) :: log_buf
269 character(len=10) :: format_str
270 logical :: found
271
272 call c%case%params%get('case.checkpoint_format', chkp_format, found)
273 call c%case%fluid%chkp%sync_host()
274 format_str = '.chkp'
275 if (found) then
276 if (chkp_format .eq. 'hdf5') then
277 format_str = '.h5'
278 end if
279 end if
280 call chkpf%init(c%case%output_directory // 'joblimit' // trim(format_str))
281 call chkpf%write(c%case%fluid%chkp, t)
282 write(log_buf, '(A)') '! saving checkpoint >>>'
283 call neko_log%message(log_buf)
284
285 end subroutine simulation_adjoint_joblimit_chkp
286
287end module simulation_adjoint
Adjoint simulation driver.
Adjoint case type. Todo: This should Ideally be a subclass of case_t, however, that is not yet suppor...