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