Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
design_brinkman.f90
Go to the documentation of this file.
1
34
35! Implements the `brinkman_design_t` type.
36module brinkman_design
37 use num_types, only: rp, sp, dp
38 use field, only: field_t
39 use json_module, only: json_file
41 use coefs, only: coef_t
43 use scratch_registry, only: neko_scratch_registry
44 use point_zone_registry, only: neko_point_zone_registry
45 use point_zone, only: point_zone_t
47 use neko_config, only: neko_bcknd_device
48 use device, only: device_memcpy, host_to_device
49 use design, only: design_t
50 use math, only: rzero
51 use simulation_m, only: simulation_t
52 use json_module, only: json_file
54 use vector, only: vector_t
55 use math, only: copy
56 use device_math, only: device_copy
57 use registry, only: neko_registry
60 use field_math, only: field_rzero
61 use json_utils, only: json_get, json_get_or_default, json_get
62 use utils, only: neko_error
63 use comm, only: neko_comm
64 implicit none
65 private
66
68 type, extends(design_t), public :: brinkman_design_t
69 private
70
71 ! TODO
72 ! in the future make this a derived type of a `design_variable`
73 ! type, public, extends(design_variable_t) :: brinkman_design_t
75 type(field_t), pointer :: design_indicator
76
78 ! TODO
79 ! NOTE: Tim, right now we only have the brinkman term
80 ! in the future we may map to other coeeficients for other equations...
81 ! in which case this should be a field list
82 !
83 ! or as I describe below, we also have multiple constraints,
84 ! so a list-of-lists may be the correct way forward
85 type(field_t), pointer :: brinkman_amplitude
86
87 ! NOTE:
88 ! again, we have to be so clear with nomenclature.
89 ! If we have an objective function F.
90 ! I also like to use \rho to describe the design_indicator
91 !
92 ! and let's say, for completness, we're doing conjugate heat transfer,
93 ! so we have to map
94 ! to 3 coefficients, \chi, C and \kappa.
95 !
96 ! then we will perform 3 mapping,
97 ! \rho -> \chi
98 ! \rho -> C
99 ! \rho -> \kappa
100 !
101 ! Then during the forward/adjoint looping there will be an additional
102 ! object, the "objective_function" object that will be responsible for
103 ! computing the sensitivity of the objective function with respect to the
104 ! coefficients
105 ! ie,
106 ! dF/d\chi, dF/dC and dF/d\kappa
107 !
108 ! What I'm calling "sensitivity" here, is the sensitivity with respect to
109 ! the design indicator
110 ! so dF/d\rho
111 !
112 ! so the proceedure "map_backwards" will take in the field list
113 ! dF/d\chi, dF/dC and dF/d\kappa
114 !and chain rule it's way back to
115 ! dF/d\rho
116 ! and store it here v
117 type(field_t), pointer :: sensitivity
118 ! have a field list here
119 ! type(filed_list_t), public :: constraint_sensitivity
120 ! HOWEVER !
121 ! What is a bit confusing to me... is how we'll deal with constraints.
122 !
123 ! Because in principle we could have constraints C1, C2, C3 etc
124 ! Implying we also want dC1/d\rho, dC2/d\rho etc
125 ! So this "sensitivity" should also be a list...
126 !
127 ! So I bet you'll have a nice abstract way of constructing this in the
128 ! future but I think a sort of list-of-lists would be nice.
129 !
130 ! For F:
131 ! \rho -> \tild{\rho} -> \chi
132 ! \rho -> \tild{\rho} -> C
133 ! \rho -> \tild{\rho} -> \kappa
134 !
135 ! For C1:
136 ! \rho -> \tild{\rho} (eg for a volume constraint or so)
137 !
138 ! For C2:
139 ! \rho -> \tild{\rho} (for something else)
140 !
141 ! etc..
142 !
143 ! So now we have multiple instances of the "objective" type,
144 ! each for F, C1, C2 etc
145 ! each of them can pass their dF\d_coefficents to the design
146 ! and we continue from there.
147 !
148 ! perhaps even the "objective" type is defined as a list of objectives.
149
150 ! Let's say way have a chain of two mappings
151
152 ! This is still up for discussion where the mapping will eventually live.
153 ! Maybe it makes more sense to keep the design as clean as possible,
154 ! and then the resposibility of converting the design into coefficients
155 ! used in the PDE being solved would be the responsibility of the module
156 ! used for solving the PDE, ie, the simulation.
157 ! Then the simulation would responsiblity for asking "do we have a scalar?"
158 ! ok, let me figure out how to map this design into a coeffient in the
159 ! scalar equation.
160 !
161 ! I guess then the mapping backwards gets confusing to me.
162 ! anyway, for now, let's keep things how they are and I suppose the
163 ! intention is to have different types of designs for different types of
164 ! problems.
165
171
173 class(point_zone_t), pointer :: optimization_domain
175 logical :: has_mask
176
177 ! TODO
178 ! you also had logicals for convergence etc,
179 ! I feel they should live in "problem"
180
181 contains
182
183 ! ----------------------------------------------------------------------- !
184 ! Initializations
185
187 generic, public :: init => init_from_json_sim, init_from_components
189 procedure, pass(this), public :: init_from_json_sim => &
190 brinkman_design_init_from_json_sim
192 procedure, pass(this), public :: init_from_components => &
193 brinkman_design_init_from_components
194
196 procedure, pass(this) :: get_values => brinkman_design_get_design
198 procedure, pass(this) :: get_sensitivity => brinkman_design_get_sensitivity
199
201 procedure, pass(this) :: design_get_x => brinkman_design_get_x
203 procedure, pass(this) :: design_get_x_i => brinkman_design_get_x_i
205 procedure, pass(this) :: design_get_y => brinkman_design_get_y
207 procedure, pass(this) :: design_get_y_i => brinkman_design_get_y_i
209 procedure, pass(this) :: design_get_z => brinkman_design_get_z
211 procedure, pass(this) :: design_get_z_i => brinkman_design_get_z_i
212
214 procedure, pass(this) :: update_design => brinkman_design_update_design
215
217 ! design_indicator -> filtering -> chi
218 ! and ultimately handle mapping different coeficients!
219 procedure, pass(this) :: map_forward => brinkman_design_map_forward
221 ! d_design_indicator <- d_filtering <- d_chi
222 ! and ultimately handle mapping different coeficients!
223 procedure, pass(this) :: map_backward => brinkman_design_map_backward
224 ! TODO
225 ! maybe it would have been smarter to have a "coeficient" type,
226 ! which is just a scalar field and set of mappings going from
227 ! design_indicator -> coeficient and their corresponding chain rules
228 ! maybe also some information about what equation they live in...
229
230 ! a writer being called from outside would be nice
231 procedure, pass(this) :: write => brinkman_design_write
233 procedure, pass(this) :: set_output_counter => &
234 brinkman_design_set_output_counter
235
237 procedure, pass(this) :: free => brinkman_design_free
238 ! TODO
239 ! I'm not sure who owns the optimizer...
240 ! but it would make sense to have it in here so you provide it
241 ! with dF/d_design_indicator and it updates itself.
242 ! procedure, pass(this) :: update => brinkman_design_update_design
243 end type brinkman_design_t
244
245contains
246
248 subroutine brinkman_design_init_from_json_sim(this, parameters, simulation)
249 class(brinkman_design_t), intent(inout) :: this
250 type(json_file), intent(inout) :: parameters
251 type(simulation_t), intent(inout) :: simulation
252 type(json_file) :: json_subdict
253 character(len=:), allocatable :: domain_name, domain_type, name
254 character(len=:), allocatable :: output_precision_str
255 logical :: dealias, verbose_design, verbose_sensitivity
256 integer :: output_precision
257
258 call json_get_or_default(parameters, 'name', name, 'Brinkman Design')
259 call json_get_or_default(parameters, 'domain.type', domain_type, 'full')
260 call json_get_or_default(parameters, 'dealias', dealias, .true.)
261 call json_get_or_default(parameters, 'verbose_design', verbose_design, &
262 .false.)
263 call json_get_or_default(parameters, 'verbose_sensitivity', &
264 verbose_sensitivity, .false.)
265 call json_get_or_default(parameters, 'output_precision', &
266 output_precision_str, 'sp')
267
268 select case (trim(output_precision_str))
269 case ('sp', 'SP')
270 output_precision = sp
271 case ('dp', 'DP')
272 output_precision = dp
273 case default
274 call neko_error('Invalid output_precision in design.brinkman. ' // &
275 'Expected ''sp'' or ''dp''.')
276 end select
277
278 select case (trim(domain_type))
279 case ('full')
280 this%has_mask = .false.
281 case ('point_zone')
282 this%has_mask = .true.
283 call json_get(parameters, 'domain.zone_name', domain_name)
284 this%optimization_domain => &
285 neko_point_zone_registry%get_point_zone(domain_name)
286
287 case default
288 call neko_error('brinkman design only supports point_zones for ' // &
289 'optimization domain types')
290
291 end select
292
293 ! Initialize and inject into the simulation
294 call this%init_from_components(name, simulation, dealias)
295
296 ! Initialize the mapper
297 associate(coef => simulation%neko_case%fluid%c_Xh, &
298 gs => simulation%neko_case%fluid%gs_Xh)
299
300 if ('mapping' .in. parameters) then
301 call this%mapping%init_base(coef)
302 call this%mapping%add(parameters, 'mapping')
303 call this%mapping%init_output_fields(this%brinkman_amplitude, &
304 this%sensitivity, verbose_design, verbose_sensitivity, &
305 output_precision)
306 end if
307
308 if ('initial_distribution' .in. parameters) then
309 call json_get(parameters, 'initial_distribution', json_subdict)
310 call set_optimization_ic(this%design_indicator, coef, gs, &
311 json_subdict)
312 else
313 call field_rzero(this%design_indicator)
314 end if
315 end associate
316
317 ! Map to the Brinkman amplitude
318 call this%map_forward()
319
320 end subroutine brinkman_design_init_from_json_sim
321
323 subroutine brinkman_design_free(this)
324 class(brinkman_design_t), intent(inout) :: this
325
326 call this%free_base()
327 nullify(this%brinkman_amplitude)
328 nullify(this%design_indicator)
329 nullify(this%sensitivity)
330
331 end subroutine brinkman_design_free
332
333 subroutine brinkman_design_init_from_components(this, name, simulation, &
334 dealias)
335 class(brinkman_design_t), intent(inout) :: this
336 character(len=*), intent(in) :: name
337 type(simulation_t), intent(inout) :: simulation
338 logical, intent(in) :: dealias
339 integer :: n
340 type(simple_brinkman_source_term_t) :: forward_brinkman, adjoint_brinkman
341
342 associate(dof => simulation%neko_case%fluid%dm_Xh)
343
344 call neko_registry%add_field(dof, "design_indicator", .true.)
345 call neko_registry%add_field(dof, "brinkman_amplitude", .true.)
346 call neko_registry%add_field(dof, "sensitivity", .true.)
347
348 end associate
349
350 this%design_indicator => neko_registry%get_field("design_indicator")
351 this%brinkman_amplitude => neko_registry%get_field("brinkman_amplitude")
352 this%sensitivity => neko_registry%get_field("sensitivity")
353
354 ! TODO
355 ! this is where we steal basically everything in
356 ! brinkman_source_term regarding loading initial fields
357 ! for now, make it a cylinder by hand
358 this%design_indicator = 0.0_rp
359 this%brinkman_amplitude = 0.0_rp
360 this%design_indicator = 0.0_rp
361
362 ! TODO
363 ! Regarding masks and filters,
364 ! I suppose there are two ways of thinking about it:
365 ! 1) Mask first, then filter
366 ! 2) filter first, then mask
367 !
368 ! Each one can be used to achieve different results, and when we do complex
369 ! non-linear filter cascades the choice here has implications for exactly
370 ! how we define "minimum size control".
371 !
372 ! On top of this, I've only ever used spatial convolution filters before,
373 ! not PDE based filters.
374 ! With these filters you need to think of how your domain is "padded" when
375 ! you move the kernel over the boundary. Again, you can achieve different
376 ! effects depending on the choice of padding.
377 ! I don't really know if there's a similar implication for PDE based
378 ! based filters, I bet there is. (We should ask Niels and Casper) I bet it
379 ! means we need to consider the boundary of the mask as the boundary we
380 ! enforce to be Nuemann, not the boundary of the computational domain.
381 ! That sounds REALLY hard to implement...I hope it doesn't come down to
382 ! that.
383 !
384 ! We can always change this decision later, but I'm going with (1)
385 ! mask first, then filter.
386 ! The reason being, one of the purposes of the filtering is to avoid sharp
387 ! interfaces, if we filter first then mask there's a chance we have a sharp
388 ! interface on the boundary of the optimization domain.
389 ! if we mask first then filter, at least all the boundaries will be smooth.
390 if (this%has_mask) then
391 call mask_exterior_const(this%design_indicator, &
392 this%optimization_domain, 0.0_rp)
393 end if
394
395 n = this%design_indicator%dof%size()
396 call this%init_base(name, n)
397
398 ! init the simple brinkman term for the forward problem
399 call forward_brinkman%init_from_components( &
400 simulation%fluid%f_x, &
401 simulation%fluid%f_y, &
402 simulation%fluid%f_z, &
403 this%brinkman_amplitude, &
404 simulation%fluid%u, &
405 simulation%fluid%v, &
406 simulation%fluid%w, &
407 simulation%fluid%c_Xh, &
408 simulation%adjoint_fluid%c_Xh_GL, &
409 simulation%adjoint_fluid%GLL_to_GL, &
410 dealias, simulation%adjoint_fluid%scratch_GL)
411 ! append brinkman source term to the forward problem
412 call simulation%fluid%source_term%add(forward_brinkman)
413
414 ! init the simple brinkman term for the adjoint
415 call adjoint_brinkman%init_from_components( &
416 simulation%adjoint_fluid%f_adj_x, &
417 simulation%adjoint_fluid%f_adj_y, &
418 simulation%adjoint_fluid%f_adj_z, &
419 this%brinkman_amplitude, &
420 simulation%adjoint_fluid%u_adj, &
421 simulation%adjoint_fluid%v_adj, &
422 simulation%adjoint_fluid%w_adj, &
423 simulation%adjoint_fluid%c_Xh, &
424 simulation%adjoint_fluid%c_Xh_GL, &
425 simulation%adjoint_fluid%GLL_to_GL, &
426 dealias, simulation%adjoint_fluid%scratch_GL)
427 ! append brinkman source term based on design
428
429 select type (f => simulation%adjoint_fluid)
430 type is (adjoint_fluid_pnpn_t)
431 call f%source_term%add(adjoint_brinkman)
432 class default
433 end select
434
435 end subroutine brinkman_design_init_from_components
436
437
438 subroutine brinkman_design_map_forward(this)
439 class(brinkman_design_t), intent(inout) :: this
440
441 ! TODO, see previous todo about mask first, then mapping
442 if (this%has_mask) then
443 call mask_exterior_const(this%design_indicator, &
444 this%optimization_domain, 0.0_rp)
445 end if
446
447 call this%mapping%apply_forward(this%brinkman_amplitude, &
448 this%design_indicator)
449
450 end subroutine brinkman_design_map_forward
451
452 subroutine brinkman_design_get_design(this, values)
453 class(brinkman_design_t), intent(in) :: this
454 type(vector_t), intent(inout) :: values
455 integer :: n
456
457 n = this%size()
458 if (n .ne. values%size()) then
459 call neko_error('Get design: size mismatch')
460 end if
461
462 if (neko_bcknd_device .eq. 1) then
463 call device_copy(values%x_d, this%design_indicator%x_d, n)
464 else
465 call copy(values%x, this%design_indicator%x, n)
466 end if
467
468 end subroutine brinkman_design_get_design
469
470 subroutine brinkman_design_get_sensitivity(this, values)
471 class(brinkman_design_t), intent(in) :: this
472 type(vector_t), intent(inout) :: values
473 integer :: n
474
475 n = this%size()
476 if (n .ne. values%size()) then
477 call neko_error('Get design: size mismatch')
478 end if
479
480 if (neko_bcknd_device .eq. 1) then
481 call device_copy(values%x_d, this%sensitivity%x_d, n)
482 else
483 call copy(values%x, this%sensitivity%x, n)
484 end if
485
486 end subroutine brinkman_design_get_sensitivity
487
488 subroutine brinkman_design_get_x(this, x)
489 class(brinkman_design_t), intent(in) :: this
490 type(vector_t), intent(inout) :: x
491 integer :: n
492
493 n = this%size()
494 if (n .ne. x%size()) then
495 call neko_error('Get x: size mismatch')
496 end if
497
498 if (neko_bcknd_device .eq. 1) then
499 call device_copy(x%x_d, this%design_indicator%dof%x_d, n)
500 else
501 call copy(x%x, this%design_indicator%dof%x, n)
502 end if
503
504 end subroutine brinkman_design_get_x
505
506 function brinkman_design_get_x_i(this, i) result(x_i)
507 class(brinkman_design_t), intent(in) :: this
508 integer, intent(in) :: i
509 real(kind=rp) :: x_i
510 integer :: n
511
512 n = this%size()
513 if (i .lt. 1 .or. i .gt. n) then
514 call neko_error('brinkman_design_get_x_i: index out of bounds')
515 end if
516
517 x_i = this%design_indicator%dof%x(i,1,1,1)
518
519 end function brinkman_design_get_x_i
520
521 subroutine brinkman_design_get_y(this, y)
522 class(brinkman_design_t), intent(in) :: this
523 type(vector_t), intent(inout) :: y
524 integer :: n
525
526 n = this%size()
527 if (n .ne. y%size()) then
528 call neko_error('Get y: size mismatch')
529 end if
530
531 if (neko_bcknd_device .eq. 1) then
532 call device_copy(y%x_d, this%design_indicator%dof%y_d, n)
533 else
534 call copy(y%x, this%design_indicator%dof%y, n)
535 end if
536
537 end subroutine brinkman_design_get_y
538
539 function brinkman_design_get_y_i(this, i) result(y_i)
540 class(brinkman_design_t), intent(in) :: this
541 integer, intent(in) :: i
542 real(kind=rp) :: y_i
543 integer :: n
544
545 n = this%size()
546 if (i .lt. 1 .or. i .gt. n) then
547 call neko_error('brinkman_design_get_y_i: index out of bounds')
548 end if
549
550 y_i = this%design_indicator%dof%y(i,1,1,1)
551
552 end function brinkman_design_get_y_i
553
554 subroutine brinkman_design_get_z(this, z)
555 class(brinkman_design_t), intent(in) :: this
556 type(vector_t), intent(inout) :: z
557 integer :: n
558
559 n = this%size()
560 if (n .ne. z%size()) then
561 call neko_error('Get z: size mismatch')
562 end if
563
564 if (neko_bcknd_device .eq. 1) then
565 call device_copy(z%x_d, this%design_indicator%dof%z_d, n)
566 else
567 call copy(z%x, this%design_indicator%dof%z, n)
568 end if
569
570 end subroutine brinkman_design_get_z
571
572 function brinkman_design_get_z_i(this, i) result(z_i)
573 class(brinkman_design_t), intent(in) :: this
574 integer, intent(in) :: i
575 real(kind=rp) :: z_i
576 integer :: n
577
578 n = this%size()
579 if (i .lt. 1 .or. i .gt. n) then
580 call neko_error('brinkman_design_get_z_i: index out of bounds')
581 end if
582
583 z_i = this%design_indicator%dof%z(i,1,1,1)
584
585 end function brinkman_design_get_z_i
586
587 subroutine brinkman_design_update_design(this, values)
588 class(brinkman_design_t), intent(inout) :: this
589 type(vector_t), intent(inout) :: values
590 integer :: n
591
592 n = this%size()
593 if (neko_bcknd_device .eq. 1) then
594 call device_copy(this%design_indicator%x_d, values%x_d, n)
595 else
596 call copy(this%design_indicator%x, values%x, n)
597 end if
598
599 call this%map_forward()
600
601 if (neko_bcknd_device .eq. 1) then
602 call device_copy(values%x_d, this%design_indicator%x_d, n)
603 else
604 call copy(values%x, this%design_indicator%x, n)
605 end if
606
607 end subroutine brinkman_design_update_design
608
609 subroutine brinkman_design_map_backward(this, sensitivity)
610 class(brinkman_design_t), intent(inout) :: this
611 type(vector_t), intent(in) :: sensitivity
612 type(field_t), pointer :: tmp_fld
613 integer :: temp_index
614
615 call neko_scratch_registry%request(tmp_fld, temp_index, .false.)
616
617 call vector_to_field(tmp_fld, sensitivity)
618
619 call this%mapping%apply_backward(this%sensitivity, tmp_fld)
620
621 if (this%has_mask) then
622 call mask_exterior_const(this%sensitivity, this%optimization_domain, &
623 0.0_rp)
624 end if
625
626 call neko_scratch_registry%relinquish_field(temp_index)
627
628 end subroutine brinkman_design_map_backward
629
630 subroutine brinkman_design_write(this, idx)
631 class(brinkman_design_t), intent(inout) :: this
632 integer, intent(in) :: idx
633
634 call this%mapping%write_design(idx)
635 call this%mapping%write_sensitivity(idx)
636
637 end subroutine brinkman_design_write
638
639 subroutine brinkman_design_set_output_counter(this, idx)
640 class(brinkman_design_t), intent(inout) :: this
641 integer, intent(in) :: idx
642
643 call this%mapping%set_output_counter(idx)
644
645 end subroutine brinkman_design_set_output_counter
646
647end module brinkman_design
Adjoint Pn/Pn formulation.
Implements the design_t.
Definition design.f90:36
Implements the mapping_handler_t type.
Mappings to be applied to a scalar field.
Definition mapping.f90:36
Some common Masking operations we may need.
Definition mask_ops.f90:36
Contains extensions to the neko library required to run the topology optimization code.
Definition neko_ext.f90:42
subroutine, public field_to_vector(vector, field)
Field to vector.
Definition neko_ext.f90:409
subroutine, public vector_to_field(field, vector)
Vector to field.
Definition neko_ext.f90:385
Optimization initial condition.
Implements the simple_brinkman_source_term_t type.
Implements the steady_problem_t type.
A topology optimization design variable.
An abstract design type.
Definition design.f90:54
Abstract class for handling mapping_cascade.