Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
simulation_adjoint.f90
Go to the documentation of this file.
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_step_controller, only: time_step_controller_t
48 implicit none
49 private
50
52
53contains
54
55 ! Compute the simcomp_test field.
56 subroutine solve_adjoint(this)
57 type(adjoint_case_t), intent(inout) :: this
58 real(kind=rp) :: t
59 integer :: tstep
60
61 real(kind=rp) :: t_adj
62 real(kind=dp) :: start_time_org, start_time, end_time
63 character(len=LOG_SIZE) :: log_buf
64 integer :: tstep_adj
65 logical :: output_at_end
66 type(time_step_controller_t) :: dt_controller
67
68 ! ------------------------------------------------------------------------ !
69 ! Computation of the adjoint field.
70 !
71 ! This is where the actual computation of the adjoint field should be
72 ! implemented. This will so far just be a steady state field based on the
73 ! current state of the fluid fields.
74
75 ! ------------------------------------------------------------------------ !
76
77 t_adj = 0d0
78 tstep_adj = 0
79 call neko_log%section('Starting adjoint')
80 write(log_buf, '(A,E15.7,A,E15.7,A)') 'T : [', 0d0, ', ', &
81 this%case%end_time, ')'
82 call neko_log%message(log_buf)
83 call dt_controller%init(this%case%params)
84 if (.not. dt_controller%if_variable_dt) then
85 write(log_buf, '(A, E15.7)') 'dt : ', this%case%dt
86 call neko_log%message(log_buf)
87 else
88 write(log_buf, '(A, E15.7)') 'CFL : ', dt_controller%set_cfl
89 call neko_log%message(log_buf)
90 end if
91
93 call neko_log%section('Postprocessing')
94 ! call this%case%q%eval(t_adj, this%case%dt, tstep_adj)
95 call this%output_controller%execute(t_adj, tstep_adj)
96
97 ! HARRY
98 ! ok this I guess this is techincally where we set the initial condition
99 ! of adjoint yeh?
100 call this%case%usr%user_init_modules(t_adj, this%scheme%u_adj, &
101 this%scheme%v_adj, this%scheme%w_adj,&
102 this%scheme%p_adj, this%scheme%c_Xh, this%case%params)
103 call neko_log%end_section()
104 call neko_log%newline()
105
106 call profiler_start
107 start_time_org = mpi_wtime()
108
109 do while (t_adj .lt. this%case%end_time .and. (.not. jobctrl_time_limit()))
110 call profiler_start_region('Time-Step')
111 tstep_adj = tstep_adj + 1
112 start_time = mpi_wtime()
113
114 call neko_log%status(t_adj, this%case%end_time)
115 write(log_buf, '(A,I6)') 'Time-step: ', tstep_adj
116 call neko_log%message(log_buf)
117 call neko_log%begin()
118
119 ! write(log_buf, '(A,E15.7,1x,A,E15.7)') 'CFL:', cfl, 'dt:', this%case%dt
120 call neko_log%message(log_buf)
121 call simulation_settime(t_adj, this%case%dt, this%case%fluid%ext_bdf, &
122 this%case%tlag, this%case%dtlag, tstep_adj)
123
124 call neko_log%section('Adjoint fluid')
125 call this%scheme%step(t_adj, tstep_adj, this%case%dt, &
126 this%case%fluid%ext_bdf, dt_controller)
127 end_time = mpi_wtime()
128 write(log_buf, '(A,E15.7,A,E15.7)') &
129 'Elapsed time (s):', end_time-start_time_org, ' Step time:', &
130 end_time-start_time
131 call neko_log%end_section(log_buf)
132
133 ! Scalar step
134 if (allocated(this%case%scalar)) then
135 start_time = mpi_wtime()
136 call neko_log%section('Scalar')
137 call this%case%scalar%step(t_adj, tstep_adj, this%case%dt, &
138 this%case%fluid%ext_bdf, dt_controller)
139 end_time = mpi_wtime()
140 write(log_buf, '(A,E15.7,A,E15.7)') &
141 'Elapsed time (s):', end_time-start_time_org, ' Step time:', &
142 end_time-start_time
143 call neko_log%end_section(log_buf)
144 end if
145
146 call neko_log%section('Postprocessing')
147
148 ! call this%case%q%eval(t_adj, this%case%dt, tstep_adj)
149 call this%output_controller%execute(t_adj, tstep_adj)
150
151 ! Update material properties
152 call this%case%usr%material_properties(t, tstep, this%case%fluid%rho, &
153 this%case%fluid%mu, &
154 this%case%scalar%cp, &
155 this%case%scalar%lambda, &
156 this%case%params)
157
158 call neko_log%end_section()
159
160 call neko_log%end()
161 call profiler_end_region
162 end do
163 call profiler_stop
164
165 call json_get_or_default(this%case%params, 'case.output_at_end',&
166 output_at_end, .true.)
167 call this%output_controller%execute(t_adj, tstep_adj, output_at_end)
168
169 if (.not. (output_at_end) .and. t_adj .lt. this%case%end_time) then
170 call simulation_joblimit_chkp(this%case, t_adj)
171 end if
172
173 call this%case%usr%user_finalize_modules(t_adj, this%case%params)
174
175 call neko_log%end_section('Normal end.')
176
177 end subroutine solve_adjoint
178
179 subroutine simulation_settime(t, dt, ext_bdf, tlag, dtlag, step)
180 real(kind=rp), intent(inout) :: t
181 real(kind=rp), intent(in) :: dt
182 type(time_scheme_controller_t), intent(inout) :: ext_bdf
183 real(kind=rp), dimension(10) :: tlag
184 real(kind=rp), dimension(10) :: dtlag
185 integer, intent(in) :: step
186 integer :: i
187
188
189 do i = 10, 2, -1
190 tlag(i) = tlag(i-1)
191 dtlag(i) = dtlag(i-1)
192 end do
193
194 dtlag(1) = dt
195 tlag(1) = t
196 if (ext_bdf%ndiff .eq. 0) then
197 dtlag(2) = dt
198 tlag(2) = t
199 end if
200
201 t = t + dt
202
203 call ext_bdf%set_coeffs(dtlag)
204
205 end subroutine simulation_settime
206
208 subroutine simulation_restart(C, t)
209 type(case_t), intent(inout) :: c
210 real(kind=rp), intent(inout) :: t
211 integer :: i
212 type(file_t) :: chkpf, previous_meshf
213 character(len=LOG_SIZE) :: log_buf
214 character(len=:), allocatable :: restart_file
215 character(len=:), allocatable :: restart_mesh_file
216 real(kind=rp) :: tol
217 logical :: found
218
219 call c%params%get('case.restart_file', restart_file, found)
220 call c%params%get('case.restart_mesh_file', restart_mesh_file, found)
221
222 if (found) then
223 previous_meshf = file_t(trim(restart_mesh_file))
224 call previous_meshf%read(c%fluid%chkp%previous_mesh)
225 end if
226
227 call c%params%get('case.mesh2mesh_tolerance', tol,&
228 found)
229
230 if (found) c%fluid%chkp%mesh2mesh_tol = tol
231
232 chkpf = file_t(trim(restart_file))
233 call chkpf%read(c%fluid%chkp)
234 c%dtlag = c%fluid%chkp%dtlag
235 c%tlag = c%fluid%chkp%tlag
236
237 !Free the previous mesh, dont need it anymore
238 call c%fluid%chkp%previous_mesh%free()
239 do i = 1, size(c%dtlag)
240 call c%fluid%ext_bdf%set_coeffs(c%dtlag)
241 end do
242
243 call c%fluid%restart(c%dtlag, c%tlag)
244 if (allocated(c%scalar)) call c%scalar%restart(c%dtlag, c%tlag)
245
246 t = c%fluid%chkp%restart_time()
247 call neko_log%section('Restarting from checkpoint')
248 write(log_buf, '(A,A)') 'File : ', trim(restart_file)
249 call neko_log%message(log_buf)
250 write(log_buf, '(A,E15.7)') 'Time : ', t
251 call neko_log%message(log_buf)
252 call neko_log%end_section()
253
254 call c%output_controller%set_counter(t)
255 end subroutine simulation_restart
256
259 type(case_t), intent(inout) :: C
260 real(kind=rp), intent(inout) :: t
261 type(file_t) :: chkpf
262 character(len=:), allocatable :: chkp_format
263 character(len=LOG_SIZE) :: log_buf
264 character(len=10) :: format_str
265 logical :: found
266
267 call c%params%get('case.checkpoint_format', chkp_format, found)
268 call c%fluid%chkp%sync_host()
269 format_str = '.chkp'
270 if (found) then
271 if (chkp_format .eq. 'hdf5') then
272 format_str = '.h5'
273 end if
274 end if
275 chkpf = file_t('joblimit'//trim(format_str))
276 call chkpf%write(c%fluid%chkp, t)
277 write(log_buf, '(A)') '! saving checkpoint >>>'
278 call neko_log%message(log_buf)
279
280 end subroutine simulation_joblimit_chkp
281
282
283end module simulation_adjoint
284
285
Adjoint simulation driver.
subroutine simulation_settime(t, dt, ext_bdf, tlag, dtlag, step)
subroutine, public simulation_restart(c, t)
Restart a case C from a given checkpoint.
subroutine, public solve_adjoint(this)
subroutine simulation_joblimit_chkp(c, t)
Write a checkpoint at joblimit.
Adjoint case type. Todo: This should Ideally be a subclass of case_t, however, that is not yet suppoe...