Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
design_simple.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 `simple_design_t` type.
34module simple_design
35 use num_types, only: rp, sp
36 use field, only: field_t
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 design, only: design_t
48 use math, only: rzero
49 use simulation_m, only: simulation_t
50 use comm, only: pe_size
51 use json_module, only: json_file
52 use json_utils, only: json_get, json_get_or_default
54 use vector, only: vector_t
55 use math, only: copy
56 use field_registry, only: neko_field_registry
57 use mpi_f08, only: mpi_comm_size, mpi_comm_world
58 use utils, only: neko_error
59 implicit none
60 private
61
63 type, extends(design_t), public :: simple_design_t
64 private
65
66 type(vector_t) :: values
67 type(vector_t) :: x_coord
68 type(vector_t) :: y_coord
69 type(vector_t) :: z_coord
70
71 ! needed to write the design into a vtk file with connectivity
72 integer :: nx, ny, nz
73
74 contains
75
76 ! ----------------------------------------------------------------------- !
77 ! Initializations
78
80 generic, public :: init => init_from_json, init_from_components
82 procedure, pass(this) :: init_from_json => &
83 design_simple_init_from_json
85 procedure, pass(this) :: init_from_components => &
86 design_simple_init_from_components
87
89 procedure, pass(this) :: add_mapping => design_simple_add_mapping
90
92 procedure, pass(this) :: get_values => design_simple_get_values
94 procedure, pass(this) :: design_get_x => design_simple_get_x
96 procedure, pass(this) :: design_get_y => design_simple_get_y
98 procedure, pass(this) :: design_get_z => design_simple_get_z
99
101 procedure, pass(this) :: update_design => design_simple_update_design
102
104 procedure, pass(this) :: map_forward => design_simple_map_forward
105
107 procedure, pass(this) :: map_backward => design_simple_map_backward
108
110 procedure, pass(this) :: write => design_simple_write
111
113 procedure, pass(this) :: free => design_simple_free
114
115 end type simple_design_t
116
117contains
118
120 subroutine design_simple_init_from_json(this, parameters)
121 class(simple_design_t), intent(inout) :: this
122 type(json_file), intent(inout) :: parameters
123 character(len=:), allocatable :: type, name
124 integer :: n, nx, ny, nz, i, j, k, index
125 real(kind=rp), dimension(:), allocatable :: limits
126 type(vector_t) :: x, y, z
127
128 call json_get(parameters, 'domain.type', type)
129 call json_get_or_default(parameters, 'name', name, 'Simple Design')
130
131 select case (trim(type))
132 case ("box")
133 call json_get(parameters, 'optimization.design.domain.nx', nx)
134 call json_get(parameters, 'optimization.design.domain.ny', ny)
135 call json_get(parameters, 'optimization.design.domain.nz', nz)
136 call json_get(parameters, 'optimization.design.domain.limits', limits)
137 n = (nx + 1) * (ny + 1) * (nz + 1)
138 this%nx = nx
139 this%ny = ny
140 this%nz = nz
141
142
143 call x%init(n)
144 call y%init(n)
145 call z%init(n)
146 index = 1
147 do k = 1, nz + 1
148 do j = 1, ny + 1
149 do i = 1, nx + 1
150 x%x(index) = limits(1) + (limits(2) - limits(1)) * &
151 real(i - 1, kind=rp) / real(nx, kind=rp)
152 y%x(index) = limits(3) + (limits(4) - limits(3)) * &
153 real(j - 1, kind=rp) / real(ny, kind=rp)
154 z%x(index) = limits(5) + (limits(6) - limits(5)) * &
155 real(k - 1, kind=rp) / real(nz, kind=rp)
156 index = index + 1
157 end do
158 end do
159 end do
160
161 if (neko_bcknd_device .eq. 1) then
162 call device_memcpy(x%x, x%x_d, n, host_to_device, sync = .false.)
163 call device_memcpy(y%x, y%x_d, n, host_to_device, sync = .false.)
164 call device_memcpy(z%x, z%x_d, n, host_to_device, sync = .true.)
165 end if
166
167 end select
168
169 call this%init_from_components(name, n, x, y, z)
170
171 end subroutine design_simple_init_from_json
172
173 subroutine design_simple_init_from_components(this, name, n, x, y, z)
174 class(simple_design_t), intent(inout) :: this
175 character(len=*), intent(in) :: name
176 integer, intent(in) :: n
177 type(vector_t), intent(in) :: x, y, z
178 integer :: nproc, ierr
179
180 ! Throw an error if we run on multiple MPI ranks
181 call mpi_comm_size(mpi_comm_world, nproc, ierr)
182 if (nproc .gt. 1) then
183 call neko_error('simple_design_t can only be used on a single MPI rank')
184 end if
185
186 if (pe_size .ne. 1) then
187 call neko_error("Simple design can only be used with a single MPI " // &
188 "process.")
189 end if
190
191 call this%init_base(name, n)
192
193 call this%values%init(n)
194 this%x_coord = x
195 this%y_coord = y
196 this%z_coord = z
197
198 end subroutine design_simple_init_from_components
199
201 subroutine design_simple_free(this)
202 class(simple_design_t), intent(inout) :: this
203
204 call this%free_base()
205 call this%values%free()
206 call this%x_coord%free()
207 call this%y_coord%free()
208 call this%z_coord%free()
209
210 end subroutine design_simple_free
211
213 subroutine design_simple_add_mapping(this, parameters, simulation)
214 class(simple_design_t), intent(inout) :: this
215 type(json_file), intent(inout) :: parameters
216 type(simulation_t), intent(inout) :: simulation
217
218 end subroutine design_simple_add_mapping
219
220
221 subroutine design_simple_map_forward(this)
222 class(simple_design_t), intent(inout) :: this
223
224
225 end subroutine design_simple_map_forward
226
227 subroutine design_simple_get_values(this, values)
228 class(simple_design_t), intent(in) :: this
229 type(vector_t), intent(inout) :: values
230
231 values = this%values
232
233 end subroutine design_simple_get_values
234
235 subroutine design_simple_get_x(this, x)
236 class(simple_design_t), intent(in) :: this
237 type(vector_t), intent(inout) :: x
238
239 x = this%x_coord
240
241 end subroutine design_simple_get_x
242
243 subroutine design_simple_get_y(this, y)
244 class(simple_design_t), intent(in) :: this
245 type(vector_t), intent(inout) :: y
246
247 y = this%y_coord
248
249 end subroutine design_simple_get_y
250
251 subroutine design_simple_get_z(this, z)
252 class(simple_design_t), intent(in) :: this
253 type(vector_t), intent(inout) :: z
254
255 z = this%z_coord
256
257 end subroutine design_simple_get_z
258
259 subroutine design_simple_update_design(this, values)
260 class(simple_design_t), intent(inout) :: this
261 type(vector_t), intent(inout) :: values
262
263 this%values = values
264
265 end subroutine design_simple_update_design
266
267 subroutine design_simple_map_backward(this, sensitivity)
268 class(simple_design_t), intent(inout) :: this
269 type(vector_t), intent(in) :: sensitivity
270
271 end subroutine design_simple_map_backward
272
273 subroutine design_simple_write(this, idx)
274 class(simple_design_t), intent(inout) :: this
275 integer, intent(in) :: idx
276
277 integer :: i
278 integer :: npts, nx, ny, nz
279 character(len=100) :: filename
280 integer :: ios, funit
281
282 nx = this%nx + 1
283 ny = this%ny + 1
284 nz = this%nz + 1
285 npts = nx * ny * nz
286
287 ! Synchronize the device memory if using a GPU backend is enabled
288 if (neko_bcknd_device .eq. 1) then
289 call device_memcpy(this%x_coord%x, this%x_coord%x_d, npts, &
290 device_to_host, sync = .false.)
291 call device_memcpy(this%y_coord%x, this%y_coord%x_d, npts, &
292 device_to_host, sync = .false.)
293 call device_memcpy(this%z_coord%x, this%z_coord%x_d, npts, &
294 device_to_host, sync = .false.)
295 call device_memcpy(this%values%x, this%values%x_d, npts, &
296 device_to_host, sync = .true.)
297 end if
298
299 write(filename, '(A,I4.4,A)') 'design_', idx, '.vtk'
300
301 open(newunit=funit, file=filename, status='replace', action='write', iostat=ios)
302 if (ios .ne. 0) then
303 call neko_error('Error opening file ' // filename)
304 end if
305
306 ! Header
307 write(funit, '(A)') '# vtk DataFile Version 3.0'
308 write(funit, '(A)') 'Simple Scalar Field'
309 write(funit, '(A)') 'ASCII'
310 write(funit, '(A)') 'DATASET STRUCTURED_GRID'
311 write(funit, '(A,3(1X,I0))') 'DIMENSIONS', nx, ny, nz
312
313 ! Points
314 write(funit, '(A,1X,I0,1X,A)') 'POINTS', npts, 'float'
315 do i = 1, npts
316 write(funit, '(3(F20.12,1X))') this%x_coord%x(i), this%y_coord%x(i), &
317 this%z_coord%x(i)
318 end do
319
320 ! Scalars
321 write(funit, '(A,1X,I0)') 'POINT_DATA', npts
322 write(funit, '(A)') 'SCALARS density float 1'
323 write(funit, '(A)') 'LOOKUP_TABLE default'
324 do i = 1, npts
325 write(funit, '(F20.12)') this%values%x(i)
326 end do
327
328 close(funit)
329 end subroutine design_simple_write
330
331end module simple_design
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 simple_brinkman_source_term_t type.
Implements the steady_problem_t type.
An abstract design type.
Definition design.f90:52
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....
A topology optimization design variable.