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
50 use json_module,
only: json_file
52 use vector,
only: vector_t
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
72 type(field_t),
pointer :: design_indicator
82 type(field_t),
pointer :: brinkman_amplitude
114 type(field_t),
pointer :: sensitivity
170 class(point_zone_t),
pointer :: optimization_domain
179 type(fld_file_output_t),
private :: output
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
196 procedure, pass(this) :: get_values => brinkman_design_get_design
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
206 procedure, pass(this) :: update_design => brinkman_design_update_design
211 procedure, pass(this) :: map_forward => brinkman_design_map_forward
215 procedure, pass(this) :: map_backward => brinkman_design_map_backward
223 procedure, pass(this) :: write => brinkman_design_write
226 procedure, pass(this) :: free => brinkman_design_free
237 subroutine brinkman_design_init_from_json_sim(this, parameters, simulation)
239 type(json_file),
intent(inout) :: parameters
241 type(json_file) :: json_subdict
242 character(len=:),
allocatable :: domain_name, domain_type
245 if (parameters%valid_path(
'optimization.domain'))
then
246 call json_get(parameters,
'optimization.domain.type', domain_type)
247 select case (trim(domain_type))
249 this%if_mask = .true.
250 call json_get(parameters,
'optimization.domain.zone_name', &
252 this%optimization_domain => &
253 neko_point_zone_registry%get_point_zone(domain_name)
256 call neko_error(
'brinkman design only supports point_zones for&
257 & optimization domain types')
261 this%if_mask = .false.
265 call this%init_from_components(simulation)
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')
273 if (parameters%valid_path(&
274 'optimization.design.initial_distribution'))
then
275 call json_extract_object(parameters, &
276 'optimization.design.initial_distribution', json_subdict)
280 call field_rzero(this%design_indicator)
285 call this%map_forward()
287 end subroutine brinkman_design_init_from_json_sim
290 subroutine brinkman_design_free(this)
293 call this%free_base()
294 call this%brinkman_amplitude%free()
295 call this%design_indicator%free()
296 call this%sensitivity%free()
298 end subroutine brinkman_design_free
300 subroutine brinkman_design_init_from_components(this, simulation)
302 type(simulation_t),
intent(inout) :: simulation
304 type(simple_brinkman_source_term_t) :: forward_brinkman, adjoint_brinkman
306 associate(dof => simulation%neko_case%fluid%dm_Xh)
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.)
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")
325 this%design_indicator = 0.0_rp
326 this%brinkman_amplitude = 0.0_rp
327 this%design_indicator%x = 0.0_rp
329 n = this%design_indicator%dof%size()
332 this%design_indicator%x(i,1,1,1) = 0.0_rp
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.)
370 if (this%if_mask)
then
371 call mask_exterior_const(this%design_indicator, &
372 this%optimization_domain, 0.0_rp)
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)
389 call this%init_base(n)
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)
402 call simulation%fluid%source_term%add(forward_brinkman)
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)
415 call simulation%adjoint_fluid%source_term%add(adjoint_brinkman)
417 end subroutine brinkman_design_init_from_components
420 subroutine brinkman_design_map_forward(this)
424 if (this%if_mask)
then
425 call mask_exterior_const(this%design_indicator, &
426 this%optimization_domain, 0.0_rp)
429 call this%mapping%apply_forward(this%brinkman_amplitude, &
430 this%design_indicator)
432 end subroutine brinkman_design_map_forward
434 function brinkman_design_get_design(this)
result(values)
436 type(vector_t) :: values
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)
446 end function brinkman_design_get_design
448 function brinkman_design_get_x(this)
result(x)
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)
460 end function brinkman_design_get_x
462 function brinkman_design_get_y(this)
result(y)
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)
474 end function brinkman_design_get_y
476 function brinkman_design_get_z(this)
result(z)
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)
488 end function brinkman_design_get_z
490 subroutine brinkman_design_update_design(this, values)
492 type(vector_t),
intent(inout) :: values
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)
501 call this%map_forward()
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)
508 end subroutine brinkman_design_update_design
510 subroutine brinkman_design_map_backward(this, sensitivity)
512 type(vector_t),
intent(in) :: sensitivity
513 type(field_t),
pointer :: tmp_fld
514 integer :: temp_indices(1)
516 call neko_scratch_registry%request_field(tmp_fld, temp_indices(1))
518 call vector_to_field(tmp_fld, sensitivity)
520 call this%mapping%apply_backward(this%sensitivity, tmp_fld)
532 if (this%if_mask)
then
533 call mask_exterior_const(this%sensitivity, this%optimization_domain, &
537 call neko_scratch_registry%relinquish_field(temp_indices)
539 end subroutine brinkman_design_map_backward
541 subroutine brinkman_design_write(this, idx)
543 integer,
intent(in) :: idx
545 call this%output%sample(real(idx, kind=rp))
547 end subroutine brinkman_design_write
549end module brinkman_design
Implements the mapping_handler_t type.
Mappings to be applied to a scalar field.
Some common Masking operations we may need.
Contains extensions to the neko library required to run the topology optimization code.
subroutine, public field_to_vector(vector, field)
Field to vector.
subroutine, public vector_to_field(field, vector)
Vector to field.
Optimization initial condition.
Implements the simple_brinkman_source_term_t type.
Implements the steady_problem_t type.
A topology optimization design variable.
Abstract class for handling mapping_cascade.
A simple Brinkman source term.