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
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 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 ! afield writer would be nice too
182 type(fld_file_output_t), private :: output
183
184 contains
185
186 ! ----------------------------------------------------------------------- !
187 ! Initializations
188
190 generic, public :: init => init_from_json_sim, init_from_components
192 procedure, pass(this), public :: init_from_json_sim => &
193 brinkman_design_init_from_json_sim
195 procedure, pass(this), public :: init_from_components => &
196 brinkman_design_init_from_components
197
199 procedure, pass(this) :: get_values => brinkman_design_get_design
201 procedure, pass(this) :: get_sensitivity => brinkman_design_get_sensitivity
202
204 procedure, pass(this) :: design_get_x => brinkman_design_get_x
206 procedure, pass(this) :: design_get_x_i => brinkman_design_get_x_i
208 procedure, pass(this) :: design_get_y => brinkman_design_get_y
210 procedure, pass(this) :: design_get_y_i => brinkman_design_get_y_i
212 procedure, pass(this) :: design_get_z => brinkman_design_get_z
214 procedure, pass(this) :: design_get_z_i => brinkman_design_get_z_i
215
217 procedure, pass(this) :: update_design => brinkman_design_update_design
218
220 ! design_indicator -> filtering -> chi
221 ! and ultimately handle mapping different coeficients!
222 procedure, pass(this) :: map_forward => brinkman_design_map_forward
224 ! d_design_indicator <- d_filtering <- d_chi
225 ! and ultimately handle mapping different coeficients!
226 procedure, pass(this) :: map_backward => brinkman_design_map_backward
227 ! TODO
228 ! maybe it would have been smarter to have a "coeficient" type,
229 ! which is just a scalar field and set of mappings going from
230 ! design_indicator -> coeficient and their corresponding chain rules
231 ! maybe also some information about what equation they live in...
232
233 ! a writer being called from outside would be nice
234 procedure, pass(this) :: write => brinkman_design_write
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 logical :: dealias
255
256 call json_get_or_default(parameters, 'name', name, 'Brinkman Design')
257 call json_get_or_default(parameters, 'domain.type', domain_type, 'full')
258 call json_get_or_default(parameters, 'dealias', dealias, .true.)
259
260 select case (trim(domain_type))
261 case ('full')
262 this%has_mask = .false.
263 case ('point_zone')
264 this%has_mask = .true.
265 call json_get(parameters, 'domain.zone_name', domain_name)
266 this%optimization_domain => &
267 neko_point_zone_registry%get_point_zone(domain_name)
268
269 case default
270 call neko_error('brinkman design only supports point_zones for ' // &
271 'optimization domain types')
272
273 end select
274
275 ! Initialize and inject into the simulation
276 call this%init_from_components(name, simulation, dealias)
277
278 ! Initialize the mapper
279 associate(coef => simulation%neko_case%fluid%c_Xh, &
280 gs => simulation%neko_case%fluid%gs_Xh)
281
282 if ('mapping' .in. parameters) then
283 call this%mapping%init_base(coef)
284 call this%mapping%add(parameters, 'mapping')
285 end if
286
287 if ('initial_distribution' .in. parameters) then
288 call json_get(parameters, 'initial_distribution', json_subdict)
289 call set_optimization_ic(this%design_indicator, coef, gs, &
290 json_subdict)
291 else
292 call field_rzero(this%design_indicator)
293 end if
294 end associate
295
296 ! Map to the Brinkman amplitude
297 call this%map_forward()
298
299 end subroutine brinkman_design_init_from_json_sim
300
302 subroutine brinkman_design_free(this)
303 class(brinkman_design_t), intent(inout) :: this
304
305 call this%free_base()
306 nullify(this%brinkman_amplitude)
307 nullify(this%design_indicator)
308 nullify(this%sensitivity)
309
310 end subroutine brinkman_design_free
311
312 subroutine brinkman_design_init_from_components(this, name, simulation, &
313 dealias)
314 class(brinkman_design_t), intent(inout) :: this
315 character(len=*), intent(in) :: name
316 type(simulation_t), intent(inout) :: simulation
317 logical, intent(in) :: dealias
318 integer :: n
319 type(simple_brinkman_source_term_t) :: forward_brinkman, adjoint_brinkman
320
321 associate(dof => simulation%neko_case%fluid%dm_Xh)
322
323 call neko_registry%add_field(dof, "design_indicator", .true.)
324 call neko_registry%add_field(dof, "brinkman_amplitude", .true.)
325 call neko_registry%add_field(dof, "sensitivity", .true.)
326
327 end associate
328
329 this%design_indicator => neko_registry%get_field("design_indicator")
330 this%brinkman_amplitude => neko_registry%get_field("brinkman_amplitude")
331 this%sensitivity => neko_registry%get_field("sensitivity")
332
333 ! TODO
334 ! this is where we steal basically everything in
335 ! brinkman_source_term regarding loading initial fields
336 ! for now, make it a cylinder by hand
337 this%design_indicator = 0.0_rp
338 this%brinkman_amplitude = 0.0_rp
339 this%design_indicator = 0.0_rp
340
341 ! TODO
342 ! Regarding masks and filters,
343 ! I suppose there are two ways of thinking about it:
344 ! 1) Mask first, then filter
345 ! 2) filter first, then mask
346 !
347 ! Each one can be used to achieve different results, and when we do complex
348 ! non-linear filter cascades the choice here has implications for exactly
349 ! how we define "minimum size control".
350 !
351 ! On top of this, I've only ever used spatial convolution filters before,
352 ! not PDE based filters.
353 ! With these filters you need to think of how your domain is "padded" when
354 ! you move the kernel over the boundary. Again, you can achieve different
355 ! effects depending on the choice of padding.
356 ! I don't really know if there's a similar implication for PDE based
357 ! based filters, I bet there is. (We should ask Niels and Casper) I bet it
358 ! means we need to consider the boundary of the mask as the boundary we
359 ! enforce to be Nuemann, not the boundary of the computational domain.
360 ! That sounds REALLY hard to implement...I hope it doesn't come down to
361 ! that.
362 !
363 ! We can always change this decision later, but I'm going with (1)
364 ! mask first, then filter.
365 ! The reason being, one of the purposes of the filtering is to avoid sharp
366 ! interfaces, if we filter first then mask there's a chance we have a sharp
367 ! interface on the boundary of the optimization domain.
368 ! if we mask first then filter, at least all the boundaries will be smooth.
369 if (this%has_mask) then
370 call mask_exterior_const(this%design_indicator, &
371 this%optimization_domain, 0.0_rp)
372 end if
373
374 ! a field writer would be nice to output
375 ! - design indicator (\rho)
376 ! - mapped design (\chi)
377 ! - sensitivity (dF/d\chi)
378 ! TODO
379 ! do this properly with JSON
380 ! TODO
381 ! obviously when we do the mappings properly, to many coeficients,
382 ! we'll also have to modify this
383 call this%output%init(sp, 'design', 3)
384 call this%output%fields%assign_to_field(1, this%design_indicator)
385 call this%output%fields%assign_to_field(2, this%brinkman_amplitude)
386 call this%output%fields%assign_to_field(3, this%sensitivity)
387
388 n = this%design_indicator%dof%size()
389 call this%init_base(name, n)
390
391 ! init the simple brinkman term for the forward problem
392 call forward_brinkman%init_from_components( &
393 simulation%fluid%f_x, &
394 simulation%fluid%f_y, &
395 simulation%fluid%f_z, &
396 this%brinkman_amplitude, &
397 simulation%fluid%u, &
398 simulation%fluid%v, &
399 simulation%fluid%w, &
400 simulation%fluid%c_Xh, &
401 simulation%adjoint_fluid%c_Xh_GL, &
402 simulation%adjoint_fluid%GLL_to_GL, &
403 dealias, simulation%adjoint_fluid%scratch_GL)
404 ! append brinkman source term to the forward problem
405 call simulation%fluid%source_term%add(forward_brinkman)
406
407 ! init the simple brinkman term for the adjoint
408 call adjoint_brinkman%init_from_components( &
409 simulation%adjoint_fluid%f_adj_x, &
410 simulation%adjoint_fluid%f_adj_y, &
411 simulation%adjoint_fluid%f_adj_z, &
412 this%brinkman_amplitude, &
413 simulation%adjoint_fluid%u_adj, &
414 simulation%adjoint_fluid%v_adj, &
415 simulation%adjoint_fluid%w_adj, &
416 simulation%adjoint_fluid%c_Xh, &
417 simulation%adjoint_fluid%c_Xh_GL, &
418 simulation%adjoint_fluid%GLL_to_GL, &
419 dealias, simulation%adjoint_fluid%scratch_GL)
420 ! append brinkman source term based on design
421
422 select type (f => simulation%adjoint_fluid)
423 type is (adjoint_fluid_pnpn_t)
424 call f%source_term%add(adjoint_brinkman)
425 class default
426 end select
427
428 end subroutine brinkman_design_init_from_components
429
430
431 subroutine brinkman_design_map_forward(this)
432 class(brinkman_design_t), intent(inout) :: this
433
434 ! TODO, see previous todo about mask first, then mapping
435 if (this%has_mask) then
436 call mask_exterior_const(this%design_indicator, &
437 this%optimization_domain, 0.0_rp)
438 end if
439
440 call this%mapping%apply_forward(this%brinkman_amplitude, &
441 this%design_indicator)
442
443 end subroutine brinkman_design_map_forward
444
445 subroutine brinkman_design_get_design(this, values)
446 class(brinkman_design_t), intent(in) :: this
447 type(vector_t), intent(inout) :: values
448 integer :: n
449
450 n = this%size()
451 if (n .ne. values%size()) then
452 call neko_error('Get design: size mismatch')
453 end if
454
455 if (neko_bcknd_device .eq. 1) then
456 call device_copy(values%x_d, this%design_indicator%x_d, n)
457 else
458 call copy(values%x, this%design_indicator%x, n)
459 end if
460
461 end subroutine brinkman_design_get_design
462
463 subroutine brinkman_design_get_sensitivity(this, values)
464 class(brinkman_design_t), intent(in) :: this
465 type(vector_t), intent(inout) :: values
466 integer :: n
467
468 n = this%size()
469 if (n .ne. values%size()) then
470 call neko_error('Get design: size mismatch')
471 end if
472
473 if (neko_bcknd_device .eq. 1) then
474 call device_copy(values%x_d, this%sensitivity%x_d, n)
475 else
476 call copy(values%x, this%sensitivity%x, n)
477 end if
478
479 end subroutine brinkman_design_get_sensitivity
480
481 subroutine brinkman_design_get_x(this, x)
482 class(brinkman_design_t), intent(in) :: this
483 type(vector_t), intent(inout) :: x
484 integer :: n
485
486 n = this%size()
487 if (n .ne. x%size()) then
488 call neko_error('Get x: size mismatch')
489 end if
490
491 if (neko_bcknd_device .eq. 1) then
492 call device_copy(x%x_d, this%design_indicator%dof%x_d, n)
493 else
494 call copy(x%x, this%design_indicator%dof%x, n)
495 end if
496
497 end subroutine brinkman_design_get_x
498
499 function brinkman_design_get_x_i(this, i) result(x_i)
500 class(brinkman_design_t), intent(in) :: this
501 integer, intent(in) :: i
502 real(kind=rp) :: x_i
503 integer :: n
504
505 n = this%size()
506 if (i .lt. 1 .or. i .gt. n) then
507 call neko_error('brinkman_design_get_x_i: index out of bounds')
508 end if
509
510 x_i = this%design_indicator%dof%x(i,1,1,1)
511
512 end function brinkman_design_get_x_i
513
514 subroutine brinkman_design_get_y(this, y)
515 class(brinkman_design_t), intent(in) :: this
516 type(vector_t), intent(inout) :: y
517 integer :: n
518
519 n = this%size()
520 if (n .ne. y%size()) then
521 call neko_error('Get y: size mismatch')
522 end if
523
524 if (neko_bcknd_device .eq. 1) then
525 call device_copy(y%x_d, this%design_indicator%dof%y_d, n)
526 else
527 call copy(y%x, this%design_indicator%dof%y, n)
528 end if
529
530 end subroutine brinkman_design_get_y
531
532 function brinkman_design_get_y_i(this, i) result(y_i)
533 class(brinkman_design_t), intent(in) :: this
534 integer, intent(in) :: i
535 real(kind=rp) :: y_i
536 integer :: n
537
538 n = this%size()
539 if (i .lt. 1 .or. i .gt. n) then
540 call neko_error('brinkman_design_get_y_i: index out of bounds')
541 end if
542
543 y_i = this%design_indicator%dof%y(i,1,1,1)
544
545 end function brinkman_design_get_y_i
546
547 subroutine brinkman_design_get_z(this, z)
548 class(brinkman_design_t), intent(in) :: this
549 type(vector_t), intent(inout) :: z
550 integer :: n
551
552 n = this%size()
553 if (n .ne. z%size()) then
554 call neko_error('Get z: size mismatch')
555 end if
556
557 if (neko_bcknd_device .eq. 1) then
558 call device_copy(z%x_d, this%design_indicator%dof%z_d, n)
559 else
560 call copy(z%x, this%design_indicator%dof%z, n)
561 end if
562
563 end subroutine brinkman_design_get_z
564
565 function brinkman_design_get_z_i(this, i) result(z_i)
566 class(brinkman_design_t), intent(in) :: this
567 integer, intent(in) :: i
568 real(kind=rp) :: z_i
569 integer :: n
570
571 n = this%size()
572 if (i .lt. 1 .or. i .gt. n) then
573 call neko_error('brinkman_design_get_z_i: index out of bounds')
574 end if
575
576 z_i = this%design_indicator%dof%z(i,1,1,1)
577
578 end function brinkman_design_get_z_i
579
580 subroutine brinkman_design_update_design(this, values)
581 class(brinkman_design_t), intent(inout) :: this
582 type(vector_t), intent(inout) :: values
583 integer :: n
584
585 n = this%size()
586 if (neko_bcknd_device .eq. 1) then
587 call device_copy(this%design_indicator%x_d, values%x_d, n)
588 else
589 call copy(this%design_indicator%x, values%x, n)
590 end if
591
592 call this%map_forward()
593
594 if (neko_bcknd_device .eq. 1) then
595 call device_copy(values%x_d, this%design_indicator%x_d, n)
596 else
597 call copy(values%x, this%design_indicator%x, n)
598 end if
599
600 end subroutine brinkman_design_update_design
601
602 subroutine brinkman_design_map_backward(this, sensitivity)
603 class(brinkman_design_t), intent(inout) :: this
604 type(vector_t), intent(in) :: sensitivity
605 type(field_t), pointer :: tmp_fld
606 integer :: temp_index
607
608 call neko_scratch_registry%request(tmp_fld, temp_index, .false.)
609
610 call vector_to_field(tmp_fld, sensitivity)
611
612 call this%mapping%apply_backward(this%sensitivity, tmp_fld)
613
614 if (this%has_mask) then
615 call mask_exterior_const(this%sensitivity, this%optimization_domain, &
616 0.0_rp)
617 end if
618
619 call neko_scratch_registry%relinquish_field(temp_index)
620
621 end subroutine brinkman_design_map_backward
622
623 subroutine brinkman_design_write(this, idx)
624 class(brinkman_design_t), intent(inout) :: this
625 integer, intent(in) :: idx
626
627 call this%output%sample(real(idx, kind=rp))
628
629 end subroutine brinkman_design_write
630
631end 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.