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 field_registry, only: neko_field_registry
59 use field_math, only: field_rzero
60 use json_utils, only: json_get, json_get_or_default, json_extract_object
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 :: if_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
198
200 procedure, pass(this) :: get_x => brinkman_design_get_x
202 procedure, pass(this) :: get_y => brinkman_design_get_y
204 procedure, pass(this) :: get_z => brinkman_design_get_z
205
207 procedure, pass(this) :: update_design => brinkman_design_update_design
208
210 ! design_indicator -> filtering -> chi
211 ! and ultimately handle mapping different coeficients!
212 procedure, pass(this) :: map_forward => brinkman_design_map_forward
214 ! d_design_indicator <- d_filtering <- d_chi
215 ! and ultimately handle mapping different coeficients!
216 procedure, pass(this) :: map_backward => brinkman_design_map_backward
217 ! TODO
218 ! maybe it would have been smarter to have a "coeficient" type,
219 ! which is just a scalar field and set of mappings going from
220 ! design_indicator -> coeficient and their corresponding chain rules
221 ! maybe also some information about what equation they live in...
222
223 ! a writer being called from outside would be nice
224 procedure, pass(this) :: write => brinkman_design_write
225
227 procedure, pass(this) :: free => brinkman_design_free
228 ! TODO
229 ! I'm not sure who owns the optimizer...
230 ! but it would make sense to have it in here so you provide it
231 ! with dF/d_design_indicator and it updates itself.
232 ! procedure, pass(this) :: update => brinkman_design_update_design
233 end type brinkman_design_t
234
235contains
236
238 subroutine brinkman_design_init_from_json_sim(this, parameters, simulation)
239 class(brinkman_design_t), intent(inout) :: this
240 type(json_file), intent(inout) :: parameters
241 type(simulation_t), intent(inout) :: simulation
242 type(json_file) :: json_subdict
243 character(len=:), allocatable :: domain_name, domain_type
244
245 ! Initialize the optimization domain
246 if (parameters%valid_path('optimization.domain')) then
247 call json_get(parameters, 'optimization.domain.type', domain_type)
248 select case (trim(domain_type))
249 case ('point_zone')
250 this%if_mask = .true.
251 call json_get(parameters, 'optimization.domain.zone_name', &
252 domain_name)
253 this%optimization_domain => &
254 neko_point_zone_registry%get_point_zone(domain_name)
255
256 case default
257 call neko_error('brinkman design only supports point_zones for&
258 & optimization domain types')
259
260 end select
261 else
262 this%if_mask = .false.
263 end if
264
265 ! Initialize and inject into the simulation
266 call this%init_from_components(simulation)
267
268 ! Initialize the mapper
269 associate(coef => simulation%neko_case%fluid%c_Xh, &
270 gs => simulation%neko_case%fluid%gs_Xh)
271 call this%mapping%init_base(coef)
272 call this%mapping%add(parameters, 'optimization.design.mapping')
273
274 if (parameters%valid_path(&
275 'optimization.design.initial_distribution')) then
276 call json_extract_object(parameters, &
277 'optimization.design.initial_distribution', json_subdict)
278 call set_optimization_ic(this%design_indicator, coef, gs, &
279 json_subdict)
280 else
281 call field_rzero(this%design_indicator)
282 end if
283 end associate
284
285 ! Map to the Brinkman amplitude
286 call this%map_forward()
287
288 end subroutine brinkman_design_init_from_json_sim
289
291 subroutine brinkman_design_free(this)
292 class(brinkman_design_t), intent(inout) :: this
293
294 call this%free_base()
295 call this%brinkman_amplitude%free()
296 call this%design_indicator%free()
297 call this%sensitivity%free()
298
299 end subroutine brinkman_design_free
300
301 subroutine brinkman_design_init_from_components(this, simulation)
302 class(brinkman_design_t), intent(inout) :: this
303 type(simulation_t), intent(inout) :: simulation
304 integer :: n, i
305 type(simple_brinkman_source_term_t) :: forward_brinkman, adjoint_brinkman
306
307 associate(dof => simulation%neko_case%fluid%dm_Xh)
308
309 call neko_field_registry%add_field(dof, "design_indicator", .true.)
310 call neko_field_registry%add_field(dof, "brinkman_amplitude", .true.)
311 call neko_field_registry%add_field(dof, "sensitivity", .true.)
312
313 end associate
314
315 this%design_indicator => &
316 neko_field_registry%get_field("design_indicator")
317 this%brinkman_amplitude => &
318 neko_field_registry%get_field("brinkman_amplitude")
319 this%sensitivity => &
320 neko_field_registry%get_field("sensitivity")
321
322 ! TODO
323 ! this is where we steal basically everything in
324 ! brinkman_source_term regarding loading initial fields
325 ! for now, make it a cylinder by hand
326 this%design_indicator = 0.0_rp
327 this%brinkman_amplitude = 0.0_rp
328 this%design_indicator%x = 0.0_rp
329
330 n = this%design_indicator%dof%size()
331 ! This is probably getting fixed in tim's PR anyway, otherwise I'll fix it.
332 do i = 1, n
333 this%design_indicator%x(i,1,1,1) = 0.0_rp
334 end do
335
336 ! again this will be handled better in the future...
337 if (neko_bcknd_device .eq. 1) then
338 call device_memcpy(this%design_indicator%x, &
339 this%design_indicator%x_d, n, &
340 host_to_device, sync = .false.)
341 end if
342
343 ! TODO
344 ! Regarding masks and filters,
345 ! I suppose there are two ways of thinking about it:
346 ! 1) Mask first, then filter
347 ! 2) filter first, then mask
348 !
349 ! Each one can be used to achieve different results, and when we do complex
350 ! non-linear filter cascades the choice here has implications for exactly
351 ! how we define "minimum size control".
352 !
353 ! On top of this, I've only ever used spatial convolution filters before,
354 ! not PDE based filters.
355 ! With these filters you need to think of how your domain is "padded" when
356 ! you move the kernel over the boundary. Again, you can achieve different
357 ! effects depending on the choice of padding.
358 ! I don't really know if there's a similar implication for PDE based
359 ! based filters, I bet there is. (We should ask Niels and Casper) I bet it
360 ! means we need to consider the boundary of the mask as the boundary we
361 ! enforce to be Nuemann, not the boundary of the computational domain.
362 ! That sounds REALLY hard to implement...I hope it doesn't come down to
363 ! that.
364 !
365 ! We can always change this decision later, but I'm going with (1)
366 ! mask first, then filter.
367 ! The reason being, one of the purposes of the filtering is to avoid sharp
368 ! interfaces, if we filter first then mask there's a chance we have a sharp
369 ! interface on the boundary of the optimization domain.
370 ! if we mask first then filter, at least all the boundaries will be smooth.
371 if (this%if_mask) then
372 call mask_exterior_const(this%design_indicator, &
373 this%optimization_domain, 0.0_rp)
374 end if
375
376 ! a field writer would be nice to output
377 ! - design indicator (\rho)
378 ! - mapped design (\chi)
379 ! - sensitivity (dF/d\chi)
380 ! TODO
381 ! do this properly with JSON
382 ! TODO
383 ! obviously when we do the mappings properly, to many coeficients,
384 ! we'll also have to modify this
385 call this%output%init(sp, 'design', 3)
386 call this%output%fields%assign_to_field(1, this%design_indicator)
387 call this%output%fields%assign_to_field(2, this%brinkman_amplitude)
388 call this%output%fields%assign_to_field(3, this%sensitivity)
389
390 call this%init_base(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 ! 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 ! append brinkman source term based on design
416
417 select type (f => simulation%adjoint_fluid)
418 type is (adjoint_fluid_pnpn_t)
419 call f%source_term%add(adjoint_brinkman)
420 class default
421 end select
422
423 end subroutine brinkman_design_init_from_components
424
425
426 subroutine brinkman_design_map_forward(this)
427 class(brinkman_design_t), intent(inout) :: this
428
429 ! TODO, see previous todo about mask first, then mapping
430 if (this%if_mask) then
431 call mask_exterior_const(this%design_indicator, &
432 this%optimization_domain, 0.0_rp)
433 end if
434
435 call this%mapping%apply_forward(this%brinkman_amplitude, &
436 this%design_indicator)
437
438 end subroutine brinkman_design_map_forward
439
440 function brinkman_design_get_design(this) result(values)
441 class(brinkman_design_t), intent(in) :: this
442 type(vector_t) :: values
443 integer :: n
444
445 n = this%size()
446 call values%init(n)
447 call copy(values%x, this%design_indicator%x, n)
448 if (neko_bcknd_device .eq. 1) then
449 call device_copy(values%x_d, this%design_indicator%x_d, n)
450 end if
451
452 end function brinkman_design_get_design
453
454 function brinkman_design_get_x(this) result(x)
455 class(brinkman_design_t), intent(in) :: this
456 type(vector_t) :: x
457 integer :: n
458
459 n = this%size()
460 call x%init(n)
461 call copy(x%x, this%design_indicator%dof%x, n)
462 if (neko_bcknd_device .eq. 1) then
463 call device_copy(x%x_d, this%design_indicator%dof%x_d, n)
464 end if
465
466 end function brinkman_design_get_x
467
468 function brinkman_design_get_y(this) result(y)
469 class(brinkman_design_t), intent(in) :: this
470 type(vector_t) :: y
471 integer :: n
472
473 n = this%size()
474 call y%init(n)
475 call copy(y%x, this%design_indicator%dof%y, n)
476 if (neko_bcknd_device .eq. 1) then
477 call device_copy(y%x_d, this%design_indicator%dof%y_d, n)
478 end if
479
480 end function brinkman_design_get_y
481
482 function brinkman_design_get_z(this) result(z)
483 class(brinkman_design_t), intent(in) :: this
484 type(vector_t) :: z
485 integer :: n
486
487 n = this%size()
488 call z%init(n)
489 call copy(z%x, this%design_indicator%dof%z, n)
490 if (neko_bcknd_device .eq. 1) then
491 call device_copy(z%x_d, this%design_indicator%dof%z_d, n)
492 end if
493
494 end function brinkman_design_get_z
495
496 subroutine brinkman_design_update_design(this, values)
497 class(brinkman_design_t), intent(inout) :: this
498 type(vector_t), intent(inout) :: values
499 integer :: n
500
501 n = this%size()
502 call copy(this%design_indicator%x, values%x, n)
503 if (neko_bcknd_device .eq. 1) then
504 call device_copy(this%design_indicator%x_d, values%x_d, n)
505 end if
506
507 call this%map_forward()
508
509 call copy(values%x, this%design_indicator%x, n)
510 if (neko_bcknd_device .eq. 1) then
511 call device_copy(values%x_d, this%design_indicator%x_d, n)
512 end if
513
514 end subroutine brinkman_design_update_design
515
516 subroutine brinkman_design_map_backward(this, sensitivity)
517 class(brinkman_design_t), intent(inout) :: this
518 type(vector_t), intent(in) :: sensitivity
519 type(field_t), pointer :: tmp_fld
520 integer :: temp_indices(1)
521
522 call neko_scratch_registry%request_field(tmp_fld, temp_indices(1))
523
524 call vector_to_field(tmp_fld, sensitivity)
525
526 call this%mapping%apply_backward(this%sensitivity, tmp_fld)
527
528 ! TODO
529 ! DELETE THIS LATER
530 !
531 ! When Abbas writes the interface for the optimization
532 ! module this may be a moot point, because we would only really collect
533 ! the sensitivity of the design variables inside the mask.
534 !
535 ! Note for Abbas,
536 ! I'm NOT doing this because I'm too lazy and I just need masks so I can
537 ! test something in the passive scalar.
538 if (this%if_mask) then
539 call mask_exterior_const(this%sensitivity, this%optimization_domain, &
540 0.0_rp)
541 end if
542
543 call neko_scratch_registry%relinquish_field(temp_indices)
544
545 end subroutine brinkman_design_map_backward
546
547 subroutine brinkman_design_write(this, idx)
548 class(brinkman_design_t), intent(inout) :: this
549 integer, intent(in) :: idx
550
551 call this%output%sample(real(idx, kind=rp))
552
553 end subroutine brinkman_design_write
554
555end 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:258
subroutine, public vector_to_field(field, vector)
Vector to field.
Definition neko_ext.f90:234
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:50
Abstract class for handling mapping_cascade.