Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
problem.f90
Go to the documentation of this file.
1
34
41module problem
42 use num_types, only: rp, dp
43 use design, only: design_t
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
52 use simulation_m, only: simulation_t
53 use logger, only: neko_log
54 use math, only: copy
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
63 implicit none
64 private
65
67 type, public :: problem_t
68 private
69
71 integer :: n_design = 0
73 integer :: n_objectives = 0
75 integer :: n_constraints = 0
76
78 class(objective_wrapper_t), allocatable, dimension(:) :: objective_list
80 class(constraint_wrapper_t), allocatable, dimension(:) :: constraint_list
81 contains
82
83 ! ----------------------------------------------------------------------- !
84 ! Interfaces
85
87 procedure, pass(this), public :: init => problem_init
89 procedure, pass(this), public :: free => problem_free
90
94 procedure, pass(this), public :: compute => problem_compute
95
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
109 ! ----------------------------------------------------------------------- !
110 ! Base class methods
111
113 procedure, pass(this), public :: read_objectives => problem_read_objectives
115 procedure, pass(this), public :: read_constraints => &
116 problem_read_constraints
117
118 ! ----------------------------------------------------------------------- !
119 ! Actual methods
120
122 procedure, pass(this), public :: write => problem_write
124 procedure, pass(this), public :: get_log_values => problem_get_log_values
125
127 procedure, pass(this), public :: add_objective => problem_add_objective
129 procedure, pass(this), public :: add_constraint => problem_add_constraint
130
131 ! ----------------------------------------------------------------------- !
132 ! Internal Updater methods
133
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
146
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
159
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
172
173 ! ----------------------------------------------------------------------- !
174 ! Public Getters
175
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
191
193 procedure, pass(this) :: get_n_objectives => problem_get_num_objectives
195 procedure, pass(this) :: get_n_constraints => problem_get_num_constraints
196
198 procedure, pass(this) :: get_log_header => problem_get_log_header
200 procedure, pass(this) :: get_log_size => problem_get_log_size
201
202 end type problem_t
203
204contains
205
206 ! ========================================================================== !
207 ! Base class methods
208
210 subroutine problem_init(this, parameters, design, simulation)
211 class(problem_t), intent(inout) :: this
212 type(json_file), intent(inout) :: parameters
213 class(design_t), intent(in) :: design
214 type(simulation_t), optional, intent(inout) :: simulation
215
216 call this%free()
217
218 this%n_design = design%size()
219
220 ! Read the objectives and constraints
221 call this%read_objectives(parameters, design, simulation)
222 call this%read_constraints(parameters, design, simulation)
223
224 end subroutine problem_init
225
227 subroutine problem_free(this)
228 class(problem_t), intent(inout) :: this
229 integer :: i
230
231 this%n_design = 0
232 this%n_objectives = 0
233 this%n_constraints = 0
234
235 ! Free the objective list
236 if (allocated(this%objective_list)) then
237 do i = 1, size(this%objective_list)
238 call this%objective_list(i)%free()
239 end do
240 deallocate(this%objective_list)
241 end if
242
243 ! Free the constraint list
244 if (allocated(this%constraint_list)) then
245 do i = 1, size(this%constraint_list)
246 call this%constraint_list(i)%free()
247 end do
248 deallocate(this%constraint_list)
249 end if
250 end subroutine problem_free
251
253 subroutine problem_write(this, idx)
254 class(problem_t), intent(inout) :: this
255 integer, intent(in) :: idx
256
257 end subroutine problem_write
258
263 subroutine problem_get_log_values(this, values, include_constraints)
264 class(problem_t), intent(in) :: this
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
271
272 if (present(include_constraints)) then
273 do_constraints = include_constraints
274 else
275 do_constraints = .true.
276 end if
277
278 call this%get_objective_value(objective_value)
279 values = 0.0_rp
280 values(1) = objective_value
281
282 offset = 2
283 do i = 1, this%n_objectives
284 n = this%objective_list(i)%objective%get_log_size()
285 if (n .gt. 0) then
286 allocate(tmp(n))
287 call this%objective_list(i)%objective%get_log_values(tmp)
288 values(offset:offset + n - 1) = tmp
289 offset = offset + n
290 deallocate(tmp)
291 end if
292 end do
293
294 if (do_constraints) then
295 do i = 1, this%n_constraints
296 n = this%constraint_list(i)%constraint%get_log_size()
297 if (n .gt. 0) then
298 allocate(tmp(n))
299 call this%constraint_list(i)%constraint%get_log_values(tmp)
300 values(offset:offset + n - 1) = tmp
301 offset = offset + n
302 deallocate(tmp)
303 end if
304 end do
305 end if
306 end subroutine problem_get_log_values
307
308 ! ========================================================================== !
309 ! Handling constraints and objectives
310
312 subroutine problem_read_objectives(this, parameters, design, simulation)
313 class(problem_t), intent(inout) :: this
314 type(json_file), intent(inout) :: parameters
315 class(design_t), intent(in) :: design
316 type(simulation_t), optional, intent(inout) :: simulation
317 class(objective_t), allocatable :: objective
318
319 ! A single objective term as its own json_file.
320 character(len=:), allocatable :: path, type
321 type(json_file) :: objective_json
322 integer :: n_objectives, i
323 logical :: dealias
324
325 call neko_log%section("Reading objectives")
326
327 ! Get the number of objectives.
328 path = "optimization.objectives"
329 if (parameters%valid_path(path)) then
330 call parameters%info(path, n_children = n_objectives)
331
332 ! Grab a single parameters entry and create a constraint from it.
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)
337
338 call objective_factory(objective, objective_json, design, simulation)
339 call this%add_objective(objective)
340 end do
341 end if
342
343 if (present(simulation)) then
344 if (allocated(objective)) deallocate(objective)
346 select type(alo => objective)
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 = "", &
352 dealias = dealias)
353 end select
354 call this%add_objective(objective)
355 end if
356
357 call neko_log%end_section()
358
359 end subroutine problem_read_objectives
360
362 subroutine problem_read_constraints(this, parameters, design, simulation)
363 class(problem_t), intent(inout) :: this
364 type(json_file), intent(inout) :: parameters
365 class(design_t), intent(in) :: design
366 class(constraint_t), allocatable :: constraint
367 type(simulation_t), optional, intent(inout) :: simulation
368
369 ! A single constraint term as its own json_file.
370 character(len=:), allocatable :: path, type
371 type(json_file) :: constraint_json
372 integer :: n_constraints, i
373
374 call neko_log%section("Reading constraints")
375
376 ! Get the number of constraints.
377 path = "optimization.constraints"
378
379 if (parameters%valid_path(path)) then
380 call parameters%info(path, n_children = n_constraints)
381
382 ! Grab a single parameters entry and create a constraint from it.
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)
387
388 call constraint_factory(constraint, constraint_json, design, simulation)
389 call this%add_constraint(constraint)
390 end do
391 end if
392
393 call neko_log%end_section()
394
395 end subroutine problem_read_constraints
396
398 subroutine problem_add_objective(this, objective)
399 class(problem_t), intent(inout) :: this
400 class(objective_t), allocatable, intent(inout) :: objective
401 class(objective_wrapper_t), allocatable, dimension(:) :: temp_list
402 integer :: i, n
403
404 n = 0
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
410 do i = 1, n
411 call move_alloc(temp_list(i)%objective, &
412 this%objective_list(i)%objective)
413 end do
414 end if
415 else
416 allocate(this%objective_list(1))
417 end if
418
419 call move_alloc(objective, this%objective_list(n + 1)%objective)
420 this%n_objectives = n + 1
421 end subroutine problem_add_objective
422
424 subroutine problem_add_constraint(this, constraint)
425 class(problem_t), intent(inout) :: this
426 class(constraint_t), allocatable, intent(inout) :: constraint
427 class(constraint_wrapper_t), allocatable, dimension(:) :: temp_list
428 integer :: i, n
429
430 n = 0
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
436 do i = 1, n
437 call move_alloc(temp_list(i)%constraint, &
438 this%constraint_list(i)%constraint)
439 end do
440 end if
441 else
442 allocate(this%constraint_list(1))
443 end if
444
445 call move_alloc(constraint, this%constraint_list(n + 1)%constraint)
446 this%n_constraints = n + 1
447 end subroutine problem_add_constraint
448
449 ! ========================================================================== !
450 ! Problem part computation
451
453 subroutine problem_compute(this, design, simulation)
454 class(problem_t), intent(inout) :: this
455 class(design_t), intent(inout) :: design
456 class(simulation_t), optional, intent(inout) :: simulation
457
458 if (present(simulation)) then
459 call simulation%reset()
460 if (simulation%unsteady) then
461 ! Objective value accumulated
462 call this%run_forward_unsteady(simulation, design)
463 else
464 call simulation%run_forward()
465 ! Compute objective value on steady field
466 call this%update_objectives(design)
467 end if
468 else
469 call this%update_objectives(design)
470 end if
471
472 call this%update_constraints(design)
473
474 end subroutine problem_compute
475
477 subroutine problem_compute_sensitivity(this, design, simulation)
478 class(problem_t), intent(inout) :: this
479 class(design_t), intent(inout) :: design
480 class(simulation_t), optional, intent(inout) :: simulation
481
482 type(vector_t) :: objective_sensitivity
483
484 if (present(simulation)) then
485 if (simulation%unsteady) then
486 ! Objective sensitivity accumulated
487 call this%run_backward_unsteady(simulation, design)
488 else
489 call simulation%run_backward()
490 ! Compute sensitivity on steady field
491 call this%update_objective_sensitivities(design)
492 end if
493 else
494 call this%update_objective_sensitivities(design)
495 end if
496
497 call this%update_constraint_sensitivities(design)
498
499 call objective_sensitivity%init(this%n_design)
500 call this%get_objective_sensitivities(objective_sensitivity)
501
502 call design%map_backward(objective_sensitivity)
503
504 call objective_sensitivity%free()
505 end subroutine problem_compute_sensitivity
506
508 subroutine problem_run_forward_unsteady(this, simulation, design)
509 class(problem_t), intent(inout) :: this
510 class(simulation_t), intent(inout) :: simulation
511 class(design_t), intent(inout) :: design
512 type(time_step_controller_t) :: dt_controller
513 real(kind=dp) :: loop_start
514
515 call dt_controller%init(simulation%neko_case%params)
516
517 call simulation%reset()
518 call simulation_init(simulation%neko_case, dt_controller)
519
520 ! Reset the objective value to zero
521 call this%reset_objectives()
522
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
528 ! step forward
529 call simulation_step(simulation%neko_case, dt_controller, loop_start)
530 ! accumulate objective value
531 call this%accumulate_objectives(design, simulation%neko_case%time)
532 ! save a checkpoint
533 call simulation%checkpoint%save(simulation%neko_case)
534 end do
535 call profiler_end_region("Forward simulation")
536
537 call simulation_finalize(simulation%neko_case)
538
539 end subroutine problem_run_forward_unsteady
540
542 subroutine problem_run_backward_unsteady(this, simulation, design)
543 class(problem_t), intent(inout) :: this
544 class(simulation_t), intent(inout) :: simulation
545 class(design_t), intent(inout) :: design
546 type(time_step_controller_t) :: dt_controller
547 real(kind=dp) :: loop_start
548 real(kind=rp) :: cfl
549 real(kind=rp) :: total_time
550 integer :: i
551 type(time_state_t) :: accumulation_time
552
553 call dt_controller%init(simulation%neko_case%params)
554
555 call simulation_adjoint_init(simulation%adjoint_case, dt_controller)
556
557 ! Reset the sensitivity value to zero
558 call this%reset_objective_sensitivities()
559
560 cfl = simulation%adjoint_case%fluid_adj%compute_cfl(simulation%adjoint_case%time%dt)
561 loop_start = mpi_wtime()
562
563 ! Total time of the forward simulation
564 total_time = simulation%n_timesteps * simulation%adjoint_case%time%dt
565
566 call profiler_start_region("Adjoint simulation")
567
568 ! TODO. IC's need to be handled rather carefully, we should take a
569 ! checkpoint at i = 0. However, 99% of the time we use an initial condition
570 ! for the fluid of u=0, so this doesn't matter. Then we use u_adj = 0 on
571 ! the other end.
572 !
573 ! we have:
574 ! - n time steps to compute
575 ! - n+1 fields to consider
576 ! - n-1 non-zero contributions to u * u_adj
577 !
578 ! non-zero
579 ! |--------|
580 ! primal o--x--x--x--x--x
581 ! x--x--x--x--x--o adjoint
582 ! ^ ^
583 ! u=0 u_adj=0
584
585 do i = simulation%n_timesteps, 1, -1
586 ! restore primal field
587 call simulation%checkpoint%restore(simulation%neko_case, i)
588 ! accumulate objective sensitivity
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)
592 ! step the adjoint backwards
593 call simulation_adjoint_step(simulation%adjoint_case, dt_controller, &
594 cfl, loop_start, total_time)
595 end do
596
597 call profiler_end_region("Adjoint simulation")
598
599 call simulation_adjoint_finalize(simulation%adjoint_case)
600
601 end subroutine problem_run_backward_unsteady
602
603 ! ========================================================================== !
604 ! Update the objectives and constraints
605
612 subroutine problem_update_objectives(this, design)
613 class(problem_t), intent(inout) :: this
614 class(design_t), intent(in) :: design
615 integer :: i
616
617 do i = 1, this%n_objectives
618 call this%objective_list(i)%objective%update_value(design)
619 end do
620 end subroutine problem_update_objectives
621
628 subroutine problem_update_constraints(this, design)
629 class(problem_t), intent(inout) :: this
630 class(design_t), intent(in) :: design
631 integer :: i
632
633 do i = 1, this%n_constraints
634 call this%constraint_list(i)%constraint%update_value(design)
635 end do
636 end subroutine problem_update_constraints
637
644 subroutine problem_update_objective_sensitivities(this, design)
645 class(problem_t), intent(inout) :: this
646 class(design_t), intent(in) :: design
647 integer :: i
648
649 do i = 1, this%n_objectives
650 call this%objective_list(i)%objective%update_sensitivity(design)
651 end do
652 end subroutine problem_update_objective_sensitivities
653
660 subroutine problem_update_constraint_sensitivities(this, design)
661 class(problem_t), intent(inout) :: this
662 class(design_t), intent(in) :: design
663 integer :: i
664
665 do i = 1, this%n_constraints
666 call this%constraint_list(i)%constraint%update_sensitivity(design)
667 end do
668 end subroutine problem_update_constraint_sensitivities
669
670 ! ========================================================================== !
671 ! Reset the objectives and constraints
672
677 subroutine problem_reset_objectives(this)
678 class(problem_t), intent(inout) :: this
679 integer :: i
680
681 do i = 1, this%n_objectives
682 call this%objective_list(i)%objective%reset_value()
683 end do
684 end subroutine problem_reset_objectives
685
690 subroutine problem_reset_constraints(this)
691 class(problem_t), intent(inout) :: this
692 integer :: i
693
694 do i = 1, this%n_constraints
695 call this%constraint_list(i)%constraint%reset_value()
696 end do
697 end subroutine problem_reset_constraints
698
703 subroutine problem_reset_objective_sensitivities(this)
704 class(problem_t), intent(inout) :: this
705 integer :: i
706
707 do i = 1, this%n_objectives
708 call this%objective_list(i)%objective%reset_sensitivity()
709 end do
710 end subroutine problem_reset_objective_sensitivities
711
716 subroutine problem_reset_constraint_sensitivities(this)
717 class(problem_t), intent(inout) :: this
718 integer :: i
719
720 do i = 1, this%n_constraints
721 call this%constraint_list(i)%constraint%reset_sensitivity()
722 end do
723 end subroutine problem_reset_constraint_sensitivities
724
725 ! ========================================================================== !
726 ! Accumulate the objectives and constraints
727
734 subroutine problem_accumulate_objectives(this, design, time)
735 class(problem_t), intent(inout) :: this
736 class(design_t), intent(in) :: design
737 type(time_state_t), intent(in) :: time
738 integer :: i
739
740 do i = 1, this%n_objectives
741 call this%objective_list(i)%objective%accumulate_value(design, time)
742 end do
743 end subroutine problem_accumulate_objectives
744
751 subroutine problem_accumulate_constraints(this, design, time)
752 class(problem_t), intent(inout) :: this
753 class(design_t), intent(in) :: design
754 type(time_state_t), intent(in) :: time
755 integer :: i
756
757 do i = 1, this%n_constraints
758 call this%constraint_list(i)%constraint%accumulate_value(design, time)
759 end do
760 end subroutine problem_accumulate_constraints
761
768 subroutine problem_accumulate_objective_sensitivities(this, design, time)
769 class(problem_t), intent(inout) :: this
770 class(design_t), intent(in) :: design
771 type(time_state_t), intent(in) :: time
772 integer :: i
773
774 do i = 1, this%n_objectives
775 call this%objective_list(i)%objective%accumulate_sensitivity(design, &
776 time)
777 end do
778 end subroutine problem_accumulate_objective_sensitivities
779
786 subroutine problem_accumulate_constraint_sensitivities(this, design, time)
787 class(problem_t), intent(inout) :: this
788 class(design_t), intent(in) :: design
789 type(time_state_t), intent(in) :: time
790 integer :: i
791
792 do i = 1, this%n_constraints
793 call this%constraint_list(i)%constraint%accumulate_sensitivity(design, &
794 time)
795 end do
796 end subroutine problem_accumulate_constraint_sensitivities
797
798 ! ========================================================================== !
799 ! Problem part getters
800
807 subroutine problem_get_objective_value(this, objective_value)
808 class(problem_t), intent(in) :: this
809 real(kind=rp), intent(out) :: objective_value
810 integer :: i
811
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()
817 end do
818
819 end subroutine problem_get_objective_value
820
827 subroutine problem_get_all_objective_values(this, all_objective_values)
828 class(problem_t), intent(in) :: this
829 type(vector_t), intent(inout) :: all_objective_values
830 integer :: i
831
832 do i = 1, this%n_objectives
833 all_objective_values%x(i) = this%objective_list(i)%objective%value
834 end do
835
836 call all_objective_values%copy_from(host_to_device, sync = .true.)
837
838 end subroutine problem_get_all_objective_values
839
846 subroutine problem_get_constraint_values(this, constraint_value)
847 class(problem_t), intent(in) :: this
848 type(vector_t), intent(inout) :: constraint_value
849 integer :: i
850
851 do i = 1, this%n_constraints
852 constraint_value%x(i) = this%constraint_list(i)%constraint%value
853 end do
854
855 call constraint_value%copy_from(host_to_device, sync = .true.)
856
857 end subroutine problem_get_constraint_values
858
865 subroutine problem_get_objective_sensitivities(this, sensitivity)
866 class(problem_t), intent(in) :: this
867 type(vector_t), intent(inout) :: sensitivity
868 integer :: i
869
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)
874 end do
875
876 end subroutine problem_get_objective_sensitivities
877
884 subroutine problem_get_constraint_sensitivities(this, sensitivity)
885 class(problem_t), intent(inout) :: this
886 type(matrix_t), target, intent(inout) :: sensitivity
887 real(kind=rp), pointer :: row(:)
888 integer :: i
889
890 ! Copy all constraint sensitivities to host, sync on last one
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)
894 end do
895
896 do i = 1, this%n_constraints
897 row(1:this%n_design) => sensitivity%x(i, :)
898
899 call copy(row, this%constraint_list(i)%constraint%sensitivity%x, &
900 this%n_design)
901 end do
902
903 call sensitivity%copy_from(host_to_device, sync = .true.)
904
905 end subroutine problem_get_constraint_sensitivities
906
907 ! ========================================================================== !
908 ! Simple getters
909
911 pure function problem_get_num_objectives(this) result(n)
912 class(problem_t), intent(in) :: this
913 integer :: n
914
915 n = this%n_objectives
916 end function problem_get_num_objectives
917
919 pure function problem_get_num_constraints(this) result(n)
920 class(problem_t), intent(in) :: this
921 integer :: n
922
923 n = this%n_constraints
924 end function problem_get_num_constraints
925
930 function problem_get_log_header(this, include_constraints) result(buff)
931 class(problem_t), intent(in) :: this
932 logical, intent(in), optional :: include_constraints
933 character(len=4096) :: buff
934 character(len=128) :: mini_buff
935 character(len=128), allocatable :: headers(:)
936 integer :: i, j, n
937 logical :: do_constraints
938
939 buff = "Total objective function"
940 if (present(include_constraints)) then
941 do_constraints = include_constraints
942 else
943 do_constraints = .true.
944 end if
945
946 do i = 1, this%get_n_objectives()
947 n = this%objective_list(i)%objective%get_log_size()
948 if (n .gt. 0) then
949 allocate(headers(n))
950 call this%objective_list(i)%objective%get_log_headers(headers)
951 do j = 1, n
952 mini_buff = ""
953 write(mini_buff, '(", ", A)') trim(headers(j))
954 buff = trim(buff) // trim(mini_buff)
955 end do
956 deallocate(headers)
957 end if
958 end do
959
960 if (do_constraints) then
961 do i = 1, this%get_n_constraints()
962 n = this%constraint_list(i)%constraint%get_log_size()
963 if (n .gt. 0) then
964 allocate(headers(n))
965 call this%constraint_list(i)%constraint%get_log_headers(headers)
966 do j = 1, n
967 mini_buff = ""
968 write(mini_buff, '(", ", A)') trim(headers(j))
969 buff = trim(buff) // trim(mini_buff)
970 end do
971 deallocate(headers)
972 end if
973 end do
974 end if
975
976 end function problem_get_log_header
977
982 function problem_get_log_size(this, include_constraints) result(n)
983 class(problem_t), intent(in) :: this
984 logical, intent(in), optional :: include_constraints
985 integer :: n, i
986 logical :: do_constraints
987
988 n = 1
989 if (present(include_constraints)) then
990 do_constraints = include_constraints
991 else
992 do_constraints = .true.
993 end if
994 do i = 1, this%get_n_objectives()
995 n = n + this%objective_list(i)%objective%get_log_size()
996 end do
997
998 if (do_constraints) then
999 do i = 1, this%get_n_constraints()
1000 n = n + this%constraint_list(i)%constraint%get_log_size()
1001 end do
1002 end if
1003
1004 end function problem_get_log_size
1005end module problem
Factory function Allocates and initializes an constraint function object.
Factory function Allocates and initializes an objective function object.
Definition objective.f90:89
Implements the augmented_lagrangian_objective_t type.
Implements the constraint_t type.
Implements the design_t.
Definition design.f90:36
Implements the objective_t type.
Definition objective.f90:36
Module for handling the optimization problem.
Definition problem.f90:41
subroutine problem_init(this, parameters, design, simulation)
The constructor for the base problem.
Definition problem.f90:211
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.
An abstract design type.
Definition design.f90:53
The abstract objective type.
Definition objective.f90:52
Wrapper for objectives for use in lists.
Definition objective.f90:74
The abstract problem type.
Definition problem.f90:67