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
51 use json_module,
only: json_file
53 use vector,
only: vector_t
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
73 type(field_t),
pointer :: design_indicator
83 type(field_t),
pointer :: brinkman_amplitude
115 type(field_t),
pointer :: sensitivity
171 class(point_zone_t),
pointer :: optimization_domain
180 type(fld_file_output_t),
private :: output
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
197 procedure, pass(this) :: get_values => brinkman_design_get_design
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
207 procedure, pass(this) :: update_design => brinkman_design_update_design
212 procedure, pass(this) :: map_forward => brinkman_design_map_forward
216 procedure, pass(this) :: map_backward => brinkman_design_map_backward
224 procedure, pass(this) :: write => brinkman_design_write
227 procedure, pass(this) :: free => brinkman_design_free
238 subroutine brinkman_design_init_from_json_sim(this, parameters, simulation)
240 type(json_file),
intent(inout) :: parameters
242 type(json_file) :: json_subdict
243 character(len=:),
allocatable :: domain_name, domain_type
246 if (parameters%valid_path(
'optimization.domain'))
then
247 call json_get(parameters,
'optimization.domain.type', domain_type)
248 select case (trim(domain_type))
250 this%if_mask = .true.
251 call json_get(parameters,
'optimization.domain.zone_name', &
253 this%optimization_domain => &
254 neko_point_zone_registry%get_point_zone(domain_name)
257 call neko_error(
'brinkman design only supports point_zones for&
258 & optimization domain types')
262 this%if_mask = .false.
266 call this%init_from_components(simulation)
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')
274 if (parameters%valid_path(&
275 'optimization.design.initial_distribution'))
then
276 call json_extract_object(parameters, &
277 'optimization.design.initial_distribution', json_subdict)
281 call field_rzero(this%design_indicator)
286 call this%map_forward()
288 end subroutine brinkman_design_init_from_json_sim
291 subroutine brinkman_design_free(this)
294 call this%free_base()
295 call this%brinkman_amplitude%free()
296 call this%design_indicator%free()
297 call this%sensitivity%free()
299 end subroutine brinkman_design_free
301 subroutine brinkman_design_init_from_components(this, simulation)
303 type(simulation_t),
intent(inout) :: simulation
305 type(simple_brinkman_source_term_t) :: forward_brinkman, adjoint_brinkman
307 associate(dof => simulation%neko_case%fluid%dm_Xh)
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.)
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")
326 this%design_indicator = 0.0_rp
327 this%brinkman_amplitude = 0.0_rp
328 this%design_indicator%x = 0.0_rp
330 n = this%design_indicator%dof%size()
333 this%design_indicator%x(i,1,1,1) = 0.0_rp
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.)
371 if (this%if_mask)
then
372 call mask_exterior_const(this%design_indicator, &
373 this%optimization_domain, 0.0_rp)
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)
390 call this%init_base(n)
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)
403 call simulation%fluid%source_term%add(forward_brinkman)
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)
417 select type (f => simulation%adjoint_fluid)
418 type is (adjoint_fluid_pnpn_t)
419 call f%source_term%add(adjoint_brinkman)
423 end subroutine brinkman_design_init_from_components
426 subroutine brinkman_design_map_forward(this)
430 if (this%if_mask)
then
431 call mask_exterior_const(this%design_indicator, &
432 this%optimization_domain, 0.0_rp)
435 call this%mapping%apply_forward(this%brinkman_amplitude, &
436 this%design_indicator)
438 end subroutine brinkman_design_map_forward
440 function brinkman_design_get_design(this)
result(values)
442 type(vector_t) :: values
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)
452 end function brinkman_design_get_design
454 function brinkman_design_get_x(this)
result(x)
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)
466 end function brinkman_design_get_x
468 function brinkman_design_get_y(this)
result(y)
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)
480 end function brinkman_design_get_y
482 function brinkman_design_get_z(this)
result(z)
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)
494 end function brinkman_design_get_z
496 subroutine brinkman_design_update_design(this, values)
498 type(vector_t),
intent(inout) :: values
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)
507 call this%map_forward()
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)
514 end subroutine brinkman_design_update_design
516 subroutine brinkman_design_map_backward(this, sensitivity)
518 type(vector_t),
intent(in) :: sensitivity
519 type(field_t),
pointer :: tmp_fld
520 integer :: temp_indices(1)
522 call neko_scratch_registry%request_field(tmp_fld, temp_indices(1))
524 call vector_to_field(tmp_fld, sensitivity)
526 call this%mapping%apply_backward(this%sensitivity, tmp_fld)
538 if (this%if_mask)
then
539 call mask_exterior_const(this%sensitivity, this%optimization_domain, &
543 call neko_scratch_registry%relinquish_field(temp_indices)
545 end subroutine brinkman_design_map_backward
547 subroutine brinkman_design_write(this, idx)
549 integer,
intent(in) :: idx
551 call this%output%sample(real(idx, kind=rp))
553 end subroutine brinkman_design_write
555end module brinkman_design
Adjoint Pn/Pn formulation.
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.