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