37 use num_types,
only: rp, sp
38 use json_module,
only: json_file
42 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
48 use device,
only: device_memcpy, host_to_device, device_to_host
49 use device_math,
only: device_copy
53 use json_module,
only: json_file
54 use json_utils,
only: json_get
55 use utils,
only: neko_error
57 use vector,
only: vector_t
60 use mpi_f08,
only: mpi_exscan, mpi_sum, mpi_integer, mpi_allreduce
61 use comm,
only: pe_rank, pe_size, neko_comm, mpi_real_precision
63 use fld_file_output,
only: fld_file_output_t
72 type(vector_t) :: values
80 generic,
public :: init => init_from_components
82 procedure, pass(this) :: init_from_components => &
83 design_3dto1d_init_from_components
86 procedure, pass(this) :: get_values => design_3dto1d_get_values
89 procedure, pass(this) :: update_design => design_3dto1d_update_design
92 procedure, pass(this) :: write => design_3dto1d_write
95 procedure, pass(this) :: free => design_3dto1d_free
98 procedure, pass(this) :: map_forward => design_3dto1d_map_forward
100 procedure, pass(this) :: map_backward => design_3dto1d_map_backward
107 subroutine design_3dto1d_init_from_components(this, n)
109 integer,
intent(in) :: n
111 call this%init_base(
'design_3dto1d', n)
113 call this%values%init(n)
115 end subroutine design_3dto1d_init_from_components
118 subroutine design_3dto1d_free(this)
121 call this%free_base()
122 call this%values%free()
123 end subroutine design_3dto1d_free
125 subroutine design_3dto1d_map_forward(this)
128 end subroutine design_3dto1d_map_forward
130 subroutine design_3dto1d_map_backward(this, sensitivity)
132 type(vector_t),
intent(in) :: sensitivity
133 end subroutine design_3dto1d_map_backward
136 subroutine design_3dto1d_get_values(this, values)
138 type(vector_t),
intent(inout) :: values
140 if (this%size() .ne. values%size())
then
141 call neko_error(
'Get design: size mismatch')
146 end subroutine design_3dto1d_get_values
148 subroutine design_3dto1d_update_design(this, values)
150 type(vector_t),
intent(inout) :: values
155 if (neko_bcknd_device .eq. 1)
then
156 call device_copy(this%values%x_d, values%x_d, n)
162 end subroutine design_3dto1d_update_design
164 subroutine design_3dto1d_write(this, idx)
166 integer,
intent(in) :: idx
168 character(len=100) :: filename
169 integer :: i, iunit, ierr
171 real(rp),
allocatable :: global_values(:)
172 real(rp),
allocatable :: global_x(:)
173 integer :: global_size, local_size
175 integer,
allocatable :: recvcounts(:), displs(:)
176 integer :: root_rank = 0
181 local_size = this%size()
184 if (pe_rank == root_rank)
then
186 global_size = this%size_global()
187 allocate(global_values(global_size))
188 allocate(global_x(global_size))
189 allocate(recvcounts(pe_size), displs(pe_size))
191 le = l_total / real(global_size, kind=rp)
192 do i = 1, global_size
193 global_x(i) = le * (i - 0.5_rp)
197 allocate(global_values(1), recvcounts(1), displs(1))
201 call mpi_gather(local_size, 1, mpi_integer, &
202 recvcounts, 1, mpi_integer, &
203 root_rank, neko_comm, ierr)
206 if (pe_rank == root_rank)
then
209 displs(i) = displs(i-1) + recvcounts(i-1)
213 if (neko_bcknd_device .eq. 1)
then
214 call device_memcpy(this%values%x, this%values%x_d, local_size, &
215 device_to_host, sync = .false.)
218 call mpi_gatherv(this%values%x, local_size, mpi_real_precision, &
219 global_values, recvcounts, displs, mpi_real_precision, &
220 root_rank, neko_comm, ierr)
222 if (pe_rank == root_rank)
then
224 write(filename,
'(A,I0.6,A)')
'design_iter_', idx,
'.txt'
226 open(newunit=iunit, file=trim(filename), status=
'replace', action=
'write')
228 write(iunit,
'(A,I0)')
'# Iteration: ', idx
229 write(iunit,
'(A)')
'# x_position height'
232 do i = 1, global_size
233 write(iunit,
'(2E16.8)') global_x(i), global_values(i)
238 deallocate(global_values, global_x, recvcounts, displs)
239 print *,
"Design written to ", trim(filename)
241 deallocate(global_values, recvcounts, displs)
243 end subroutine design_3dto1d_write
245end module design_3dto1d
Mappings to be applied to a scalar field.
Some common Masking operations we may need.
A RAMP mapping of coefficients.
Implements the steady_problem_t type.
A topology optimization design variable.
Base abstract class for mapping.
A PDE based filter mapping , see Lazarov & O. Sigmund 2010, by solving an equation of the form .
A RAMP mapping of coefficients This is the standard RAMP described in https://doi....