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