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 vector_math, only: vector_add2, vector_cfill
56 use time_step_controller, only: time_step_controller_t
59 use simulation, only: simulation_init, simulation_step, simulation_finalize
60 use mpi_f08, only: mpi_wtime
61 use profiler, only: profiler_start_region, profiler_end_region
62 implicit none
63 private
64
66 type, public :: problem_t
67 private
68
70 integer :: n_design = 0
72 integer :: n_objectives = 0
74 integer :: n_constraints = 0
75
77 class(objective_wrapper_t), allocatable, dimension(:) :: objective_list
79 class(constraint_wrapper_t), allocatable, dimension(:) :: constraint_list
80 contains
81
82 ! ----------------------------------------------------------------------- !
83 ! Interfaces
84
86 procedure, pass(this), public :: init => problem_init
88 procedure, pass(this), public :: free => problem_free
89
93 procedure, pass(this), public :: compute => problem_compute
94
98 procedure, pass(this), public :: compute_sensitivity => &
99 problem_compute_sensitivity
102 procedure, pass(this), public :: run_forward_unsteady => &
103 problem_run_forward_unsteady
106 procedure, pass(this), public :: run_backward_unsteady => &
107 problem_run_backward_unsteady
108 ! ----------------------------------------------------------------------- !
109 ! Base class methods
110
112 procedure, pass(this), public :: read_objectives => problem_read_objectives
114 procedure, pass(this), public :: read_constraints => &
115 problem_read_constraints
116
117 ! ----------------------------------------------------------------------- !
118 ! Actual methods
119
121 procedure, pass(this), public :: write => problem_write
123 procedure, pass(this), public :: get_log_values => problem_get_log_values
124
126 procedure, pass(this), public :: add_objective => problem_add_objective
128 procedure, pass(this), public :: add_constraint => problem_add_constraint
129
130 ! ----------------------------------------------------------------------- !
131 ! Internal Updater methods
132
134 procedure, pass(this) :: update_objectives => &
135 problem_update_objectives
137 procedure, pass(this) :: update_constraints => &
138 problem_update_constraints
140 procedure, pass(this) :: update_objective_sensitivities => &
141 problem_update_objective_sensitivities
143 procedure, pass(this) :: update_constraint_sensitivities => &
144 problem_update_constraint_sensitivities
145
147 procedure, pass(this) :: reset_objectives => &
148 problem_reset_objectives
150 procedure, pass(this) :: reset_constraints => &
151 problem_reset_constraints
153 procedure, pass(this) :: reset_objective_sensitivities => &
154 problem_reset_objective_sensitivities
156 procedure, pass(this) :: reset_constraint_sensitivities => &
157 problem_reset_constraint_sensitivities
158
160 procedure, pass(this) :: accumulate_objectives => &
161 problem_accumulate_objectives
163 procedure, pass(this) :: accumulate_constraints => &
164 problem_accumulate_constraints
166 procedure, pass(this) :: accumulate_objective_sensitivities => &
167 problem_accumulate_objective_sensitivities
169 procedure, pass(this) :: accumulate_constraint_sensitivities => &
170 problem_accumulate_constraint_sensitivities
171
172 ! ----------------------------------------------------------------------- !
173 ! Public Getters
174
176 procedure, pass(this), public :: get_objective_value => &
177 problem_get_objective_value
179 procedure, pass(this), public :: get_all_objective_values => &
180 problem_get_all_objective_values
182 procedure, pass(this), public :: get_constraint_values => &
183 problem_get_constraint_values
185 procedure, pass(this), public :: get_objective_sensitivities => &
186 problem_get_objective_sensitivities
188 procedure, pass(this), public :: get_constraint_sensitivities => &
189 problem_get_constraint_sensitivities
190
192 procedure, pass(this) :: get_n_objectives => problem_get_num_objectives
194 procedure, pass(this) :: get_n_constraints => problem_get_num_constraints
195
197 procedure, pass(this) :: get_log_header => problem_get_log_header
199 procedure, pass(this) :: get_log_size => problem_get_log_size
200
201 end type problem_t
202
203contains
204
205 ! ========================================================================== !
206 ! Base class methods
207
209 subroutine problem_init(this, parameters, design, simulation)
210 class(problem_t), intent(inout) :: this
211 type(json_file), intent(inout) :: parameters
212 class(design_t), intent(in) :: design
213 type(simulation_t), optional, intent(inout) :: simulation
214
215 call this%free()
216
217 this%n_design = design%size()
218
219 ! Read the objectives and constraints
220 call this%read_objectives(parameters, design, simulation)
221 call this%read_constraints(parameters, design, simulation)
222
223 end subroutine problem_init
224
226 subroutine problem_free(this)
227 class(problem_t), intent(inout) :: this
228 integer :: i
229
230 this%n_design = 0
231 this%n_objectives = 0
232 this%n_constraints = 0
233
234 ! Free the objective list
235 if (allocated(this%objective_list)) then
236 do i = 1, size(this%objective_list)
237 call this%objective_list(i)%free()
238 end do
239 deallocate(this%objective_list)
240 end if
241
242 ! Free the constraint list
243 if (allocated(this%constraint_list)) then
244 do i = 1, size(this%constraint_list)
245 call this%constraint_list(i)%free()
246 end do
247 deallocate(this%constraint_list)
248 end if
249 end subroutine problem_free
250
252 subroutine problem_write(this, idx)
253 class(problem_t), intent(inout) :: this
254 integer, intent(in) :: idx
255
256 end subroutine problem_write
257
262 subroutine problem_get_log_values(this, values, include_constraints)
263 class(problem_t), intent(in) :: this
264 real(kind=rp), intent(out) :: values(:)
265 logical, intent(in), optional :: include_constraints
266 integer :: i, n, offset
267 real(kind=rp) :: objective_value
268 real(kind=rp), allocatable :: tmp(:)
269 logical :: do_constraints
270
271 if (present(include_constraints)) then
272 do_constraints = include_constraints
273 else
274 do_constraints = .true.
275 end if
276
277 call this%get_objective_value(objective_value)
278 values = 0.0_rp
279 values(1) = objective_value
280
281 offset = 2
282 do i = 1, this%n_objectives
283 n = this%objective_list(i)%objective%get_log_size()
284 if (n .gt. 0) then
285 allocate(tmp(n))
286 call this%objective_list(i)%objective%get_log_values(tmp)
287 values(offset:offset + n - 1) = tmp
288 offset = offset + n
289 deallocate(tmp)
290 end if
291 end do
292
293 if (do_constraints) then
294 do i = 1, this%n_constraints
295 n = this%constraint_list(i)%constraint%get_log_size()
296 if (n .gt. 0) then
297 allocate(tmp(n))
298 call this%constraint_list(i)%constraint%get_log_values(tmp)
299 values(offset:offset + n - 1) = tmp
300 offset = offset + n
301 deallocate(tmp)
302 end if
303 end do
304 end if
305 end subroutine problem_get_log_values
306
307 ! ========================================================================== !
308 ! Handling constraints and objectives
309
311 subroutine problem_read_objectives(this, parameters, design, simulation)
312 class(problem_t), intent(inout) :: this
313 type(json_file), intent(inout) :: parameters
314 class(design_t), intent(in) :: design
315 type(simulation_t), optional, intent(inout) :: simulation
316 class(objective_t), allocatable :: objective
317
318 ! A single objective term as its own json_file.
319 character(len=:), allocatable :: path, type
320 type(json_file) :: objective_json
321 integer :: n_objectives, i
322 logical :: dealias
323
324 call neko_log%section("Reading objectives")
325
326 ! Get the number of objectives.
327 path = "optimization.objectives"
328 if (parameters%valid_path(path)) then
329 call parameters%info(path, n_children = n_objectives)
330
331 ! Grab a single parameters entry and create a constraint from it.
332 do i = 1, n_objectives
333 call json_extract_item(parameters, path, i, objective_json)
334 call json_get(objective_json, "type", type)
335 call neko_log%message(type)
336
337 call objective_factory(objective, objective_json, design, simulation)
338 call this%add_objective(objective)
339 end do
340 end if
341
342 if (present(simulation)) then
343 if (allocated(objective)) deallocate(objective)
345 select type(alo => objective)
347 call json_get_or_default(parameters, &
348 "adjoint_fluid.dealias_sensitivity", dealias, .true.)
349 call alo%init_from_attributes(design, simulation, weight = 1.0_rp, &
350 name = "Augmented Lagrangian", mask_name = "", &
351 dealias = dealias)
352 end select
353 call this%add_objective(objective)
354 end if
355
356 call neko_log%end_section()
357
358 end subroutine problem_read_objectives
359
361 subroutine problem_read_constraints(this, parameters, design, simulation)
362 class(problem_t), intent(inout) :: this
363 type(json_file), intent(inout) :: parameters
364 class(design_t), intent(in) :: design
365 class(constraint_t), allocatable :: constraint
366 type(simulation_t), optional, intent(inout) :: simulation
367
368 ! A single constraint term as its own json_file.
369 character(len=:), allocatable :: path, type
370 type(json_file) :: constraint_json
371 integer :: n_constraints, i
372
373 call neko_log%section("Reading constraints")
374
375 ! Get the number of constraints.
376 path = "optimization.constraints"
377
378 if (parameters%valid_path(path)) then
379 call parameters%info(path, n_children = n_constraints)
380
381 ! Grab a single parameters entry and create a constraint from it.
382 do i = 1, n_constraints
383 call json_extract_item(parameters, path, i, constraint_json)
384 call json_get(constraint_json, "type", type)
385 call neko_log%message(type)
386
387 call constraint_factory(constraint, constraint_json, design, simulation)
388 call this%add_constraint(constraint)
389 end do
390 end if
391
392 call neko_log%end_section()
393
394 end subroutine problem_read_constraints
395
397 subroutine problem_add_objective(this, objective)
398 class(problem_t), intent(inout) :: this
399 class(objective_t), allocatable, intent(inout) :: objective
400 class(objective_wrapper_t), allocatable, dimension(:) :: temp_list
401 integer :: i, n
402
403 n = 0
404 if (allocated(this%objective_list)) then
405 n = size(this%objective_list)
406 call move_alloc(this%objective_list, temp_list)
407 allocate(this%objective_list(n + 1))
408 if (allocated(temp_list)) then
409 do i = 1, n
410 call move_alloc(temp_list(i)%objective, &
411 this%objective_list(i)%objective)
412 end do
413 end if
414 else
415 allocate(this%objective_list(1))
416 end if
417
418 call move_alloc(objective, this%objective_list(n + 1)%objective)
419 this%n_objectives = n + 1
420 end subroutine problem_add_objective
421
423 subroutine problem_add_constraint(this, constraint)
424 class(problem_t), intent(inout) :: this
425 class(constraint_t), allocatable, intent(inout) :: constraint
426 class(constraint_wrapper_t), allocatable, dimension(:) :: temp_list
427 integer :: i, n
428
429 n = 0
430 if (allocated(this%constraint_list)) then
431 n = size(this%constraint_list)
432 call move_alloc(this%constraint_list, temp_list)
433 allocate(this%constraint_list(n + 1))
434 if (allocated(temp_list)) then
435 do i = 1, n
436 call move_alloc(temp_list(i)%constraint, &
437 this%constraint_list(i)%constraint)
438 end do
439 end if
440 else
441 allocate(this%constraint_list(1))
442 end if
443
444 call move_alloc(constraint, this%constraint_list(n + 1)%constraint)
445 this%n_constraints = n + 1
446 end subroutine problem_add_constraint
447
448 ! ========================================================================== !
449 ! Problem part computation
450
452 subroutine problem_compute(this, design, simulation)
453 class(problem_t), intent(inout) :: this
454 class(design_t), intent(inout) :: design
455 class(simulation_t), optional, intent(inout) :: simulation
456
457 if (present(simulation)) then
458 call simulation%reset()
459 if (simulation%unsteady) then
460 ! Objective value accumulated
461 call this%run_forward_unsteady(simulation, design)
462 else
463 call simulation%run_forward()
464 ! Compute objective value on steady field
465 call this%update_objectives(design)
466 end if
467 else
468 call this%update_objectives(design)
469 end if
470
471 call this%update_constraints(design)
472
473 end subroutine problem_compute
474
476 subroutine problem_compute_sensitivity(this, design, simulation)
477 class(problem_t), intent(inout) :: this
478 class(design_t), intent(inout) :: design
479 class(simulation_t), optional, intent(inout) :: simulation
480
481 type(vector_t) :: objective_sensitivity
482
483 if (present(simulation)) then
484 if (simulation%unsteady) then
485 ! Objective sensitivity accumulated
486 call this%run_backward_unsteady(simulation, design)
487 else
488 call simulation%run_backward()
489 ! Compute sensitivity on steady field
490 call this%update_objective_sensitivities(design)
491 end if
492 else
493 call this%update_objective_sensitivities(design)
494 end if
495
496 call this%update_constraint_sensitivities(design)
497
498 call objective_sensitivity%init(this%n_design)
499 call this%get_objective_sensitivities(objective_sensitivity)
500
501 call design%map_backward(objective_sensitivity)
502
503 call objective_sensitivity%free()
504 end subroutine problem_compute_sensitivity
505
507 subroutine problem_run_forward_unsteady(this, simulation, design)
508 class(problem_t), intent(inout) :: this
509 class(simulation_t), intent(inout) :: simulation
510 class(design_t), intent(inout) :: design
511 type(time_step_controller_t) :: dt_controller
512 real(kind=dp) :: loop_start
513
514 call dt_controller%init(simulation%neko_case%params)
515
516 call simulation%reset()
517 call simulation_init(simulation%neko_case, dt_controller)
518
519 ! Reset the objective value to zero
520 call this%reset_objectives()
521
522 call profiler_start_region("Forward simulation")
523 loop_start = mpi_wtime()
524 simulation%n_timesteps = 0
525 do while (simulation%neko_case%time%t .lt. simulation%neko_case%time%end_time)
526 simulation%n_timesteps = simulation%n_timesteps + 1
527 ! step forward
528 call simulation_step(simulation%neko_case, dt_controller, loop_start)
529 ! accumulate objective value
530 call this%accumulate_objectives(design, simulation%adjoint_case%time%dt)
531 ! save a checkpoint
532 call simulation%checkpoint%save(simulation%neko_case)
533 end do
534 call profiler_end_region("Forward simulation")
535
536 call simulation_finalize(simulation%neko_case)
537
538 end subroutine problem_run_forward_unsteady
539
541 subroutine problem_run_backward_unsteady(this, simulation, design)
542 class(problem_t), intent(inout) :: this
543 class(simulation_t), intent(inout) :: simulation
544 class(design_t), intent(inout) :: design
545 type(time_step_controller_t) :: dt_controller
546 real(kind=dp) :: loop_start
547 real(kind=rp) :: cfl
548 real(kind=rp) :: total_time
549 integer :: i
550
551 call dt_controller%init(simulation%neko_case%params)
552
553 call simulation_adjoint_init(simulation%adjoint_case, dt_controller)
554
555 ! Reset the sensitivity value to zero
556 call this%reset_objective_sensitivities()
557
558 cfl = simulation%adjoint_case%fluid_adj%compute_cfl(simulation%adjoint_case%time%dt)
559 loop_start = mpi_wtime()
560
561 total_time = simulation%n_timesteps * simulation%adjoint_case%time%dt
562
563 call profiler_start_region("Adjoint simulation")
564
565 ! TODO. IC's need to be handled rather carefully, we should take a
566 ! checkpoint at i = 0. However, 99% of the time we use an initial condition
567 ! for the fluid of u=0, so this doesn't matter. Then we use u_adj = 0 on
568 ! the other end.
569 !
570 ! we have:
571 ! - n time steps to compute
572 ! - n+1 fields to consider
573 ! - n-1 non-zero contributions to u * u_adj
574 !
575 ! non-zero
576 ! |--------|
577 ! primal o--x--x--x--x--x
578 ! x--x--x--x--x--o adjoint
579 ! ^ ^
580 ! u=0 u_adj=0
581
582 do i = simulation%n_timesteps, 1, -1
583 ! restore primal field
584 call simulation%checkpoint%restore(simulation%neko_case, i)
585 ! accumulate objective sensitivity
586 call this%accumulate_objective_sensitivities(design, &
587 simulation%adjoint_case%time%dt)
588 ! step the adjoint backwards
589 call simulation_adjoint_step(simulation%adjoint_case, dt_controller, &
590 cfl, loop_start, total_time)
591 end do
592
593 call profiler_end_region("Forward simulation")
594
595 call simulation_adjoint_finalize(simulation%adjoint_case)
596
597 end subroutine problem_run_backward_unsteady
598
599 ! ========================================================================== !
600 ! Update the objectives and constraints
601
608 subroutine problem_update_objectives(this, design)
609 class(problem_t), intent(inout) :: this
610 class(design_t), intent(in) :: design
611 integer :: i
612
613 do i = 1, this%n_objectives
614 call this%objective_list(i)%objective%update_value(design)
615 end do
616 end subroutine problem_update_objectives
617
624 subroutine problem_update_constraints(this, design)
625 class(problem_t), intent(inout) :: this
626 class(design_t), intent(in) :: design
627 integer :: i
628
629 do i = 1, this%n_constraints
630 call this%constraint_list(i)%constraint%update_value(design)
631 end do
632 end subroutine problem_update_constraints
633
640 subroutine problem_update_objective_sensitivities(this, design)
641 class(problem_t), intent(inout) :: this
642 class(design_t), intent(in) :: design
643 integer :: i
644
645 do i = 1, this%n_objectives
646 call this%objective_list(i)%objective%update_sensitivity(design)
647 end do
648 end subroutine problem_update_objective_sensitivities
649
656 subroutine problem_update_constraint_sensitivities(this, design)
657 class(problem_t), intent(inout) :: this
658 class(design_t), intent(in) :: design
659 integer :: i
660
661 do i = 1, this%n_constraints
662 call this%constraint_list(i)%constraint%update_sensitivity(design)
663 end do
664 end subroutine problem_update_constraint_sensitivities
665
666 ! ========================================================================== !
667 ! Reset the objectives and constraints
668
673 subroutine problem_reset_objectives(this)
674 class(problem_t), intent(inout) :: this
675 integer :: i
676
677 do i = 1, this%n_objectives
678 call this%objective_list(i)%objective%reset_value()
679 end do
680 end subroutine problem_reset_objectives
681
686 subroutine problem_reset_constraints(this)
687 class(problem_t), intent(inout) :: this
688 integer :: i
689
690 do i = 1, this%n_constraints
691 call this%constraint_list(i)%constraint%reset_value()
692 end do
693 end subroutine problem_reset_constraints
694
699 subroutine problem_reset_objective_sensitivities(this)
700 class(problem_t), intent(inout) :: this
701 integer :: i
702
703 do i = 1, this%n_objectives
704 call this%objective_list(i)%objective%reset_sensitivity()
705 end do
706 end subroutine problem_reset_objective_sensitivities
707
712 subroutine problem_reset_constraint_sensitivities(this)
713 class(problem_t), intent(inout) :: this
714 integer :: i
715
716 do i = 1, this%n_constraints
717 call this%constraint_list(i)%constraint%reset_sensitivity()
718 end do
719 end subroutine problem_reset_constraint_sensitivities
720
721 ! ========================================================================== !
722 ! Accumulate the objectives and constraints
723
730 subroutine problem_accumulate_objectives(this, design, dt)
731 class(problem_t), intent(inout) :: this
732 class(design_t), intent(in) :: design
733 real(kind=rp), intent(in) :: dt
734 integer :: i
735
736 do i = 1, this%n_objectives
737 call this%objective_list(i)%objective%accumulate_value(design, dt)
738 end do
739 end subroutine problem_accumulate_objectives
740
747 subroutine problem_accumulate_constraints(this, design, dt)
748 class(problem_t), intent(inout) :: this
749 class(design_t), intent(in) :: design
750 real(kind=rp), intent(in) :: dt
751 integer :: i
752
753 do i = 1, this%n_constraints
754 call this%constraint_list(i)%constraint%accumulate_value(design, dt)
755 end do
756 end subroutine problem_accumulate_constraints
757
764 subroutine problem_accumulate_objective_sensitivities(this, design, dt)
765 class(problem_t), intent(inout) :: this
766 class(design_t), intent(in) :: design
767 real(kind=rp), intent(in) :: dt
768 integer :: i
769
770 do i = 1, this%n_objectives
771 call this%objective_list(i)%objective%accumulate_sensitivity(design, dt)
772 end do
773 end subroutine problem_accumulate_objective_sensitivities
774
781 subroutine problem_accumulate_constraint_sensitivities(this, design, dt)
782 class(problem_t), intent(inout) :: this
783 class(design_t), intent(in) :: design
784 real(kind=rp), intent(in) :: dt
785 integer :: i
786
787 do i = 1, this%n_constraints
788 call this%constraint_list(i)%constraint%accumulate_sensitivity(design, dt)
789 end do
790 end subroutine problem_accumulate_constraint_sensitivities
791
792 ! ========================================================================== !
793 ! Problem part getters
794
801 subroutine problem_get_objective_value(this, objective_value)
802 class(problem_t), intent(in) :: this
803 real(kind=rp), intent(out) :: objective_value
804 integer :: i
805
806 objective_value = 0.0_rp
807 do i = 1, this%n_objectives
808 objective_value = objective_value + &
809 this%objective_list(i)%objective%get_weight() * &
810 this%objective_list(i)%objective%get_value()
811 end do
812
813 end subroutine problem_get_objective_value
814
821 subroutine problem_get_all_objective_values(this, all_objective_values)
822 class(problem_t), intent(in) :: this
823 type(vector_t), intent(inout) :: all_objective_values
824 integer :: i
825
826 do i = 1, this%n_objectives
827 all_objective_values%x(i) = this%objective_list(i)%objective%value
828 end do
829
830 call all_objective_values%copy_from(host_to_device, sync = .true.)
831
832 end subroutine problem_get_all_objective_values
833
840 subroutine problem_get_constraint_values(this, constraint_value)
841 class(problem_t), intent(in) :: this
842 type(vector_t), intent(inout) :: constraint_value
843 integer :: i
844
845 do i = 1, this%n_constraints
846 constraint_value%x(i) = this%constraint_list(i)%constraint%value
847 end do
848
849 call constraint_value%copy_from(host_to_device, sync = .true.)
850
851 end subroutine problem_get_constraint_values
852
859 subroutine problem_get_objective_sensitivities(this, sensitivity)
860 class(problem_t), intent(in) :: this
861 type(vector_t), intent(inout) :: sensitivity
862 integer :: i
863
864 call vector_cfill(sensitivity, 0.0_rp)
865 do i = 1, this%n_objectives
866 call vector_add2(sensitivity, &
867 this%objective_list(i)%objective%sensitivity)
868 end do
869
870 end subroutine problem_get_objective_sensitivities
871
878 subroutine problem_get_constraint_sensitivities(this, sensitivity)
879 class(problem_t), intent(inout) :: this
880 type(matrix_t), target, intent(inout) :: sensitivity
881 real(kind=rp), pointer :: row(:)
882 integer :: i
883
884 ! Copy all constraint sensitivities to host, sync on last one
885 do i = 1, this%n_constraints
886 call this%constraint_list(i)%constraint%sensitivity%copy_from( &
887 device_to_host, sync = i .eq. this%n_constraints)
888 end do
889
890 do i = 1, this%n_constraints
891 row(1:this%n_design) => sensitivity%x(i, :)
892
893 call copy(row, this%constraint_list(i)%constraint%sensitivity%x, &
894 this%n_design)
895 end do
896
897 call sensitivity%copy_from(host_to_device, sync = .true.)
898
899 end subroutine problem_get_constraint_sensitivities
900
901 ! ========================================================================== !
902 ! Simple getters
903
905 pure function problem_get_num_objectives(this) result(n)
906 class(problem_t), intent(in) :: this
907 integer :: n
908
909 n = this%n_objectives
910 end function problem_get_num_objectives
911
913 pure function problem_get_num_constraints(this) result(n)
914 class(problem_t), intent(in) :: this
915 integer :: n
916
917 n = this%n_constraints
918 end function problem_get_num_constraints
919
924 function problem_get_log_header(this, include_constraints) result(buff)
925 class(problem_t), intent(in) :: this
926 logical, intent(in), optional :: include_constraints
927 character(len=4096) :: buff
928 character(len=128) :: mini_buff
929 character(len=128), allocatable :: headers(:)
930 integer :: i, j, n
931 logical :: do_constraints
932
933 buff = "Total objective function"
934 if (present(include_constraints)) then
935 do_constraints = include_constraints
936 else
937 do_constraints = .true.
938 end if
939
940 do i = 1, this%get_n_objectives()
941 n = this%objective_list(i)%objective%get_log_size()
942 if (n .gt. 0) then
943 allocate(headers(n))
944 call this%objective_list(i)%objective%get_log_headers(headers)
945 do j = 1, n
946 mini_buff = ""
947 write(mini_buff, '(", ", A)') trim(headers(j))
948 buff = trim(buff) // trim(mini_buff)
949 end do
950 deallocate(headers)
951 end if
952 end do
953
954 if (do_constraints) then
955 do i = 1, this%get_n_constraints()
956 n = this%constraint_list(i)%constraint%get_log_size()
957 if (n .gt. 0) then
958 allocate(headers(n))
959 call this%constraint_list(i)%constraint%get_log_headers(headers)
960 do j = 1, n
961 mini_buff = ""
962 write(mini_buff, '(", ", A)') trim(headers(j))
963 buff = trim(buff) // trim(mini_buff)
964 end do
965 deallocate(headers)
966 end if
967 end do
968 end if
969
970 end function problem_get_log_header
971
976 function problem_get_log_size(this, include_constraints) result(n)
977 class(problem_t), intent(in) :: this
978 logical, intent(in), optional :: include_constraints
979 integer :: n, i
980 logical :: do_constraints
981
982 n = 1
983 if (present(include_constraints)) then
984 do_constraints = include_constraints
985 else
986 do_constraints = .true.
987 end if
988 do i = 1, this%get_n_objectives()
989 n = n + this%objective_list(i)%objective%get_log_size()
990 end do
991
992 if (do_constraints) then
993 do i = 1, this%get_n_constraints()
994 n = n + this%constraint_list(i)%constraint%get_log_size()
995 end do
996 end if
997
998 end function problem_get_log_size
999end 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:210
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:54
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:66