37 use num_types,
only: rp, sp
38 use field,
only: field_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
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
58 use mpi_f08,
only: mpi_comm_size, mpi_comm_world
59 use utils,
only: neko_error
67 type(vector_t) :: values
68 type(vector_t) :: x_coord
69 type(vector_t) :: y_coord
70 type(vector_t) :: z_coord
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
90 procedure, pass(this) :: add_mapping => design_simple_add_mapping
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
102 procedure, pass(this) :: update_design => design_simple_update_design
105 procedure, pass(this) :: map_forward => design_simple_map_forward
108 procedure, pass(this) :: map_backward => design_simple_map_backward
111 procedure, pass(this) :: write => design_simple_write
114 procedure, pass(this) :: free => design_simple_free
121 subroutine design_simple_init_from_json(this, parameters)
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
129 call json_get(parameters,
'domain.type', type)
130 call json_get_or_default(parameters,
'name', name,
'Simple Design')
132 select case (trim(type))
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)
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)
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.)
170 call this%init_from_components(name, n, x, y, z)
172 end subroutine design_simple_init_from_json
174 subroutine design_simple_init_from_components(this, name, n, x, y, z)
176 character(len=*),
intent(in) :: name
177 integer,
intent(in) :: n
178 type(vector_t),
intent(in) :: x, y, z
179 integer :: nproc, ierr
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')
187 if (pe_size .ne. 1)
then
188 call neko_error(
"Simple design can only be used with a single MPI " // &
192 call this%init_base(name, n)
194 call this%values%init(n)
199 end subroutine design_simple_init_from_components
202 subroutine design_simple_free(this)
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()
211 end subroutine design_simple_free
214 subroutine design_simple_add_mapping(this, parameters, simulation)
216 type(json_file),
intent(inout) :: parameters
219 end subroutine design_simple_add_mapping
222 subroutine design_simple_map_forward(this)
226 end subroutine design_simple_map_forward
228 subroutine design_simple_get_values(this, values)
230 type(vector_t),
intent(inout) :: values
232 if (this%size() .ne. values%size())
then
233 call neko_error(
'Get design: size mismatch')
237 end subroutine design_simple_get_values
239 subroutine design_simple_get_x(this, x)
241 type(vector_t),
intent(inout) :: x
243 if (this%size() .ne. x%size())
then
244 call neko_error(
'Get x: size mismatch')
248 end subroutine design_simple_get_x
250 subroutine design_simple_get_y(this, y)
252 type(vector_t),
intent(inout) :: y
254 if (this%size() .ne. y%size())
then
255 call neko_error(
'Get y: size mismatch')
259 end subroutine design_simple_get_y
261 subroutine design_simple_get_z(this, z)
263 type(vector_t),
intent(inout) :: z
265 if (this%size() .ne. z%size())
then
266 call neko_error(
'Get z: size mismatch')
270 end subroutine design_simple_get_z
272 subroutine design_simple_update_design(this, values)
274 type(vector_t),
intent(inout) :: values
276 if (this%size() .ne. values%size())
then
277 call neko_error(
'Update design: size mismatch')
281 end subroutine design_simple_update_design
283 subroutine design_simple_map_backward(this, sensitivity)
285 type(vector_t),
intent(in) :: sensitivity
287 end subroutine design_simple_map_backward
289 subroutine design_simple_write(this, idx)
291 integer,
intent(in) :: idx
294 integer :: npts, nx, ny, nz
295 character(len=100) :: filename
296 integer :: ios, funit
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.)
315 write(filename,
'(A,I4.4,A)')
'design_', idx,
'.vtk'
317 open(newunit=funit, file=filename, status=
'replace', action=
'write', iostat=ios)
319 call neko_error(
'Error opening file ' // filename)
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
330 write(funit,
'(A,1X,I0,1X,A)')
'POINTS', npts,
'float'
332 write(funit,
'(3(F20.12,1X))') this%x_coord%x(i), this%y_coord%x(i), &
337 write(funit,
'(A,1X,I0)')
'POINT_DATA', npts
338 write(funit,
'(A)')
'SCALARS density float 1'
339 write(funit,
'(A)')
'LOOKUP_TABLE default'
341 write(funit,
'(F20.12)') this%values%x(i)
345 end subroutine design_simple_write
347end module simple_design
Mappings to be applied to a scalar field.
Some common Masking operations we may need.
A RAMP mapping of coefficients.
Implements the simple_brinkman_source_term_t type.
Implements the steady_problem_t type.
Base abstract class for mapping.
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 simple Brinkman source term.
A topology optimization design variable.