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
40 use scratch_registry, only: neko_scratch_registry
41 use fld_file_output, only: fld_file_output_t
42 use point_zone_registry, only: neko_point_zone_registry
43 use point_zone, only: point_zone_t
45 use neko_config, only: neko_bcknd_device
46 use device, only: device_memcpy, host_to_device
47 use design, only: design_t
48 use math, only: rzero
49 use simulation_m, only: simulation_t
50 use json_module, only: json_file
52 use vector, only: vector_t
53 use math, only: copy
54 use device_math, only: device_copy
55 use field_registry, only: neko_field_registry
58 use field_math, only: field_rzero
59 use json_utils, only: json_get, json_get_or_default, json_extract_object
60 use utils, only: neko_error
61 implicit none
62 private
63
65 type, extends(design_t), public :: brinkman_design_t
66 private
67
68 ! TODO
69 ! in the future make this a derived type of a `design_variable`
70 ! type, public, extends(design_variable_t) :: brinkman_design_t
72 type(field_t), pointer :: design_indicator
73
75 ! TODO
76 ! NOTE: Tim, right now we only have the brinkman term
77 ! in the future we may map to other coeeficients for other equations...
78 ! in which case this should be a field list
79 !
80 ! or as I describe below, we also have multiple constraints,
81 ! so a list-of-lists may be the correct way forward
82 type(field_t), pointer :: brinkman_amplitude
83
84 ! NOTE:
85 ! again, we have to be so clear with nomenclature.
86 ! If we have an objective function F.
87 ! I also like to use \rho to describe the design_indicator
88 !
89 ! and let's say, for completness, we're doing conjugate heat transfer,
90 ! so we have to map
91 ! to 3 coefficients, \chi, C and \kappa.
92 !
93 ! then we will perform 3 mapping,
94 ! \rho -> \chi
95 ! \rho -> C
96 ! \rho -> \kappa
97 !
98 ! Then during the forward/adjoint looping there will be an additional
99 ! object, the "objective_function" object that will be responsible for
100 ! computing the sensitivity of the objective function with respect to the
101 ! coefficients
102 ! ie,
103 ! dF/d\chi, dF/dC and dF/d\kappa
104 !
105 ! What I'm calling "sensitivity" here, is the sensitivity with respect to
106 ! the design indicator
107 ! so dF/d\rho
108 !
109 ! so the proceedure "map_backwards" will take in the field list
110 ! dF/d\chi, dF/dC and dF/d\kappa
111 !and chain rule it's way back to
112 ! dF/d\rho
113 ! and store it here v
114 type(field_t), pointer :: sensitivity
115 ! have a field list here
116 ! type(filed_list_t), public :: constraint_sensitivity
117 ! HOWEVER !
118 ! What is a bit confusing to me... is how we'll deal with constraints.
119 !
120 ! Because in principle we could have constraints C1, C2, C3 etc
121 ! Implying we also want dC1/d\rho, dC2/d\rho etc
122 ! So this "sensitivity" should also be a list...
123 !
124 ! So I bet you'll have a nice abstract way of constructing this in the
125 ! future but I think a sort of list-of-lists would be nice.
126 !
127 ! For F:
128 ! \rho -> \tild{\rho} -> \chi
129 ! \rho -> \tild{\rho} -> C
130 ! \rho -> \tild{\rho} -> \kappa
131 !
132 ! For C1:
133 ! \rho -> \tild{\rho} (eg for a volume constraint or so)
134 !
135 ! For C2:
136 ! \rho -> \tild{\rho} (for something else)
137 !
138 ! etc..
139 !
140 ! So now we have multiple instances of the "objective" type,
141 ! each for F, C1, C2 etc
142 ! each of them can pass their dF\d_coefficents to the design
143 ! and we continue from there.
144 !
145 ! perhaps even the "objective" type is defined as a list of objectives.
146
147 ! Let's say way have a chain of two mappings
148
149 ! This is still up for discussion where the mapping will eventually live.
150 ! Maybe it makes more sense to keep the design as clean as possible,
151 ! and then the resposibility of converting the design into coefficients
152 ! used in the PDE being solved would be the responsibility of the module
153 ! used for solving the PDE, ie, the simulation.
154 ! Then the simulation would responsiblity for asking "do we have a scalar?"
155 ! ok, let me figure out how to map this design into a coeffient in the
156 ! scalar equation.
157 !
158 ! I guess then the mapping backwards gets confusing to me.
159 ! anyway, for now, let's keep things how they are and I suppose the
160 ! intention is to have different types of designs for different types of
161 ! problems.
162
168
170 class(point_zone_t), pointer :: optimization_domain
172 logical :: if_mask
173
174 ! TODO
175 ! you also had logicals for convergence etc,
176 ! I feel they should live in "problem"
177
178 ! afield writer would be nice too
179 type(fld_file_output_t), private :: output
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
197
199 procedure, pass(this) :: get_x => brinkman_design_get_x
201 procedure, pass(this) :: get_y => brinkman_design_get_y
203 procedure, pass(this) :: get_z => brinkman_design_get_z
204
206 procedure, pass(this) :: update_design => brinkman_design_update_design
207
209 ! design_indicator -> filtering -> chi
210 ! and ultimately handle mapping different coeficients!
211 procedure, pass(this) :: map_forward => brinkman_design_map_forward
213 ! d_design_indicator <- d_filtering <- d_chi
214 ! and ultimately handle mapping different coeficients!
215 procedure, pass(this) :: map_backward => brinkman_design_map_backward
216 ! TODO
217 ! maybe it would have been smarter to have a "coeficient" type,
218 ! which is just a scalar field and set of mappings going from
219 ! design_indicator -> coeficient and their corresponding chain rules
220 ! maybe also some information about what equation they live in...
221
222 ! a writer being called from outside would be nice
223 procedure, pass(this) :: write => brinkman_design_write
224
226 procedure, pass(this) :: free => brinkman_design_free
227 ! TODO
228 ! I'm not sure who owns the optimizer...
229 ! but it would make sense to have it in here so you provide it
230 ! with dF/d_design_indicator and it updates itself.
231 ! procedure, pass(this) :: update => brinkman_design_update_design
232 end type brinkman_design_t
233
234contains
235
237 subroutine brinkman_design_init_from_json_sim(this, parameters, simulation)
238 class(brinkman_design_t), intent(inout) :: this
239 type(json_file), intent(inout) :: parameters
240 type(simulation_t), intent(inout) :: simulation
241 type(json_file) :: json_subdict
242 character(len=:), allocatable :: domain_name, domain_type
243
244 ! Initialize the optimization domain
245 if (parameters%valid_path('optimization.domain')) then
246 call json_get(parameters, 'optimization.domain.type', domain_type)
247 select case (trim(domain_type))
248 case ('point_zone')
249 this%if_mask = .true.
250 call json_get(parameters, 'optimization.domain.zone_name', &
251 domain_name)
252 this%optimization_domain => &
253 neko_point_zone_registry%get_point_zone(domain_name)
254
255 case default
256 call neko_error('brinkman design only supports point_zones for&
257 & optimization domain types')
258
259 end select
260 else
261 this%if_mask = .false.
262 end if
263
264 ! Initialize and inject into the simulation
265 call this%init_from_components(simulation)
266
267 ! Initialize the mapper
268 associate(coef => simulation%neko_case%fluid%c_Xh, &
269 gs => simulation%neko_case%fluid%gs_Xh)
270 call this%mapping%init_base(coef)
271 call this%mapping%add(parameters, 'optimization.design.mapping')
272
273 if (parameters%valid_path(&
274 'optimization.design.initial_distribution')) then
275 call json_extract_object(parameters, &
276 'optimization.design.initial_distribution', json_subdict)
277 call set_optimization_ic(this%design_indicator, coef, gs, &
278 json_subdict)
279 else
280 call field_rzero(this%design_indicator)
281 end if
282 end associate
283
284 ! Map to the Brinkman amplitude
285 call this%map_forward()
286
287 end subroutine brinkman_design_init_from_json_sim
288
290 subroutine brinkman_design_free(this)
291 class(brinkman_design_t), intent(inout) :: this
292
293 call this%free_base()
294 call this%brinkman_amplitude%free()
295 call this%design_indicator%free()
296 call this%sensitivity%free()
297
298 end subroutine brinkman_design_free
299
300 subroutine brinkman_design_init_from_components(this, simulation)
301 class(brinkman_design_t), intent(inout) :: this
302 type(simulation_t), intent(inout) :: simulation
303 integer :: n, i
304 type(simple_brinkman_source_term_t) :: forward_brinkman, adjoint_brinkman
305
306 associate(dof => simulation%neko_case%fluid%dm_Xh)
307
308 call neko_field_registry%add_field(dof, "design_indicator", .true.)
309 call neko_field_registry%add_field(dof, "brinkman_amplitude", .true.)
310 call neko_field_registry%add_field(dof, "sensitivity", .true.)
311
312 end associate
313
314 this%design_indicator => &
315 neko_field_registry%get_field("design_indicator")
316 this%brinkman_amplitude => &
317 neko_field_registry%get_field("brinkman_amplitude")
318 this%sensitivity => &
319 neko_field_registry%get_field("sensitivity")
320
321 ! TODO
322 ! this is where we steal basically everything in
323 ! brinkman_source_term regarding loading initial fields
324 ! for now, make it a cylinder by hand
325 this%design_indicator = 0.0_rp
326 this%brinkman_amplitude = 0.0_rp
327 this%design_indicator%x = 0.0_rp
328
329 n = this%design_indicator%dof%size()
330 ! This is probably getting fixed in tim's PR anyway, otherwise I'll fix it.
331 do i = 1, n
332 this%design_indicator%x(i,1,1,1) = 0.0_rp
333 end do
334
335 ! again this will be handled better in the future...
336 if (neko_bcknd_device .eq. 1) then
337 call device_memcpy(this%design_indicator%x, &
338 this%design_indicator%x_d, n, &
339 host_to_device, sync = .false.)
340 end if
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%if_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 call this%init_base(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 ! append brinkman source term to the forward problem
402 call simulation%fluid%source_term%add(forward_brinkman)
403
404 ! init the simple brinkman term for the adjoint
405 call adjoint_brinkman%init_from_components( &
406 simulation%adjoint_fluid%f_adj_x, &
407 simulation%adjoint_fluid%f_adj_y, &
408 simulation%adjoint_fluid%f_adj_z, &
409 this%brinkman_amplitude, &
410 simulation%adjoint_fluid%u_adj, &
411 simulation%adjoint_fluid%v_adj, &
412 simulation%adjoint_fluid%w_adj, &
413 simulation%adjoint_fluid%c_Xh)
414 ! append brinkman source term based on design
415 call simulation%adjoint_fluid%source_term%add(adjoint_brinkman)
416
417 end subroutine brinkman_design_init_from_components
418
419
420 subroutine brinkman_design_map_forward(this)
421 class(brinkman_design_t), intent(inout) :: this
422
423 ! TODO, see previous todo about mask first, then mapping
424 if (this%if_mask) then
425 call mask_exterior_const(this%design_indicator, &
426 this%optimization_domain, 0.0_rp)
427 end if
428
429 call this%mapping%apply_forward(this%brinkman_amplitude, &
430 this%design_indicator)
431
432 end subroutine brinkman_design_map_forward
433
434 function brinkman_design_get_design(this) result(values)
435 class(brinkman_design_t), intent(in) :: this
436 type(vector_t) :: values
437 integer :: n
438
439 n = this%size()
440 call values%init(n)
441 call copy(values%x, this%design_indicator%x, n)
442 if (neko_bcknd_device .eq. 1) then
443 call device_copy(values%x_d, this%design_indicator%x_d, n)
444 end if
445
446 end function brinkman_design_get_design
447
448 function brinkman_design_get_x(this) result(x)
449 class(brinkman_design_t), intent(in) :: this
450 type(vector_t) :: x
451 integer :: n
452
453 n = this%size()
454 call x%init(n)
455 call copy(x%x, this%design_indicator%dof%x, n)
456 if (neko_bcknd_device .eq. 1) then
457 call device_copy(x%x_d, this%design_indicator%dof%x_d, n)
458 end if
459
460 end function brinkman_design_get_x
461
462 function brinkman_design_get_y(this) result(y)
463 class(brinkman_design_t), intent(in) :: this
464 type(vector_t) :: y
465 integer :: n
466
467 n = this%size()
468 call y%init(n)
469 call copy(y%x, this%design_indicator%dof%y, n)
470 if (neko_bcknd_device .eq. 1) then
471 call device_copy(y%x_d, this%design_indicator%dof%y_d, n)
472 end if
473
474 end function brinkman_design_get_y
475
476 function brinkman_design_get_z(this) result(z)
477 class(brinkman_design_t), intent(in) :: this
478 type(vector_t) :: z
479 integer :: n
480
481 n = this%size()
482 call z%init(n)
483 call copy(z%x, this%design_indicator%dof%z, n)
484 if (neko_bcknd_device .eq. 1) then
485 call device_copy(z%x_d, this%design_indicator%dof%z_d, n)
486 end if
487
488 end function brinkman_design_get_z
489
490 subroutine brinkman_design_update_design(this, values)
491 class(brinkman_design_t), intent(inout) :: this
492 type(vector_t), intent(inout) :: values
493 integer :: n
494
495 n = this%size()
496 call copy(this%design_indicator%x, values%x, n)
497 if (neko_bcknd_device .eq. 1) then
498 call device_copy(this%design_indicator%x_d, values%x_d, n)
499 end if
500
501 call this%map_forward()
502
503 call copy(values%x, this%design_indicator%x, n)
504 if (neko_bcknd_device .eq. 1) then
505 call device_copy(values%x_d, this%design_indicator%x_d, n)
506 end if
507
508 end subroutine brinkman_design_update_design
509
510 subroutine brinkman_design_map_backward(this, sensitivity)
511 class(brinkman_design_t), intent(inout) :: this
512 type(vector_t), intent(in) :: sensitivity
513 type(field_t), pointer :: tmp_fld
514 integer :: temp_indices(1)
515
516 call neko_scratch_registry%request_field(tmp_fld, temp_indices(1))
517
518 call vector_to_field(tmp_fld, sensitivity)
519
520 call this%mapping%apply_backward(this%sensitivity, tmp_fld)
521
522 ! TODO
523 ! DELETE THIS LATER
524 !
525 ! When Abbas writes the interface for the optimization
526 ! module this may be a moot point, because we would only really collect
527 ! the sensitivity of the design variables inside the mask.
528 !
529 ! Note for Abbas,
530 ! I'm NOT doing this because I'm too lazy and I just need masks so I can
531 ! test something in the passive scalar.
532 if (this%if_mask) then
533 call mask_exterior_const(this%sensitivity, this%optimization_domain, &
534 0.0_rp)
535 end if
536
537 call neko_scratch_registry%relinquish_field(temp_indices)
538
539 end subroutine brinkman_design_map_backward
540
541 subroutine brinkman_design_write(this, idx)
542 class(brinkman_design_t), intent(inout) :: this
543 integer, intent(in) :: idx
544
545 call this%output%sample(real(idx, kind=rp))
546
547 end subroutine brinkman_design_write
548
549end module brinkman_design
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:222
subroutine, public vector_to_field(field, vector)
Vector to field.
Definition neko_ext.f90:198
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.