42 use num_types,
only: rp, dp
47 use vector,
only: vector_t
48 use matrix,
only: matrix_t
49 use device,
only: host_to_device, device_to_host
50 use json_module,
only: json_file
51 use json_utils,
only: json_extract_item, json_get, json_get_or_default
53 use logger,
only: neko_log
55 use time_state,
only: time_state_t
56 use vector_math,
only: vector_add2, vector_cfill
57 use time_step_controller,
only: time_step_controller_t
60 use simulation,
only: simulation_init, simulation_step, simulation_finalize
61 use mpi_f08,
only: mpi_wtime
62 use profiler,
only: profiler_start_region, profiler_end_region
71 integer :: n_design = 0
73 integer :: n_objectives = 0
75 integer :: n_constraints = 0
89 procedure, pass(this),
public :: free => problem_free
94 procedure, pass(this),
public :: compute => problem_compute
99 procedure, pass(this),
public :: compute_sensitivity => &
100 problem_compute_sensitivity
103 procedure, pass(this),
public :: run_forward_unsteady => &
104 problem_run_forward_unsteady
107 procedure, pass(this),
public :: run_backward_unsteady => &
108 problem_run_backward_unsteady
113 procedure, pass(this),
public :: read_objectives => problem_read_objectives
115 procedure, pass(this),
public :: read_constraints => &
116 problem_read_constraints
122 procedure, pass(this),
public :: write => problem_write
124 procedure, pass(this),
public :: get_log_values => problem_get_log_values
127 procedure, pass(this),
public :: add_objective => problem_add_objective
129 procedure, pass(this),
public :: add_constraint => problem_add_constraint
135 procedure, pass(this) :: update_objectives => &
136 problem_update_objectives
138 procedure, pass(this) :: update_constraints => &
139 problem_update_constraints
141 procedure, pass(this) :: update_objective_sensitivities => &
142 problem_update_objective_sensitivities
144 procedure, pass(this) :: update_constraint_sensitivities => &
145 problem_update_constraint_sensitivities
148 procedure, pass(this) :: reset_objectives => &
149 problem_reset_objectives
151 procedure, pass(this) :: reset_constraints => &
152 problem_reset_constraints
154 procedure, pass(this) :: reset_objective_sensitivities => &
155 problem_reset_objective_sensitivities
157 procedure, pass(this) :: reset_constraint_sensitivities => &
158 problem_reset_constraint_sensitivities
161 procedure, pass(this) :: accumulate_objectives => &
162 problem_accumulate_objectives
164 procedure, pass(this) :: accumulate_constraints => &
165 problem_accumulate_constraints
167 procedure, pass(this) :: accumulate_objective_sensitivities => &
168 problem_accumulate_objective_sensitivities
170 procedure, pass(this) :: accumulate_constraint_sensitivities => &
171 problem_accumulate_constraint_sensitivities
177 procedure, pass(this),
public :: get_objective_value => &
178 problem_get_objective_value
180 procedure, pass(this),
public :: get_all_objective_values => &
181 problem_get_all_objective_values
183 procedure, pass(this),
public :: get_constraint_values => &
184 problem_get_constraint_values
186 procedure, pass(this),
public :: get_objective_sensitivities => &
187 problem_get_objective_sensitivities
189 procedure, pass(this),
public :: get_constraint_sensitivities => &
190 problem_get_constraint_sensitivities
193 procedure, pass(this) :: get_n_objectives => problem_get_num_objectives
195 procedure, pass(this) :: get_n_constraints => problem_get_num_constraints
198 procedure, pass(this) :: get_log_header => problem_get_log_header
200 procedure, pass(this) :: get_log_size => problem_get_log_size
212 type(json_file),
intent(inout) :: parameters
213 class(
design_t),
intent(in) :: design
214 type(
simulation_t),
optional,
intent(inout) :: simulation
218 this%n_design =
design%size()
221 call this%read_objectives(parameters,
design, simulation)
222 call this%read_constraints(parameters,
design, simulation)
227 subroutine problem_free(this)
232 this%n_objectives = 0
233 this%n_constraints = 0
236 if (
allocated(this%objective_list))
then
237 do i = 1,
size(this%objective_list)
238 call this%objective_list(i)%free()
240 deallocate(this%objective_list)
244 if (
allocated(this%constraint_list))
then
245 do i = 1,
size(this%constraint_list)
246 call this%constraint_list(i)%free()
248 deallocate(this%constraint_list)
250 end subroutine problem_free
253 subroutine problem_write(this, idx)
255 integer,
intent(in) :: idx
257 end subroutine problem_write
263 subroutine problem_get_log_values(this, values, include_constraints)
265 real(kind=rp),
intent(out) :: values(:)
266 logical,
intent(in),
optional :: include_constraints
267 integer :: i, n, offset
268 real(kind=rp) :: objective_value
269 real(kind=rp),
allocatable :: tmp(:)
270 logical :: do_constraints
272 if (
present(include_constraints))
then
273 do_constraints = include_constraints
275 do_constraints = .true.
278 call this%get_objective_value(objective_value)
280 values(1) = objective_value
283 do i = 1, this%n_objectives
284 n = this%objective_list(i)%objective%get_log_size()
287 call this%objective_list(i)%objective%get_log_values(tmp)
288 values(offset:offset + n - 1) = tmp
294 if (do_constraints)
then
295 do i = 1, this%n_constraints
296 n = this%constraint_list(i)%constraint%get_log_size()
299 call this%constraint_list(i)%constraint%get_log_values(tmp)
300 values(offset:offset + n - 1) = tmp
306 end subroutine problem_get_log_values
312 subroutine problem_read_objectives(this, parameters, design, simulation)
314 type(json_file),
intent(inout) :: parameters
315 class(
design_t),
intent(in) :: design
316 type(
simulation_t),
optional,
intent(inout) :: simulation
320 character(len=:),
allocatable :: path, type
321 type(json_file) :: objective_json
322 integer :: n_objectives, i
325 call neko_log%section(
"Reading objectives")
328 path =
"optimization.objectives"
329 if (parameters%valid_path(path))
then
330 call parameters%info(path, n_children = n_objectives)
333 do i = 1, n_objectives
334 call json_extract_item(parameters, path, i, objective_json)
335 call json_get(objective_json,
"type", type)
336 call neko_log%message(type)
343 if (
present(simulation))
then
348 call json_get_or_default(parameters, &
349 "adjoint_fluid.dealias_sensitivity", dealias, .true.)
350 call alo%init_from_attributes(
design, simulation, weight = 1.0_rp, &
351 name =
"Augmented Lagrangian", mask_name =
"", &
357 call neko_log%end_section()
359 end subroutine problem_read_objectives
362 subroutine problem_read_constraints(this, parameters, design, simulation)
364 type(json_file),
intent(inout) :: parameters
365 class(
design_t),
intent(in) :: design
367 type(
simulation_t),
optional,
intent(inout) :: simulation
370 character(len=:),
allocatable :: path, type
371 type(json_file) :: constraint_json
372 integer :: n_constraints, i
374 call neko_log%section(
"Reading constraints")
377 path =
"optimization.constraints"
379 if (parameters%valid_path(path))
then
380 call parameters%info(path, n_children = n_constraints)
383 do i = 1, n_constraints
384 call json_extract_item(parameters, path, i, constraint_json)
385 call json_get(constraint_json,
"type", type)
386 call neko_log%message(type)
393 call neko_log%end_section()
395 end subroutine problem_read_constraints
398 subroutine problem_add_objective(this, objective)
400 class(
objective_t),
allocatable,
intent(inout) :: objective
405 if (
allocated(this%objective_list))
then
406 n =
size(this%objective_list)
407 call move_alloc(this%objective_list, temp_list)
408 allocate(this%objective_list(n + 1))
409 if (
allocated(temp_list))
then
411 call move_alloc(temp_list(i)%objective, &
412 this%objective_list(i)%objective)
416 allocate(this%objective_list(1))
419 call move_alloc(
objective, this%objective_list(n + 1)%objective)
420 this%n_objectives = n + 1
421 end subroutine problem_add_objective
424 subroutine problem_add_constraint(this, constraint)
426 class(
constraint_t),
allocatable,
intent(inout) :: constraint
431 if (
allocated(this%constraint_list))
then
432 n =
size(this%constraint_list)
433 call move_alloc(this%constraint_list, temp_list)
434 allocate(this%constraint_list(n + 1))
435 if (
allocated(temp_list))
then
437 call move_alloc(temp_list(i)%constraint, &
438 this%constraint_list(i)%constraint)
442 allocate(this%constraint_list(1))
445 call move_alloc(
constraint, this%constraint_list(n + 1)%constraint)
446 this%n_constraints = n + 1
447 end subroutine problem_add_constraint
453 subroutine problem_compute(this, design, simulation)
455 class(
design_t),
intent(inout) :: design
456 class(
simulation_t),
optional,
intent(inout) :: simulation
458 if (
present(simulation))
then
459 call simulation%reset()
460 if (simulation%unsteady)
then
462 call this%run_forward_unsteady(simulation,
design)
464 call simulation%run_forward()
466 call this%update_objectives(
design)
469 call this%update_objectives(
design)
472 call this%update_constraints(
design)
474 end subroutine problem_compute
477 subroutine problem_compute_sensitivity(this, design, simulation)
479 class(
design_t),
intent(inout) :: design
480 class(
simulation_t),
optional,
intent(inout) :: simulation
482 type(vector_t) :: objective_sensitivity
484 if (
present(simulation))
then
485 if (simulation%unsteady)
then
487 call this%run_backward_unsteady(simulation,
design)
489 call simulation%run_backward()
491 call this%update_objective_sensitivities(
design)
494 call this%update_objective_sensitivities(
design)
497 call this%update_constraint_sensitivities(
design)
499 call objective_sensitivity%init(this%n_design)
500 call this%get_objective_sensitivities(objective_sensitivity)
502 call design%map_backward(objective_sensitivity)
504 call objective_sensitivity%free()
505 end subroutine problem_compute_sensitivity
508 subroutine problem_run_forward_unsteady(this, simulation, design)
511 class(
design_t),
intent(inout) :: design
512 type(time_step_controller_t) :: dt_controller
513 real(kind=dp) :: loop_start
515 call dt_controller%init(simulation%neko_case%params)
517 call simulation%reset()
518 call simulation_init(simulation%neko_case, dt_controller)
521 call this%reset_objectives()
523 call profiler_start_region(
"Forward simulation")
524 loop_start = mpi_wtime()
525 simulation%n_timesteps = 0
526 do while (simulation%neko_case%time%t .lt. simulation%neko_case%time%end_time)
527 simulation%n_timesteps = simulation%n_timesteps + 1
529 call simulation_step(simulation%neko_case, dt_controller, loop_start)
531 call this%accumulate_objectives(
design, simulation%neko_case%time)
533 call simulation%checkpoint%save(simulation%neko_case)
535 call profiler_end_region(
"Forward simulation")
537 call simulation_finalize(simulation%neko_case)
539 end subroutine problem_run_forward_unsteady
542 subroutine problem_run_backward_unsteady(this, simulation, design)
545 class(
design_t),
intent(inout) :: design
546 type(time_step_controller_t) :: dt_controller
547 real(kind=dp) :: loop_start
549 real(kind=rp) :: total_time
551 type(time_state_t) :: accumulation_time
553 call dt_controller%init(simulation%neko_case%params)
558 call this%reset_objective_sensitivities()
560 cfl = simulation%adjoint_case%fluid_adj%compute_cfl(simulation%adjoint_case%time%dt)
561 loop_start = mpi_wtime()
564 total_time = simulation%n_timesteps * simulation%adjoint_case%time%dt
566 call profiler_start_region(
"Adjoint simulation")
585 do i = simulation%n_timesteps, 1, -1
587 call simulation%checkpoint%restore(simulation%neko_case, i)
589 accumulation_time = simulation%adjoint_case%time
590 accumulation_time%t = total_time - simulation%adjoint_case%time%t
591 call this%accumulate_objective_sensitivities(
design, accumulation_time)
594 cfl, loop_start, total_time)
597 call profiler_end_region(
"Adjoint simulation")
601 end subroutine problem_run_backward_unsteady
612 subroutine problem_update_objectives(this, design)
614 class(
design_t),
intent(in) :: design
617 do i = 1, this%n_objectives
618 call this%objective_list(i)%objective%update_value(
design)
620 end subroutine problem_update_objectives
628 subroutine problem_update_constraints(this, design)
630 class(
design_t),
intent(in) :: design
633 do i = 1, this%n_constraints
634 call this%constraint_list(i)%constraint%update_value(
design)
636 end subroutine problem_update_constraints
644 subroutine problem_update_objective_sensitivities(this, design)
646 class(
design_t),
intent(in) :: design
649 do i = 1, this%n_objectives
650 call this%objective_list(i)%objective%update_sensitivity(
design)
652 end subroutine problem_update_objective_sensitivities
660 subroutine problem_update_constraint_sensitivities(this, design)
662 class(
design_t),
intent(in) :: design
665 do i = 1, this%n_constraints
666 call this%constraint_list(i)%constraint%update_sensitivity(
design)
668 end subroutine problem_update_constraint_sensitivities
677 subroutine problem_reset_objectives(this)
681 do i = 1, this%n_objectives
682 call this%objective_list(i)%objective%reset_value()
684 end subroutine problem_reset_objectives
690 subroutine problem_reset_constraints(this)
694 do i = 1, this%n_constraints
695 call this%constraint_list(i)%constraint%reset_value()
697 end subroutine problem_reset_constraints
703 subroutine problem_reset_objective_sensitivities(this)
707 do i = 1, this%n_objectives
708 call this%objective_list(i)%objective%reset_sensitivity()
710 end subroutine problem_reset_objective_sensitivities
716 subroutine problem_reset_constraint_sensitivities(this)
720 do i = 1, this%n_constraints
721 call this%constraint_list(i)%constraint%reset_sensitivity()
723 end subroutine problem_reset_constraint_sensitivities
734 subroutine problem_accumulate_objectives(this, design, time)
736 class(
design_t),
intent(in) :: design
737 type(time_state_t),
intent(in) :: time
740 do i = 1, this%n_objectives
741 call this%objective_list(i)%objective%accumulate_value(
design, time)
743 end subroutine problem_accumulate_objectives
751 subroutine problem_accumulate_constraints(this, design, time)
753 class(
design_t),
intent(in) :: design
754 type(time_state_t),
intent(in) :: time
757 do i = 1, this%n_constraints
758 call this%constraint_list(i)%constraint%accumulate_value(
design, time)
760 end subroutine problem_accumulate_constraints
768 subroutine problem_accumulate_objective_sensitivities(this, design, time)
770 class(
design_t),
intent(in) :: design
771 type(time_state_t),
intent(in) :: time
774 do i = 1, this%n_objectives
775 call this%objective_list(i)%objective%accumulate_sensitivity(
design, &
778 end subroutine problem_accumulate_objective_sensitivities
786 subroutine problem_accumulate_constraint_sensitivities(this, design, time)
788 class(
design_t),
intent(in) :: design
789 type(time_state_t),
intent(in) :: time
792 do i = 1, this%n_constraints
793 call this%constraint_list(i)%constraint%accumulate_sensitivity(
design, &
796 end subroutine problem_accumulate_constraint_sensitivities
807 subroutine problem_get_objective_value(this, objective_value)
809 real(kind=rp),
intent(out) :: objective_value
812 objective_value = 0.0_rp
813 do i = 1, this%n_objectives
814 objective_value = objective_value + &
815 this%objective_list(i)%objective%get_weight() * &
816 this%objective_list(i)%objective%get_value()
819 end subroutine problem_get_objective_value
827 subroutine problem_get_all_objective_values(this, all_objective_values)
829 type(vector_t),
intent(inout) :: all_objective_values
832 do i = 1, this%n_objectives
833 all_objective_values%x(i) = this%objective_list(i)%objective%value
836 call all_objective_values%copy_from(host_to_device, sync = .true.)
838 end subroutine problem_get_all_objective_values
846 subroutine problem_get_constraint_values(this, constraint_value)
848 type(vector_t),
intent(inout) :: constraint_value
851 do i = 1, this%n_constraints
852 constraint_value%x(i) = this%constraint_list(i)%constraint%value
855 call constraint_value%copy_from(host_to_device, sync = .true.)
857 end subroutine problem_get_constraint_values
865 subroutine problem_get_objective_sensitivities(this, sensitivity)
867 type(vector_t),
intent(inout) :: sensitivity
870 call vector_cfill(sensitivity, 0.0_rp)
871 do i = 1, this%n_objectives
872 call vector_add2(sensitivity, &
873 this%objective_list(i)%objective%sensitivity)
876 end subroutine problem_get_objective_sensitivities
884 subroutine problem_get_constraint_sensitivities(this, sensitivity)
886 type(matrix_t),
target,
intent(inout) :: sensitivity
887 real(kind=rp),
pointer :: row(:)
891 do i = 1, this%n_constraints
892 call this%constraint_list(i)%constraint%sensitivity%copy_from( &
893 device_to_host, sync = i .eq. this%n_constraints)
896 do i = 1, this%n_constraints
897 row(1:this%n_design) => sensitivity%x(i, :)
899 call copy(row, this%constraint_list(i)%constraint%sensitivity%x, &
903 call sensitivity%copy_from(host_to_device, sync = .true.)
905 end subroutine problem_get_constraint_sensitivities
911 pure function problem_get_num_objectives(this)
result(n)
915 n = this%n_objectives
916 end function problem_get_num_objectives
919 pure function problem_get_num_constraints(this)
result(n)
923 n = this%n_constraints
924 end function problem_get_num_constraints
930 function problem_get_log_header(this, include_constraints)
result(buff)
932 logical,
intent(in),
optional :: include_constraints
933 character(len=4096) :: buff
934 character(len=128) :: mini_buff
935 character(len=128),
allocatable :: headers(:)
937 logical :: do_constraints
939 buff =
"Total objective function"
940 if (
present(include_constraints))
then
941 do_constraints = include_constraints
943 do_constraints = .true.
946 do i = 1, this%get_n_objectives()
947 n = this%objective_list(i)%objective%get_log_size()
950 call this%objective_list(i)%objective%get_log_headers(headers)
953 write(mini_buff,
'(", ", A)') trim(headers(j))
954 buff = trim(buff) // trim(mini_buff)
960 if (do_constraints)
then
961 do i = 1, this%get_n_constraints()
962 n = this%constraint_list(i)%constraint%get_log_size()
965 call this%constraint_list(i)%constraint%get_log_headers(headers)
968 write(mini_buff,
'(", ", A)') trim(headers(j))
969 buff = trim(buff) // trim(mini_buff)
976 end function problem_get_log_header
982 function problem_get_log_size(this, include_constraints)
result(n)
984 logical,
intent(in),
optional :: include_constraints
986 logical :: do_constraints
989 if (
present(include_constraints))
then
990 do_constraints = include_constraints
992 do_constraints = .true.
994 do i = 1, this%get_n_objectives()
995 n = n + this%objective_list(i)%objective%get_log_size()
998 if (do_constraints)
then
999 do i = 1, this%get_n_constraints()
1000 n = n + this%constraint_list(i)%constraint%get_log_size()
1004 end function problem_get_log_size
Factory function Allocates and initializes an constraint function object.
Factory function Allocates and initializes an objective function object.
Implements the augmented_lagrangian_objective_t type.
Implements the constraint_t type.
Implements the objective_t type.
Module for handling the optimization problem.
subroutine problem_init(this, parameters, design, simulation)
The constructor for the base problem.
Adjoint simulation driver.
subroutine, public simulation_adjoint_init(c, dt_controller)
Initialise a simulation_adjoint of a case.
subroutine, public simulation_adjoint_step(c, dt_controller, cfl, tstep_loop_start_time, final_time)
Compute a single time-step of an adjoint case.
subroutine, public simulation_adjoint_finalize(c)
Finalize a simulation of a case.
Implements the steady_problem_t type.
An objective function implementing our augmented lagrangian sensitivity contribution.
The abstract constraint type.
Wrapper for constraints for use in lists.
The abstract objective type.
Wrapper for objectives for use in lists.
The abstract problem type.