37 use num_types,
only: rp, sp, dp
38 use field,
only: field_t
39 use json_module,
only: json_file
41 use coefs,
only: coef_t
43 use scratch_registry,
only: neko_scratch_registry
44 use point_zone_registry,
only: neko_point_zone_registry
45 use point_zone,
only: point_zone_t
47 use neko_config,
only: neko_bcknd_device
51 use vector,
only: vector_t
52 use matrix,
only: matrix_t
53 use vector_math,
only: vector_cmult
54 use matrix_math,
only: matrix_cmult
55 use math,
only: copy, col2
56 use device_math,
only: device_copy, device_col2
57 use registry,
only: neko_registry
60 use field_math,
only: field_rzero
61 use json_utils,
only: json_get, json_get_or_default
62 use utils,
only: neko_error
74 type(field_t),
pointer :: design_indicator
84 type(field_t),
pointer :: brinkman_amplitude
116 type(field_t),
pointer :: sensitivity
172 class(point_zone_t),
pointer :: optimization_domain
176 class(coef_t),
pointer :: coef
178 real(kind=rp) :: avg_b
190 generic,
public :: init => init_from_json_sim, init_from_components
192 procedure, pass(this),
public :: init_from_json_sim => &
193 brinkman_design_init_from_json_sim
195 procedure, pass(this),
public :: init_from_components => &
196 brinkman_design_init_from_components
199 procedure, pass(this) :: get_values => brinkman_design_get_design
201 procedure, pass(this) :: get_sensitivity => brinkman_design_get_sensitivity
204 procedure, pass(this) :: design_get_x => brinkman_design_get_x
206 procedure, pass(this) :: design_get_x_i => brinkman_design_get_x_i
208 procedure, pass(this) :: design_get_y => brinkman_design_get_y
210 procedure, pass(this) :: design_get_y_i => brinkman_design_get_y_i
212 procedure, pass(this) :: design_get_z => brinkman_design_get_z
214 procedure, pass(this) :: design_get_z_i => brinkman_design_get_z_i
217 procedure, pass(this) :: update_design => brinkman_design_update_design
222 procedure, pass(this) :: map_forward => brinkman_design_map_forward
226 procedure, pass(this) :: map_backward => brinkman_design_map_backward
229 procedure, pass(this) :: convert_to_directional_derivative => &
230 brinkman_design_convert_to_directional_derivative
232 procedure, pass(this) :: project_sensitivity_vector => &
233 brinkman_design_project_sensitivity_vector
235 procedure, pass(this) :: project_sensitivity_matrix => &
236 brinkman_design_project_sensitivity_matrix
237 generic,
public :: project_sensitivity => &
238 project_sensitivity_vector, project_sensitivity_matrix
246 procedure, pass(this) :: write => brinkman_design_write
248 procedure, pass(this) :: set_output_counter => &
249 brinkman_design_set_output_counter
252 procedure, pass(this) :: free => brinkman_design_free
263 subroutine brinkman_design_init_from_json_sim(this, parameters, simulation)
265 type(json_file),
intent(inout) :: parameters
267 type(json_file) :: json_subdict
268 character(len=:),
allocatable :: domain_name, domain_type, name
269 character(len=:),
allocatable :: output_format_str, output_precision_str
270 logical :: dealias, verbose_design, verbose_sensitivity
271 integer :: output_precision
273 call json_get_or_default(parameters,
'name', name,
'Brinkman Design')
274 call json_get_or_default(parameters,
'domain.type', domain_type,
'full')
275 call json_get_or_default(parameters,
'dealias', dealias, .true.)
276 call json_get_or_default(parameters,
'verbose_design', verbose_design, &
278 call json_get_or_default(parameters,
'verbose_sensitivity', &
279 verbose_sensitivity, .false.)
280 call json_get_or_default(parameters,
'output_precision', &
281 output_precision_str,
'sp')
282 call json_get_or_default(parameters,
'output_format', &
283 output_format_str,
'fld')
285 select case (trim(output_precision_str))
286 case (
'sp',
'SP',
'single')
287 output_precision = sp
288 case (
'dp',
'DP',
'double')
289 output_precision = dp
291 call neko_error(
'Invalid output_precision in design.brinkman. ' // &
292 'Expected "sp" or "dp".')
295 select case (trim(domain_type))
297 this%has_mask = .false.
299 this%has_mask = .true.
300 call json_get(parameters,
'domain.zone_name', domain_name)
301 this%optimization_domain => &
302 neko_point_zone_registry%get_point_zone(domain_name)
305 call neko_error(
'brinkman design only supports point_zones for ' // &
306 'optimization domain types')
311 call this%init_from_components(name, simulation, dealias)
314 associate(coef => simulation%neko_case%fluid%c_Xh, &
315 gs => simulation%neko_case%fluid%gs_Xh)
317 if (
'mapping' .in. parameters)
then
318 call this%mapping%init_base(coef)
319 call this%mapping%add(parameters,
'mapping')
320 call this%mapping%init_output_fields(this%design_indicator, &
321 this%brinkman_amplitude, this%sensitivity, verbose_design, &
322 verbose_sensitivity, output_precision, output_format_str, &
323 'design',
'sensitivity')
326 if (
'initial_distribution' .in. parameters)
then
327 call json_get(parameters,
'initial_distribution', json_subdict)
331 call field_rzero(this%design_indicator)
336 call this%map_forward()
338 end subroutine brinkman_design_init_from_json_sim
341 subroutine brinkman_design_free(this)
344 call this%free_base()
345 nullify(this%brinkman_amplitude)
346 nullify(this%design_indicator)
347 nullify(this%sensitivity)
349 end subroutine brinkman_design_free
351 subroutine brinkman_design_init_from_components(this, name, simulation, &
354 character(len=*),
intent(in) :: name
355 type(simulation_t),
intent(inout) :: simulation
356 logical,
intent(in) :: dealias
358 type(simple_brinkman_source_term_t) :: forward_brinkman, adjoint_brinkman
360 associate(dof => simulation%neko_case%fluid%dm_Xh)
362 call neko_registry%add_field(dof,
"design_indicator", .true.)
363 call neko_registry%add_field(dof,
"brinkman_amplitude", .true.)
364 call neko_registry%add_field(dof,
"sensitivity", .true.)
368 this%design_indicator => neko_registry%get_field(
"design_indicator")
369 this%brinkman_amplitude => neko_registry%get_field(
"brinkman_amplitude")
370 this%sensitivity => neko_registry%get_field(
"sensitivity")
372 this%coef => simulation%fluid%c_Xh
378 this%design_indicator = 0.0_rp
379 this%brinkman_amplitude = 0.0_rp
380 this%design_indicator = 0.0_rp
410 if (this%has_mask)
then
411 call mask_exterior_const(this%design_indicator, &
412 this%optimization_domain, 0.0_rp)
415 n = this%design_indicator%dof%size()
416 call this%init_base(name, n)
419 this%avg_B = this%coef%volume / real(simulation%fluid%glb_unique_points)
422 call forward_brinkman%init_from_components( &
423 simulation%fluid%f_x, &
424 simulation%fluid%f_y, &
425 simulation%fluid%f_z, &
426 this%brinkman_amplitude, &
427 simulation%fluid%u, &
428 simulation%fluid%v, &
429 simulation%fluid%w, &
430 simulation%fluid%c_Xh, &
431 simulation%adjoint_fluid%c_Xh_GL, &
432 simulation%adjoint_fluid%GLL_to_GL, &
433 dealias, simulation%adjoint_fluid%scratch_GL)
435 call simulation%fluid%source_term%add(forward_brinkman)
438 call adjoint_brinkman%init_from_components( &
439 simulation%adjoint_fluid%f_adj_x, &
440 simulation%adjoint_fluid%f_adj_y, &
441 simulation%adjoint_fluid%f_adj_z, &
442 this%brinkman_amplitude, &
443 simulation%adjoint_fluid%u_adj, &
444 simulation%adjoint_fluid%v_adj, &
445 simulation%adjoint_fluid%w_adj, &
446 simulation%adjoint_fluid%c_Xh, &
447 simulation%adjoint_fluid%c_Xh_GL, &
448 simulation%adjoint_fluid%GLL_to_GL, &
449 dealias, simulation%adjoint_fluid%scratch_GL)
452 select type (f => simulation%adjoint_fluid)
453 type is (adjoint_fluid_pnpn_t)
454 call f%source_term%add(adjoint_brinkman)
458 end subroutine brinkman_design_init_from_components
461 subroutine brinkman_design_map_forward(this)
465 if (this%has_mask)
then
466 call mask_exterior_const(this%design_indicator, &
467 this%optimization_domain, 0.0_rp)
470 call this%mapping%apply_forward(this%brinkman_amplitude, &
471 this%design_indicator)
473 end subroutine brinkman_design_map_forward
475 subroutine brinkman_design_get_design(this, values)
477 type(vector_t),
intent(inout) :: values
481 if (n .ne. values%size())
then
482 call neko_error(
'Get design: size mismatch')
485 if (neko_bcknd_device .eq. 1)
then
486 call device_copy(values%x_d, this%design_indicator%x_d, n)
488 call copy(values%x, this%design_indicator%x, n)
491 end subroutine brinkman_design_get_design
493 subroutine brinkman_design_get_sensitivity(this, values)
495 type(vector_t),
intent(inout) :: values
499 if (n .ne. values%size())
then
500 call neko_error(
'Get sensitivity: size mismatch')
503 if (neko_bcknd_device .eq. 1)
then
504 call device_copy(values%x_d, this%sensitivity%x_d, n)
506 call copy(values%x, this%sensitivity%x, n)
509 end subroutine brinkman_design_get_sensitivity
515 subroutine brinkman_design_convert_to_directional_derivative(this, vec_in)
517 type(vector_t),
intent(inout) :: vec_in
521 if (n .ne. vec_in%size())
then
522 call neko_error(
'convert_to_directional_derivative: size mismatch')
525 if (neko_bcknd_device .eq. 1)
then
526 call device_col2(vec_in%x_d, this%coef%B_d, n)
528 call col2(vec_in%x, this%coef%B, n)
531 end subroutine brinkman_design_convert_to_directional_derivative
536 subroutine brinkman_design_project_sensitivity_vector(this, sensitivity)
538 type(vector_t),
intent(inout) :: sensitivity
540 call vector_cmult(sensitivity, this%avg_B)
541 end subroutine brinkman_design_project_sensitivity_vector
543 subroutine brinkman_design_project_sensitivity_matrix(this, sensitivity)
545 type(matrix_t),
intent(inout) :: sensitivity
547 call matrix_cmult(sensitivity, this%avg_B)
548 end subroutine brinkman_design_project_sensitivity_matrix
550 subroutine brinkman_design_get_x(this, x)
552 type(vector_t),
intent(inout) :: x
556 if (n .ne. x%size())
then
557 call neko_error(
'Get x: size mismatch')
560 if (neko_bcknd_device .eq. 1)
then
561 call device_copy(x%x_d, this%design_indicator%dof%x_d, n)
563 call copy(x%x, this%design_indicator%dof%x, n)
566 end subroutine brinkman_design_get_x
568 function brinkman_design_get_x_i(this, i)
result(x_i)
570 integer,
intent(in) :: i
575 if (i .lt. 1 .or. i .gt. n)
then
576 call neko_error(
'brinkman_design_get_x_i: index out of bounds')
579 x_i = this%design_indicator%dof%x(i,1,1,1)
581 end function brinkman_design_get_x_i
583 subroutine brinkman_design_get_y(this, y)
585 type(vector_t),
intent(inout) :: y
589 if (n .ne. y%size())
then
590 call neko_error(
'Get y: size mismatch')
593 if (neko_bcknd_device .eq. 1)
then
594 call device_copy(y%x_d, this%design_indicator%dof%y_d, n)
596 call copy(y%x, this%design_indicator%dof%y, n)
599 end subroutine brinkman_design_get_y
601 function brinkman_design_get_y_i(this, i)
result(y_i)
603 integer,
intent(in) :: i
608 if (i .lt. 1 .or. i .gt. n)
then
609 call neko_error(
'brinkman_design_get_y_i: index out of bounds')
612 y_i = this%design_indicator%dof%y(i,1,1,1)
614 end function brinkman_design_get_y_i
616 subroutine brinkman_design_get_z(this, z)
618 type(vector_t),
intent(inout) :: z
622 if (n .ne. z%size())
then
623 call neko_error(
'Get z: size mismatch')
626 if (neko_bcknd_device .eq. 1)
then
627 call device_copy(z%x_d, this%design_indicator%dof%z_d, n)
629 call copy(z%x, this%design_indicator%dof%z, n)
632 end subroutine brinkman_design_get_z
634 function brinkman_design_get_z_i(this, i)
result(z_i)
636 integer,
intent(in) :: i
641 if (i .lt. 1 .or. i .gt. n)
then
642 call neko_error(
'brinkman_design_get_z_i: index out of bounds')
645 z_i = this%design_indicator%dof%z(i,1,1,1)
647 end function brinkman_design_get_z_i
649 subroutine brinkman_design_update_design(this, values)
651 type(vector_t),
intent(inout) :: values
655 if (neko_bcknd_device .eq. 1)
then
656 call device_copy(this%design_indicator%x_d, values%x_d, n)
658 call copy(this%design_indicator%x, values%x, n)
661 call this%map_forward()
663 if (neko_bcknd_device .eq. 1)
then
664 call device_copy(values%x_d, this%design_indicator%x_d, n)
666 call copy(values%x, this%design_indicator%x, n)
669 end subroutine brinkman_design_update_design
671 subroutine brinkman_design_map_backward(this, sensitivity)
673 type(vector_t),
intent(in) :: sensitivity
674 type(field_t),
pointer :: tmp_fld
675 integer :: temp_index
677 call neko_scratch_registry%request(tmp_fld, temp_index, .false.)
679 call vector_to_field(tmp_fld, sensitivity)
681 call this%mapping%apply_backward(this%sensitivity, tmp_fld)
683 if (this%has_mask)
then
684 call mask_exterior_const(this%sensitivity, this%optimization_domain, &
688 call neko_scratch_registry%relinquish_field(temp_index)
690 end subroutine brinkman_design_map_backward
692 subroutine brinkman_design_write(this, idx)
694 integer,
intent(in) :: idx
696 call this%mapping%write_design(idx)
697 call this%mapping%write_sensitivity(idx)
699 end subroutine brinkman_design_write
701 subroutine brinkman_design_set_output_counter(this, idx)
703 integer,
intent(in) :: idx
705 call this%mapping%set_output_counter(idx)
707 end subroutine brinkman_design_set_output_counter
709end 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.