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