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