Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_case.f90
1! Copyright (c) 2024, 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
33! Implements the `adjoint_case_t` type.
34module adjoint_case
35 use num_types, only: rp, dp, sp
36 use case, only: case_t
37 use adjoint_fluid_scheme, only: adjoint_fluid_scheme_t
38 use adjoint_fluid_fctry, only: adjoint_fluid_scheme_factory
41 use scalar_ic, only: set_scalar_ic
42 use flow_ic, only: set_flow_ic
43 use output_controller, only: output_controller_t
44 use file, only: file_t
45 use json_module, only: json_file
46 use json_utils, only: json_get, json_get_or_default, json_extract_object, &
47 json_extract_item
50 use logger, only : neko_log
51 use time_state, only : time_state_t
52 use utils, only: neko_error
55 use json_utils_ext, only: json_key_fallback
57 implicit none
58 private
59 public :: adjoint_case_t, adjoint_init, adjoint_free
60
66 class(adjoint_fluid_scheme_t), allocatable :: fluid_adj
68 type(adjoint_scalars_t), allocatable :: adjoint_scalars
71 adjoint_convection_term
72 type(case_t), pointer :: case
73 type(time_state_t) :: time
74
75 ! Fields
76 real(kind=rp) :: tol
77 type(adjoint_output_t) :: f_out
78 type(output_controller_t) :: output_controller
79
80 logical :: have_scalar = .false.
81
82 end type adjoint_case_t
83
84 interface adjoint_init
85 module procedure adjoint_init_from_json ! todo, init from file
86 end interface adjoint_init
87
88contains
89
90 ! Constructor from json.
91 subroutine adjoint_init_from_json(this, neko_case)
92 class(adjoint_case_t), intent(inout) :: this
93 type(case_t), target, intent(inout) :: neko_case
94
95 this%case => neko_case
96 call adjoint_case_init_common(this, neko_case)
97
98 end subroutine adjoint_init_from_json
99
101 subroutine adjoint_case_init_common(this, neko_case)
102 class(adjoint_case_t), intent(inout) :: this
103 type(case_t), intent(inout) :: neko_case
104 integer :: lx = 0
105 real(kind=rp) :: real_val = 0.0_rp
106 character(len=:), allocatable :: string_val
107 integer :: precision
108 integer :: n_scalars_primal, n_scalars_adjoint, i
109 logical :: scalar = .false.
110
111 ! extra things for json
112 type(json_file) :: ic_json, numerics_params
113 type(json_file) :: scalar_params_primal, scalar_params_adjoint, json_subdict
114 character(len=:), allocatable :: json_key
115
116 !
117 ! Setup adjoint fluid
118 !
119 call json_get(neko_case%params, 'case.fluid.scheme', string_val)
120 call adjoint_fluid_scheme_factory(this%fluid_adj, trim(string_val))
121
122 call json_get(neko_case%params, 'case.numerics.polynomial_order', lx)
123 lx = lx + 1 ! add 1 to get number of gll points
124
125 select type (f => this%fluid_adj)
126 type is (adjoint_fluid_pnpn_t)
127 call f%init(neko_case%msh, lx, neko_case%params, &
128 neko_case%user, neko_case%chkp)
129 end select
130 !
131 ! Setup adjoint scalar
132 !
133 ! @todo no adjoint_scalars factory for now, probably not needed
134
135 ! check how many adjoint scalars
136 scalar = .false.
137 n_scalars_adjoint = 0
138 if (neko_case%params%valid_path('case.adjoint_scalar')) then
139 call json_get_or_default(neko_case%params, &
140 'case.adjoint_scalar.enabled', scalar, .true.)
141 n_scalars_adjoint = 1
142 n_scalars_primal = 1
143 else if (neko_case%params%valid_path('case.adjoint_scalars')) then
144 call neko_case%params%info('case.adjoint_scalars', &
145 n_children = n_scalars_adjoint)
146 call neko_case%params%info('case.scalars', n_children = n_scalars_primal)
147 if (n_scalars_adjoint > 0) then
148 scalar = .true.
149 end if
150 end if
151
152 this%have_scalar = scalar
153
154
155
156
157
158
159 if (this%have_scalar) then
160 allocate(this%adjoint_scalars)
161 call json_extract_object(neko_case%params, 'case.numerics', &
162 numerics_params)
163 if (neko_case%params%valid_path('case.adjoint_scalar')) then
164 ! For backward compatibility
165 call json_extract_object(neko_case%params, 'case.adjoint_scalar', &
166 scalar_params_adjoint)
167 call json_extract_object(neko_case%params, 'case.scalar', &
168 scalar_params_primal)
169 call this%adjoint_scalars%init(neko_case%msh, neko_case%fluid%c_Xh, &
170 neko_case%fluid%gs_Xh, scalar_params_adjoint, &
171 scalar_params_primal, numerics_params, neko_case%user, &
172 neko_case%chkp, neko_case%fluid%ulag, neko_case%fluid%vlag, &
173 neko_case%fluid%wlag, neko_case%fluid%ext_bdf, &
174 neko_case%fluid%rho)
175 ! allocate the coupling term
176 allocate(this%adjoint_convection_term)
177 ! initialize the coupling term
178 call this%adjoint_convection_term%init_from_components( &
179 this%fluid_adj%f_adj_x, this%fluid_adj%f_adj_y, &
180 this%fluid_adj%f_adj_z, this%case%scalars%scalar_fields(1)%s, &
181 this%adjoint_scalars%adjoint_scalar_fields(1)%s_adj, &
182 this%fluid_adj%c_Xh)
183
184 select type (f => this%fluid_adj)
185 type is (adjoint_fluid_pnpn_t)
186 ! append the coupling term to the adjoint velocity equation
187 call f%source_term%add(this%adjoint_convection_term)
188 end select
189 else
190 ! Multiple scalars
191
192 call json_extract_object(this%case%params, &
193 'case.adjoint_scalars', scalar_params_adjoint)
194 call json_extract_object(this%case%params, &
195 'case.scalars', scalar_params_primal)
196 call this%adjoint_scalars%init(n_scalars_adjoint, n_scalars_primal, &
197 neko_case%msh, neko_case%fluid%c_Xh, neko_case%fluid%gs_Xh, &
198 scalar_params_adjoint, scalar_params_primal, numerics_params, &
199 neko_case%user, neko_case%chkp, neko_case%fluid%ulag, &
200 neko_case%fluid%vlag, neko_case%fluid%wlag, &
201 neko_case%fluid%ext_bdf, neko_case%fluid%rho)
202 call neko_error('The adjoint scaling coupling term have not yet' // &
203 'been implemented for multiple scalars')
204 end if
205 end if
206
207 !
208 ! Time control
209 !
210 call json_extract_object(this%case%params, 'case.time', json_subdict)
211 call this%time%init(json_subdict)
212
213 !
214 ! Setup user defined conditions
215 !
216 ! if (neko_case%params%valid_path('case.fluid.inflow_condition')) then
217 ! call json_get(neko_case%params, 'case.fluid.inflow_condition.type', &
218 ! string_val)
219 ! if (trim(string_val) .eq. 'user') then
220 ! call neko_case%fluid%set_usr_inflow(neko_case%user%fluid_user_if)
221 ! end if
222 ! end if
223
224 ! Setup user boundary conditions for the scalar.
225 ! if (adjoint_scalars) then
226 ! call neko_case%adjoint_scalars%set_user_bc(&
227 ! neko_case%user%scalar_user_bc)
228 ! end if
229
230 !
231 ! Setup initial conditions
232 !
233
234 call neko_log%section("Adjoint initial condition")
235 json_key = json_key_fallback(neko_case%params, &
236 'case.adjoint_fluid.initial_condition', 'case.fluid.initial_condition')
237
238 call json_extract_object(neko_case%params, json_key, ic_json)
239 call json_get(ic_json, 'type', string_val)
240
241 if (trim(string_val) .ne. 'user') then
242 call set_flow_ic( &
243 this%fluid_adj%u_adj, this%fluid_adj%v_adj, this%fluid_adj%w_adj, &
244 this%fluid_adj%p_adj, this%fluid_adj%c_Xh, this%fluid_adj%gs_Xh, &
245 string_val, ic_json)
246 else
247 call set_flow_ic( &
248 this%fluid_adj%u_adj, this%fluid_adj%v_adj, this%fluid_adj%w_adj, &
249 this%fluid_adj%p_adj, this%fluid_adj%c_Xh, this%fluid_adj%gs_Xh, &
250 neko_case%user%fluid_user_ic, neko_case%params)
251 end if
252
253 call neko_log%end_section()
254
255 if (this%have_scalar) then
256
257 if (neko_case%params%valid_path('case.adjoint_scalar')) then
258 ! we shouldn't fallback to the primal here.
259 call json_get(neko_case%params, &
260 'case.adjoint_scalar.initial_condition.type', string_val)
261 call json_extract_object(neko_case%params, &
262 'case.adjoint_scalar.initial_condition', ic_json)
263
264 !call neko_log%section("Adjoint scalar initial condition ")
265
266 if (trim(string_val) .ne. 'user') then
267 call set_scalar_ic(&
268 this%adjoint_scalars%adjoint_scalar_fields(1)%s_adj, &
269 this%adjoint_scalars%adjoint_scalar_fields(1)%c_Xh, &
270 this%adjoint_scalars%adjoint_scalar_fields(1)%gs_Xh, &
271 string_val, ic_json)
272 else
273 call neko_error("user ICs not implemented for adjoint scalar")
274 ! call set_scalar_ic(this%adjoint_scalars%s_adj, &
275 ! this%adjoint_scalars%c_Xh, this%adjoint_scalars%gs_Xh, &
276 ! this%usr%scalar_user_ic, neko_case%params)
277 end if
278
279 ! call neko_log%end_section()
280 else
281
282 ! Handle multiple scalars
283 do i = 1, n_scalars_adjoint
284 call json_extract_item(neko_case%params, 'case.adjoint_scalars', &
285 i, scalar_params_adjoint)
286 call json_get(scalar_params_adjoint, &
287 'initial_condition.type', string_val)
288 call json_extract_object(scalar_params_adjoint, &
289 'initial_condition', json_subdict)
290
291 if (trim(string_val) .ne. 'user') then
292 call set_scalar_ic(&
293 this%adjoint_scalars%adjoint_scalar_fields(i)%s_adj, &
294 this%adjoint_scalars%adjoint_scalar_fields(i)%c_Xh, &
295 this%adjoint_scalars%adjoint_scalar_fields(i)%gs_Xh, &
296 string_val, json_subdict)
297 else
298 call neko_error("user ICs not implemented for adjoint scalar")
299 end if
300 end do
301 end if
302 end if
303
304 ! Add initial conditions to BDF fluid_adj (if present)
305 select type (f => this%fluid_adj)
306 type is (adjoint_fluid_pnpn_t)
307 call f%ulag%set(f%u_adj)
308 call f%vlag%set(f%v_adj)
309 call f%wlag%set(f%w_adj)
310 end select
311
312 !
313 ! Validate that the neko_case is properly setup for time-stepping
314 !
315 call this%fluid_adj%validate
316
317 if (this%have_scalar) then
318 call this%adjoint_scalars%validate()
319 end if
320
321 !
322 ! Setup output precision of the field files
323 !
324 call json_get_or_default(neko_case%params, 'case.output_precision', &
325 string_val, 'single')
326
327 if (trim(string_val) .eq. 'double') then
328 precision = dp
329 else
330 precision = sp
331 end if
332
333 !
334 ! Setup output_controller
335 !
336 call this%output_controller%init(neko_case%time%end_time)
337 if (this%have_scalar) then
338 this%f_out = adjoint_output_t(precision, this%fluid_adj, &
339 this%adjoint_scalars, path = trim(neko_case%output_directory))
340 else
341 this%f_out = adjoint_output_t(precision, this%fluid_adj, &
342 path = trim(neko_case%output_directory))
343 end if
344
345 call json_get_or_default(neko_case%params, 'case.fluid.output_control',&
346 string_val, 'org')
347
348 if (trim(string_val) .eq. 'org') then
349 ! yes, it should be real_val below for type compatibility
350 call json_get(neko_case%params, 'case.nsamples', real_val)
351 call this%output_controller%add(this%f_out, real_val, 'nsamples')
352 else if (trim(string_val) .eq. 'never') then
353 ! Fix a dummy 0.0 output_value
354 call json_get_or_default(neko_case%params, 'case.fluid.output_value', &
355 real_val, 0.0_rp)
356 call this%output_controller%add(this%f_out, 0.0_rp, string_val)
357 else
358 call json_get(neko_case%params, 'case.fluid.output_value', real_val)
359 call this%output_controller%add(this%f_out, real_val, string_val)
360 end if
361
362 ! !
363 ! ! Save checkpoints (if nothing specified, default to saving at end of sim)
364 ! !
365 ! call json_get_or_default(neko_case%params, 'case.output_checkpoints',&
366 ! logical_val, .true.)
367 ! if (logical_val) then
368 ! call json_get_or_default(neko_case%params, 'case.checkpoint_format', &
369 ! string_val, "chkp")
370 ! neko_case%f_chkp = chkp_output_t(this%fluid_adj%chkp, &
371 ! path = output_directory, &
372 ! ! fmt = trim(string_val))
373 ! call json_get_or_default(neko_case%params, 'case.checkpoint_control', &
374 ! string_val, "simulationtime")
375 ! call json_get_or_default(neko_case%params, 'case.checkpoint_value', &
376 ! real_val,&
377 ! 1e10_rp)
378 ! call this%output_controller%add(neko_case%f_chkp, real_val, string_val)
379 ! end if
380
381 !
382 ! Initialize time and step
383 !
384 this%time%t = 0d0
385 this%time%tstep = 0
386
387 end subroutine adjoint_case_init_common
388
389 ! Destructor.
390 subroutine adjoint_free(this)
391 class(adjoint_case_t), intent(inout) :: this
392
393 nullify(this%case)
394 if (allocated(this%adjoint_scalars)) then
395 call this%adjoint_scalars%free()
396 deallocate(this%adjoint_scalars)
397 end if
398
399 if (allocated(this%fluid_adj)) then
400 call this%fluid_adj%free()
401 deallocate(this%fluid_adj)
402 end if
403 call this%output_controller%free()
404
405 end subroutine adjoint_free
406
407end module adjoint_case
408
Factory for all adjoint fluid schemes.
Adjoint Pn/Pn formulation.
Defines an output for a adjoint.
Implements the adjoint_scalar_convection_source_term type.
Contains the adjoint_scalar_pnpn_t type.
Contains the adjoint_scalar_scheme_t type.
Contains the adjoint_scalars_t type that manages multiple scalar fields.
Adjoint case type. Todo: This should Ideally be a subclass of case_t, however, that is not yet suppor...
Base type for a scalar advection-diffusion solver.
Type to manage multiple adjoint scalar transport equations.