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 time_based_controller, only: time_based_controller_t
49 use file, only: file_t
50 use json_module, only: json_file
51 use json_utils, only: json_get, json_get_or_default, json_extract_item, &
52 json_get_or_lookup
55 use logger, only : neko_log
56 use time_state, only : time_state_t
57 use utils, only: neko_error
60 use json_utils_ext, only: json_key_fallback
62 implicit none
63 private
64 public :: adjoint_case_t
65
71 class(adjoint_fluid_scheme_t), allocatable :: fluid_adj
73 type(adjoint_scalars_t), allocatable :: adjoint_scalars
76 adjoint_convection_term
77 type(case_t), pointer :: case
78 type(time_state_t) :: time
79 type(chkp_t) :: chkp
80 type(chkp_output_t) :: chkp_out
81
82 ! Fields
83 type(adjoint_output_t) :: f_out
84 type(output_controller_t) :: output_controller
85 type(time_based_controller_t) :: norm_output_ctrl
86 type(file_t) :: norm_output_file
87 logical :: norm_output_enabled = .false.
88
89 logical :: have_scalar = .false.
90
91 contains
92 procedure, pass(this) :: init => adjoint_init_from_json
93 procedure, pass(this) :: free => adjoint_free
94 end type adjoint_case_t
95
96contains
97
98 ! Constructor from json.
99 subroutine adjoint_init_from_json(this, neko_case)
100 class(adjoint_case_t), intent(inout) :: this
101 type(case_t), target, intent(inout) :: neko_case
102
103 this%case => neko_case
104 call adjoint_case_init_common(this, neko_case)
105
106 end subroutine adjoint_init_from_json
107
109 subroutine adjoint_case_init_common(this, neko_case)
110 class(adjoint_case_t), target, intent(inout) :: this
111 type(case_t), intent(inout) :: neko_case
112 integer :: lx = 0
113 real(kind=rp) :: real_val = 0.0_rp
114 character(len=:), allocatable :: string_val, file_format, name
115 character(len=:), allocatable :: norm_control, norm_file
116 integer :: precision, integer_val, layout
117 integer :: n_scalars_primal, n_scalars_adjoint, i
118 logical :: scalar = .false.
119 logical :: temperature_found = .false.
120 logical :: logical_val
121
122 ! extra things for json
123 type(json_file) :: ic_json, numerics_params
124 type(json_file) :: scalar_params_primal, scalar_params_adjoint, json_subdict
125 character(len=:), allocatable :: json_key
126 logical :: dealias_adjoint_scalar_convection
127
128 !
129 ! Setup adjoint fluid
130 !
131 call json_get(neko_case%params, 'case.fluid.scheme', string_val)
132 call adjoint_fluid_scheme_factory(this%fluid_adj, trim(string_val))
133
134 call json_get(neko_case%params, 'case.numerics.polynomial_order', lx)
135 lx = lx + 1 ! add 1 to get number of gll points
136
137 call this%chkp%init()
138 this%chkp%tlag => this%time%tlag
139 this%chkp%dtlag => this%time%dtlag
140
141 select type (f => this%fluid_adj)
142 type is (adjoint_fluid_pnpn_t)
143 call f%init(neko_case%msh, lx, neko_case%params, &
144 neko_case%user, this%chkp)
145 end select
146 !
147 ! Setup adjoint scalar
148 !
149 ! @todo no adjoint_scalars factory for now, probably not needed
150
151 ! check how many adjoint scalars
152 scalar = .false.
153 n_scalars_adjoint = 0
154 if (neko_case%params%valid_path('case.adjoint_scalar')) then
155 call json_get_or_default(neko_case%params, &
156 'case.adjoint_scalar.enabled', scalar, .true.)
157 n_scalars_adjoint = 1
158 n_scalars_primal = 1
159 else if (neko_case%params%valid_path('case.adjoint_scalars')) then
160 call neko_case%params%info('case.adjoint_scalars', &
161 n_children = n_scalars_adjoint)
162 call neko_case%params%info('case.scalars', n_children = n_scalars_primal)
163 if (n_scalars_adjoint > 0) then
164 scalar = .true.
165 end if
166 end if
167
168 this%have_scalar = scalar
169
170
171
172
173
174
175 if (this%have_scalar) then
176 allocate(this%adjoint_scalars)
177 call json_get(neko_case%params, 'case.numerics', &
178 numerics_params)
179 if (neko_case%params%valid_path('case.adjoint_scalar')) then
180 ! For backward compatibility
181 call json_get(neko_case%params, 'case.adjoint_scalar', &
182 scalar_params_adjoint)
183 call json_get(neko_case%params, 'case.scalar', &
184 scalar_params_primal)
185 call this%adjoint_scalars%init(neko_case%msh, neko_case%fluid%c_Xh, &
186 neko_case%fluid%gs_Xh, scalar_params_adjoint, &
187 scalar_params_primal, numerics_params, neko_case%user, &
188 neko_case%chkp, neko_case%fluid%ulag, neko_case%fluid%vlag, &
189 neko_case%fluid%wlag, neko_case%fluid%ext_bdf, &
190 neko_case%fluid%rho)
191 ! allocate the coupling term
192 allocate(this%adjoint_convection_term)
193 ! initialize the coupling term
194 call json_get_or_default(neko_case%params, &
195 'case.adjoint_scalar.dealias_coupling_term', &
196 dealias_adjoint_scalar_convection, .true.)
197 call this%adjoint_convection_term%init_from_components( &
198 this%fluid_adj%f_adj_x, this%fluid_adj%f_adj_y, &
199 this%fluid_adj%f_adj_z, this%case%scalars%scalar_fields(1)%scalar%s, &
200 this%adjoint_scalars%adjoint_scalar_fields(1)%s_adj, &
201 this%fluid_adj%c_Xh, this%fluid_adj%c_Xh_GL, &
202 this%fluid_adj%GLL_to_GL, &
203 dealias_adjoint_scalar_convection, this%fluid_adj%scratch_GL)
204
205 select type (f => this%fluid_adj)
206 type is (adjoint_fluid_pnpn_t)
207 ! append the coupling term to the adjoint velocity equation
208 call f%source_term%add(this%adjoint_convection_term)
209 end select
210 else
211 ! Multiple scalars
212
213 call json_get(this%case%params, &
214 'case.adjoint_scalars', scalar_params_adjoint)
215 call json_get(this%case%params, &
216 'case.scalars', scalar_params_primal)
217 call this%adjoint_scalars%init(n_scalars_adjoint, n_scalars_primal, &
218 neko_case%msh, neko_case%fluid%c_Xh, neko_case%fluid%gs_Xh, &
219 scalar_params_adjoint, scalar_params_primal, numerics_params, &
220 neko_case%user, neko_case%chkp, neko_case%fluid%ulag, &
221 neko_case%fluid%vlag, neko_case%fluid%wlag, &
222 neko_case%fluid%ext_bdf, neko_case%fluid%rho)
223 call neko_error('The adjoint scaling coupling term have not yet' // &
224 'been implemented for multiple scalars')
225 end if
226 end if
227
228 !
229 ! Time control
230 !
231 call json_get(this%case%params, 'case.time', json_subdict)
232 call this%time%init(json_subdict)
233
234 !
235 ! Setup user defined conditions
236 !
237 ! if (neko_case%params%valid_path('case.fluid.inflow_condition')) then
238 ! call json_get(neko_case%params, 'case.fluid.inflow_condition.type', &
239 ! string_val)
240 ! if (trim(string_val) .eq. 'user') then
241 ! call neko_case%fluid%set_usr_inflow(neko_case%user%fluid_user_if)
242 ! end if
243 ! end if
244
245 ! Setup user boundary conditions for the scalar.
246 ! if (adjoint_scalars) then
247 ! call neko_case%adjoint_scalars%set_user_bc(&
248 ! neko_case%user%scalar_user_bc)
249 ! end if
250
251 !
252 ! Setup initial conditions
253 !
254
255 call neko_log%section("Adjoint initial condition")
256 json_key = json_key_fallback(neko_case%params, &
257 'case.adjoint_fluid.initial_condition', 'case.fluid.initial_condition')
258
259 call json_get(neko_case%params, json_key, ic_json)
260 call json_get(ic_json, 'type', string_val)
261
262 if (trim(string_val) .ne. 'user') then
263 call set_flow_ic( &
264 this%fluid_adj%u_adj, this%fluid_adj%v_adj, this%fluid_adj%w_adj, &
265 this%fluid_adj%p_adj, this%fluid_adj%c_Xh, this%fluid_adj%gs_Xh, &
266 string_val, ic_json)
267 else
268 call set_flow_ic( &
269 this%fluid_adj%u_adj, this%fluid_adj%v_adj, this%fluid_adj%w_adj, &
270 this%fluid_adj%p_adj, this%fluid_adj%c_Xh, this%fluid_adj%gs_Xh, &
271 neko_case%user%initial_conditions, neko_case%fluid%name)
272 end if
273
274 call neko_log%end_section()
275
276 if (this%have_scalar) then
277
278 if (neko_case%params%valid_path('case.adjoint_scalar')) then
279 ! we shouldn't fallback to the primal here.
280 call json_get(neko_case%params, &
281 'case.adjoint_scalar.initial_condition.type', string_val)
282 call json_get(neko_case%params, &
283 'case.adjoint_scalar.initial_condition', ic_json)
284
285 !call neko_log%section("Adjoint scalar initial condition ")
286
287 if (trim(string_val) .ne. 'user') then
288 if (trim(neko_case%scalars%scalar_fields(1)%scalar%name) .eq. &
289 'temperature') then
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, 0)
295 else
296 call set_scalar_ic(&
297 this%adjoint_scalars%adjoint_scalar_fields(1)%s_adj, &
298 this%adjoint_scalars%adjoint_scalar_fields(1)%c_Xh, &
299 this%adjoint_scalars%adjoint_scalar_fields(1)%gs_Xh, &
300 string_val, ic_json, 1)
301 end if
302 else
303 call neko_error("user ICs not implemented for adjoint scalar")
304 ! call set_scalar_ic(this%adjoint_scalars%s_adj, &
305 ! this%adjoint_scalars%c_Xh, this%adjoint_scalars%gs_Xh, &
306 ! this%usr%scalar_user_ic, neko_case%params)
307 end if
308
309 ! call neko_log%end_section()
310 else
311
312 ! Handle multiple scalars
313 do i = 1, n_scalars_adjoint
314 call json_extract_item(neko_case%params, 'case.adjoint_scalars', &
315 i, scalar_params_adjoint)
316 call json_get(scalar_params_adjoint, &
317 'initial_condition.type', string_val)
318 call json_get(scalar_params_adjoint, &
319 'initial_condition', json_subdict)
320
321 if (trim(string_val) .ne. 'user') then
322 if (trim(neko_case%scalars%scalar_fields(i)%scalar%name) .eq. &
323 'temperature') then
324 call set_scalar_ic(&
325 this%adjoint_scalars%adjoint_scalar_fields(i)%s_adj, &
326 this%adjoint_scalars%adjoint_scalar_fields(i)%c_Xh, &
327 this%adjoint_scalars%adjoint_scalar_fields(i)%gs_Xh, &
328 string_val, json_subdict, 0)
329 temperature_found = .true.
330 else
331 if (temperature_found) then
332 ! if temperature is found, scalars start from index 1
333 call set_scalar_ic(&
334 this%adjoint_scalars%adjoint_scalar_fields(i)%s_adj, &
335 this%adjoint_scalars%adjoint_scalar_fields(i)%c_Xh, &
336 this%adjoint_scalars%adjoint_scalar_fields(i)%gs_Xh, &
337 string_val, json_subdict, i - 1)
338 else
339 ! if temperature is not found, scalars start from index 0
340 call set_scalar_ic(&
341 this%adjoint_scalars%adjoint_scalar_fields(i)%s_adj, &
342 this%adjoint_scalars%adjoint_scalar_fields(i)%c_Xh, &
343 this%adjoint_scalars%adjoint_scalar_fields(i)%gs_Xh, &
344 string_val, json_subdict, i)
345 end if
346 end if
347 else
348 call neko_error("user ICs not implemented for adjoint scalar")
349 end if
350 end do
351 end if
352 end if
353
354 ! Add initial conditions to BDF fluid_adj (if present)
355 select type (f => this%fluid_adj)
356 type is (adjoint_fluid_pnpn_t)
357 call f%ulag%set(f%u_adj)
358 call f%vlag%set(f%v_adj)
359 call f%wlag%set(f%w_adj)
360 end select
361
362 !
363 ! Validate that the neko_case is properly setup for time-stepping
364 !
365 call this%fluid_adj%validate
366
367 if (this%have_scalar) then
368 call this%adjoint_scalars%validate()
369 end if
370
371 !
372 ! Setup output precision of the field files
373 !
374 call json_get_or_default(neko_case%params, 'case.output_precision', &
375 string_val, 'single')
376
377 if (trim(string_val) .eq. 'double') then
378 precision = dp
379 else
380 precision = sp
381 end if
382
383 !
384 ! Setup output_controller
385 !
386 call json_get_or_default(neko_case%params, &
387 'case.adjoint_fluid.output_filename', name, "adjoint")
388 call json_get_or_default(neko_case%params, &
389 'case.adjoint_fluid.output_format', file_format, 'fld')
390 call json_get_or_default(neko_case%params, &
391 'case.adjoint_fluid.output_mesh_in_all_files', &
392 logical_val, .false.)
393 call this%output_controller%init(this%time%end_time)
394 if (scalar) then
395 call this%f_out%init(precision, this%fluid_adj, &
396 this%adjoint_scalars, name = name, &
397 path = trim(neko_case%output_directory), &
398 fmt = trim(file_format), layout = layout, &
399 always_write_mesh = logical_val)
400 else
401 call this%f_out%init(precision, this%fluid_adj, name = name, &
402 path = trim(neko_case%output_directory), &
403 fmt = trim(file_format), layout = layout, &
404 always_write_mesh = logical_val)
405 end if
406
407 call json_get_or_default(neko_case%params, &
408 'case.adjoint_fluid.output_subdivide', logical_val, .false.)
409 call this%f_out%file_%set_subdivide(logical_val)
410
411 call json_get_or_default(neko_case%params, &
412 'case.adjoint_fluid.output_control', string_val, 'never')
413
414 if (trim(string_val) .eq. 'org') then
415 ! yes, it should be real_val below for type compatibility
416 call json_get_or_lookup(neko_case%params, 'case.nsamples', integer_val)
417 real_val = real(integer_val, kind=rp)
418 call this%output_controller%add(this%f_out, real_val, 'nsamples')
419 else if (trim(string_val) .eq. 'never') then
420 call this%output_controller%add(this%f_out, 0.0_rp, 'never')
421 else if (trim(string_val) .eq. 'tsteps' .or. &
422 trim(string_val) .eq. 'nsamples') then
423 call json_get_or_lookup(neko_case%params, &
424 'case.adjoint_fluid.output_value', integer_val)
425 real_val = real(integer_val, kind=rp)
426 call this%output_controller%add(this%f_out, real_val, string_val)
427 else if (trim(string_val) .eq. 'simulationtime') then
428 call json_get_or_lookup(neko_case%params, &
429 'case.adjoint_fluid.output_value', real_val)
430 call this%output_controller%add(this%f_out, real_val, string_val)
431 else
432 call neko_log%error('Unknown output control type for the fluid: ' // &
433 trim(string_val))
434 end if
435
436 !
437 ! Setup adjoint norm output
438 !
439 call json_get_or_default(neko_case%params, &
440 'case.adjoint_fluid.norm_output_control', norm_control, 'never')
441 if (trim(norm_control) .ne. 'never') then
442 call json_get_or_default(neko_case%params, &
443 'case.adjoint_fluid.norm_output_value', real_val, 1.0_rp)
444 call json_get_or_default(neko_case%params, &
445 'case.adjoint_fluid.norm_output_file', norm_file, &
446 'adjoint_norm.csv')
447 call this%norm_output_file%init(trim(neko_case%output_directory) // &
448 trim(norm_file))
449 call this%norm_output_file%set_header('Time, Norm')
450 call this%norm_output_file%set_overwrite(.true.)
451 call this%norm_output_ctrl%init(this%time%start_time, &
452 this%time%end_time, trim(norm_control), real_val)
453 this%norm_output_enabled = .true.
454 end if
455
456 ! !
457 ! ! Save checkpoints (if nothing specified, default to saving at end of sim)
458 ! !
459 ! call json_get_or_default(neko_case%params, 'case.output_checkpoints',&
460 ! logical_val, .true.)
461 ! if (logical_val) then
462 ! call json_get_or_default(neko_case%params, 'case.checkpoint_format', &
463 ! string_val, "chkp")
464 ! neko_case%f_chkp = chkp_output_t(this%fluid_adj%chkp, &
465 ! path = output_directory, &
466 ! ! fmt = trim(string_val))
467 ! call json_get_or_default(neko_case%params, 'case.checkpoint_control', &
468 ! string_val, "simulationtime")
469 ! call json_get_or_default(neko_case%params, 'case.checkpoint_value', &
470 ! real_val,&
471 ! 1e10_rp)
472 ! call this%output_controller%add(this%f_chkp, real_val, string_val)
473 ! end if
474
475 !
476 ! Initialize time and step
477 !
478 this%time%t = 0d0
479 this%time%tstep = 0
480
481 end subroutine adjoint_case_init_common
482
483 ! Destructor.
484 subroutine adjoint_free(this)
485 class(adjoint_case_t), intent(inout) :: this
486
487 if (allocated(this%fluid_adj)) then
488 call this%fluid_adj%free()
489 deallocate(this%fluid_adj)
490 end if
491
492 if (allocated(this%adjoint_scalars)) then
493 call this%adjoint_scalars%free()
494 deallocate(this%adjoint_scalars)
495 end if
496
497 if (allocated(this%adjoint_convection_term)) then
498 call this%adjoint_convection_term%free()
499 deallocate(this%adjoint_convection_term)
500 end if
501
502 nullify(this%case)
503
504 ! call this%time%free()
505 call this%chkp%free()
506 call this%chkp_out%free()
507
508 ! Fields
509 call this%f_out%free()
510 call this%output_controller%free()
511 call this%norm_output_ctrl%free()
512 call this%norm_output_file%free()
513
514 this%have_scalar = .false.
515
516 end subroutine adjoint_free
517
518end module adjoint_case
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.