Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
design_3dto1d.f90
Go to the documentation of this file.
1
34
35! Implements the `design_3dto1d_t` type.
36module design_3dto1d
37 use num_types, only: rp, sp
38 use json_module, only: json_file
39 use mapping, only: mapping_t
40 use pde_filter, only: pde_filter_t
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
50 use design, only: design_t
51 use math, only: rzero
52 use simulation_m, only: simulation_t
53 use json_module, only: json_file
54 use json_utils, only: json_get
55 use utils, only: neko_error
56
57 use vector, only: vector_t
58 use math, only: copy
59
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
62
63 use fld_file_output, only: fld_file_output_t
64
65 implicit none
66 private
67
69 type, extends(design_t), public :: design_3dto1d_t
70 private
71
72 type(vector_t) :: values
73
74 contains
75
76 ! ----------------------------------------------------------------------- !
77 ! Initializations
78
80 generic, public :: init => init_from_components
82 procedure, pass(this) :: init_from_components => &
83 design_3dto1d_init_from_components
84
86 procedure, pass(this) :: get_values => design_3dto1d_get_values
87
89 procedure, pass(this) :: update_design => design_3dto1d_update_design
90
92 procedure, pass(this) :: write => design_3dto1d_write
93
95 procedure, pass(this) :: free => design_3dto1d_free
96
98 procedure, pass(this) :: map_forward => design_3dto1d_map_forward
100 procedure, pass(this) :: map_backward => design_3dto1d_map_backward
101
102 end type design_3dto1d_t
103
104contains
105
106
107 subroutine design_3dto1d_init_from_components(this, n)
108 class(design_3dto1d_t), intent(inout) :: this
109 integer, intent(in) :: n
110
111 call this%init_base('design_3dto1d', n)
112
113 call this%values%init(n)
114
115 end subroutine design_3dto1d_init_from_components
116
118 subroutine design_3dto1d_free(this)
119 class(design_3dto1d_t), intent(inout) :: this
120
121 call this%free_base()
122 call this%values%free()
123 end subroutine design_3dto1d_free
124
125 subroutine design_3dto1d_map_forward(this)
126 class(design_3dto1d_t), intent(inout) :: this
127
128 end subroutine design_3dto1d_map_forward
129
130 subroutine design_3dto1d_map_backward(this, sensitivity)
131 class(design_3dto1d_t), intent(inout) :: this
132 type(vector_t), intent(in) :: sensitivity
133 end subroutine design_3dto1d_map_backward
134
135
136 subroutine design_3dto1d_get_values(this, values)
137 class(design_3dto1d_t), intent(in) :: this
138 type(vector_t), intent(inout) :: values
139
140 if (this%size() .ne. values%size()) then
141 call neko_error('Get design: size mismatch')
142 end if
143
144 values = this%values
145
146 end subroutine design_3dto1d_get_values
147
148 subroutine design_3dto1d_update_design(this, values)
149 class(design_3dto1d_t), intent(inout) :: this
150 type(vector_t), intent(inout) :: values
151 integer :: n
152
153 n = this%size()
154
155 if (neko_bcknd_device .eq. 1) then
156 call device_copy(this%values%x_d, values%x_d, n)
157 else
158 this%values = values
159 end if
160
161
162 end subroutine design_3dto1d_update_design
163
164 subroutine design_3dto1d_write(this, idx)
165 class(design_3dto1d_t), intent(inout) :: this
166 integer, intent(in) :: idx
167
168 character(len=100) :: filename
169 integer :: i, iunit, ierr
170 real(rp) :: le
171 real(rp), allocatable :: global_values(:)
172 real(rp), allocatable :: global_x(:)
173 integer :: global_size, local_size
174 real(rp) :: l_total
175 integer, allocatable :: recvcounts(:), displs(:)
176 integer :: root_rank = 0
177
178 l_total = 2.0_rp
179
180 ! Get local size on all ranks
181 local_size = this%size()
182
183 ! Only rank 0 handles global arrays
184 if (pe_rank == root_rank) then
185 ! Get global size
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))
190 ! Calculate element length and x positions
191 le = l_total / real(global_size, kind=rp)
192 do i = 1, global_size
193 global_x(i) = le * (i - 0.5_rp) ! Center of each element
194 end do
195 else
196 ! Non-root ranks: minimal dummy allocations
197 allocate(global_values(1), recvcounts(1), displs(1))
198 endif
199
200 ! First, gather all the local sizes to rank 0
201 call mpi_gather(local_size, 1, mpi_integer, &
202 recvcounts, 1, mpi_integer, &
203 root_rank, neko_comm, ierr)
204
205 ! Calculate displacements on rank 0
206 if (pe_rank == root_rank) then
207 displs(1) = 0
208 do i = 2, pe_size
209 displs(i) = displs(i-1) + recvcounts(i-1)
210 end do
211 endif
212
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.)
216 end if
217 ! Now gather the actual data with proper displacement handling
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)
221
222 if (pe_rank == root_rank) then
223 ! Create filename with iteration index
224 write(filename, '(A,I0.6,A)') 'design_iter_', idx, '.txt'
225 ! Open file for writing
226 open(newunit=iunit, file=trim(filename), status='replace', action='write')
227 ! Write header
228 write(iunit, '(A,I0)') '# Iteration: ', idx
229 write(iunit, '(A)') '# x_position height'
230
231 ! Write data
232 do i = 1, global_size
233 write(iunit, '(2E16.8)') global_x(i), global_values(i)
234 end do
235
236 close(iunit)
237
238 deallocate(global_values, global_x, recvcounts, displs)
239 print *, "Design written to ", trim(filename)
240 else
241 deallocate(global_values, recvcounts, displs)
242 endif
243 end subroutine design_3dto1d_write
244
245end module design_3dto1d
Implements the design_t.
Definition design.f90:36
Mappings to be applied to a scalar field.
Definition mapping.f90:36
Some common Masking operations we may need.
Definition mask_ops.f90:36
A PDE based filter.
A RAMP mapping of coefficients.
Implements the steady_problem_t type.
An abstract design type.
Definition design.f90:54
A topology optimization design variable.
Base abstract class for mapping.
Definition mapping.f90:46
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....