35 use num_types,
only: rp, sp
36 use json_module,
only: json_file
40 use coefs,
only: coef_t
41 use scratch_registry,
only: neko_scratch_registry
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, device_to_host
47 use device_math,
only: device_copy
51 use json_module,
only: json_file
52 use json_utils,
only: json_get
54 use vector,
only: vector_t
57 use mpi_f08,
only: mpi_exscan, mpi_sum, mpi_integer, mpi_allreduce
58 use comm,
only: pe_rank, pe_size, neko_comm, mpi_real_precision
60 use fld_file_output,
only: fld_file_output_t
69 type(vector_t) :: values
77 generic,
public :: init => init_from_components
79 procedure, pass(this) :: init_from_components => &
80 design_3dto1d_init_from_components
83 procedure, pass(this) :: get_values => design_3dto1d_get_values
86 procedure, pass(this) :: update_design => design_3dto1d_update_design
89 procedure, pass(this) :: write => design_3dto1d_write
92 procedure, pass(this) :: free => design_3dto1d_free
96 procedure, pass(this) :: map_forward => design_3dto1d_map_forward
98 procedure, pass(this) :: map_backward => design_3dto1d_map_backward
105 subroutine design_3dto1d_init_from_components(this, n)
107 integer,
intent(in) :: n
109 call this%init_base(
'design_3dto1d', n)
111 call this%values%init(n)
113 end subroutine design_3dto1d_init_from_components
116 subroutine design_3dto1d_free(this)
119 call this%free_base()
120 call this%values%free()
121 end subroutine design_3dto1d_free
123 subroutine design_3dto1d_map_forward(this)
126 end subroutine design_3dto1d_map_forward
128 subroutine design_3dto1d_map_backward(this, sensitivity)
130 type(vector_t),
intent(in) :: sensitivity
131 end subroutine design_3dto1d_map_backward
134 subroutine design_3dto1d_get_values(this, values)
136 type(vector_t),
intent(inout) :: values
140 end subroutine design_3dto1d_get_values
142 subroutine design_3dto1d_update_design(this, values)
144 type(vector_t),
intent(inout) :: values
149 if (neko_bcknd_device .eq. 1)
then
150 call device_copy(this%values%x_d, values%x_d, n)
156 end subroutine design_3dto1d_update_design
158 subroutine design_3dto1d_write(this, idx)
160 integer,
intent(in) :: idx
162 character(len=100) :: filename
163 integer :: i, iunit, ierr
165 real(rp),
allocatable :: global_values(:)
166 real(rp),
allocatable :: global_x(:)
167 integer :: global_size, local_size
169 integer,
allocatable :: recvcounts(:), displs(:)
170 integer :: root_rank = 0
175 local_size = this%size()
178 if (pe_rank == root_rank)
then
180 global_size = this%size_global()
181 allocate(global_values(global_size))
182 allocate(global_x(global_size))
183 allocate(recvcounts(pe_size), displs(pe_size))
185 le = l_total / real(global_size, kind=rp)
186 do i = 1, global_size
187 global_x(i) = le * (i - 0.5_rp)
191 allocate(global_values(1), recvcounts(1), displs(1))
195 call mpi_gather(local_size, 1, mpi_integer, &
196 recvcounts, 1, mpi_integer, &
197 root_rank, neko_comm, ierr)
200 if (pe_rank == root_rank)
then
203 displs(i) = displs(i-1) + recvcounts(i-1)
207 if (neko_bcknd_device .eq. 1)
then
208 call device_memcpy(this%values%x, this%values%x_d, local_size, &
209 device_to_host, sync = .false.)
212 call mpi_gatherv(this%values%x, local_size, mpi_real_precision, &
213 global_values, recvcounts, displs, mpi_real_precision, &
214 root_rank, neko_comm, ierr)
216 if (pe_rank == root_rank)
then
218 write(filename,
'(A,I0.6,A)')
'design_iter_', idx,
'.txt'
220 open(newunit=iunit, file=trim(filename), status=
'replace', action=
'write')
222 write(iunit,
'(A,I0)')
'# Iteration: ', idx
223 write(iunit,
'(A)')
'# x_position height'
226 do i = 1, global_size
227 write(iunit,
'(2E16.8)') global_x(i), global_values(i)
232 deallocate(global_values, global_x, recvcounts, displs)
233 print *,
"Design written to ", trim(filename)
235 deallocate(global_values, recvcounts, displs)
237 end subroutine design_3dto1d_write
239end 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....