Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
design_3dto1d.f90
1! Copyright (c) 2024, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32
33! Implements the `design_3dto1d_t` type.
34module design_3dto1d
35 use num_types, only: rp, sp
36 use json_module, only: json_file
37 use mapping, only: mapping_t
38 use pde_filter, only: pde_filter_t
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
48 use design, only: design_t
49 use math, only: rzero
50 use simulation_m, only: simulation_t
51 use json_module, only: json_file
52 use json_utils, only: json_get
53
54 use vector, only: vector_t
55 use math, only: copy
56
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
59
60 use fld_file_output, only: fld_file_output_t
61
62 implicit none
63 private
64
66 type, extends(design_t), public :: design_3dto1d_t
67 private
68
69 type(vector_t) :: values
70
71 contains
72
73 ! ----------------------------------------------------------------------- !
74 ! Initializations
75
77 generic, public :: init => init_from_components
79 procedure, pass(this) :: init_from_components => &
80 design_3dto1d_init_from_components
81
83 procedure, pass(this) :: get_values => design_3dto1d_get_values
84
86 procedure, pass(this) :: update_design => design_3dto1d_update_design
87
89 procedure, pass(this) :: write => design_3dto1d_write
90
92 procedure, pass(this) :: free => design_3dto1d_free
93
94
96 procedure, pass(this) :: map_forward => design_3dto1d_map_forward
98 procedure, pass(this) :: map_backward => design_3dto1d_map_backward
99
100 end type design_3dto1d_t
101
102contains
103
104
105 subroutine design_3dto1d_init_from_components(this, n)
106 class(design_3dto1d_t), intent(inout) :: this
107 integer, intent(in) :: n
108
109 call this%init_base('design_3dto1d', n)
110
111 call this%values%init(n)
112
113 end subroutine design_3dto1d_init_from_components
114
116 subroutine design_3dto1d_free(this)
117 class(design_3dto1d_t), intent(inout) :: this
118
119 call this%free_base()
120 call this%values%free()
121 end subroutine design_3dto1d_free
122
123 subroutine design_3dto1d_map_forward(this)
124 class(design_3dto1d_t), intent(inout) :: this
125
126 end subroutine design_3dto1d_map_forward
127
128 subroutine design_3dto1d_map_backward(this, sensitivity)
129 class(design_3dto1d_t), intent(inout) :: this
130 type(vector_t), intent(in) :: sensitivity
131 end subroutine design_3dto1d_map_backward
132
133
134 subroutine design_3dto1d_get_values(this, values)
135 class(design_3dto1d_t), intent(in) :: this
136 type(vector_t), intent(inout) :: values
137
138 values = this%values
139
140 end subroutine design_3dto1d_get_values
141
142 subroutine design_3dto1d_update_design(this, values)
143 class(design_3dto1d_t), intent(inout) :: this
144 type(vector_t), intent(inout) :: values
145 integer :: n
146
147 n = this%size()
148
149 if (neko_bcknd_device .eq. 1) then
150 call device_copy(this%values%x_d, values%x_d, n)
151 else
152 this%values = values
153 end if
154
155
156 end subroutine design_3dto1d_update_design
157
158 subroutine design_3dto1d_write(this, idx)
159 class(design_3dto1d_t), intent(inout) :: this
160 integer, intent(in) :: idx
161
162 character(len=100) :: filename
163 integer :: i, iunit, ierr
164 real(rp) :: le
165 real(rp), allocatable :: global_values(:)
166 real(rp), allocatable :: global_x(:)
167 integer :: global_size, local_size
168 real(rp) :: l_total
169 integer, allocatable :: recvcounts(:), displs(:)
170 integer :: root_rank = 0
171
172 l_total = 2.0_rp
173
174 ! Get local size on all ranks
175 local_size = this%size()
176
177 ! Only rank 0 handles global arrays
178 if (pe_rank == root_rank) then
179 ! Get global size
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))
184 ! Calculate element length and x positions
185 le = l_total / real(global_size, kind=rp)
186 do i = 1, global_size
187 global_x(i) = le * (i - 0.5_rp) ! Center of each element
188 end do
189 else
190 ! Non-root ranks: minimal dummy allocations
191 allocate(global_values(1), recvcounts(1), displs(1))
192 endif
193
194 ! First, gather all the local sizes to rank 0
195 call mpi_gather(local_size, 1, mpi_integer, &
196 recvcounts, 1, mpi_integer, &
197 root_rank, neko_comm, ierr)
198
199 ! Calculate displacements on rank 0
200 if (pe_rank == root_rank) then
201 displs(1) = 0
202 do i = 2, pe_size
203 displs(i) = displs(i-1) + recvcounts(i-1)
204 end do
205 endif
206
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.)
210 end if
211 ! Now gather the actual data with proper displacement handling
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)
215
216 if (pe_rank == root_rank) then
217 ! Create filename with iteration index
218 write(filename, '(A,I0.6,A)') 'design_iter_', idx, '.txt'
219 ! Open file for writing
220 open(newunit=iunit, file=trim(filename), status='replace', action='write')
221 ! Write header
222 write(iunit, '(A,I0)') '# Iteration: ', idx
223 write(iunit, '(A)') '# x_position height'
224
225 ! Write data
226 do i = 1, global_size
227 write(iunit, '(2E16.8)') global_x(i), global_values(i)
228 end do
229
230 close(iunit)
231
232 deallocate(global_values, global_x, recvcounts, displs)
233 print *, "Design written to ", trim(filename)
234 else
235 deallocate(global_values, recvcounts, displs)
236 endif
237 end subroutine design_3dto1d_write
238
239end module design_3dto1d
Implements the design_t.
Definition design.f90:34
Mappings to be applied to a scalar field.
Definition mapping.f90:35
Some common Masking operations we may need.
Definition mask_ops.f90:34
A PDE based filter.
A RAMP mapping of coefficients.
Implements the steady_problem_t type.
An abstract design type.
Definition design.f90:52
A topology optimization design variable.
Base abstract class for mapping.
Definition mapping.f90:45
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....