35 use num_types,
only: rp, sp
36 use field,
only: field_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
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
56 use mpi_f08,
only: mpi_comm_size, mpi_comm_world
57 use utils,
only: neko_error
65 type(vector_t) :: values
66 type(vector_t) :: x_coord
67 type(vector_t) :: y_coord
68 type(vector_t) :: z_coord
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
88 procedure, pass(this) :: add_mapping => design_simple_add_mapping
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
100 procedure, pass(this) :: update_design => design_simple_update_design
103 procedure, pass(this) :: map_forward => design_simple_map_forward
106 procedure, pass(this) :: map_backward => design_simple_map_backward
109 procedure, pass(this) :: write => design_simple_write
112 procedure, pass(this) :: free => design_simple_free
119 subroutine design_simple_init_from_json(this, parameters)
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
127 call json_get(parameters,
'domain.type', type)
128 call json_get_or_default(parameters,
'name', name,
'Simple Design')
130 select case (trim(type))
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)
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)
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.)
168 call this%init_from_components(name, n, x, y, z)
170 end subroutine design_simple_init_from_json
172 subroutine design_simple_init_from_components(this, name, n, x, y, z)
174 character(len=*),
intent(in) :: name
175 integer,
intent(in) :: n
176 type(vector_t),
intent(in) :: x, y, z
177 integer :: nproc, ierr
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')
185 if (pe_size .ne. 1)
then
186 call neko_error(
"Simple design can only be used with a single MPI " // &
190 call this%init_base(name, n)
192 call this%values%init(n)
197 end subroutine design_simple_init_from_components
200 subroutine design_simple_free(this)
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()
209 end subroutine design_simple_free
212 subroutine design_simple_add_mapping(this, parameters, simulation)
214 type(json_file),
intent(inout) :: parameters
217 end subroutine design_simple_add_mapping
220 subroutine design_simple_map_forward(this)
224 end subroutine design_simple_map_forward
226 subroutine design_simple_get_values(this, values)
228 type(vector_t),
intent(inout) :: values
232 end subroutine design_simple_get_values
234 subroutine design_simple_get_x(this, x)
236 type(vector_t),
intent(inout) :: x
240 end subroutine design_simple_get_x
242 subroutine design_simple_get_y(this, y)
244 type(vector_t),
intent(inout) :: y
248 end subroutine design_simple_get_y
250 subroutine design_simple_get_z(this, z)
252 type(vector_t),
intent(inout) :: z
256 end subroutine design_simple_get_z
258 subroutine design_simple_update_design(this, values)
260 type(vector_t),
intent(inout) :: values
264 end subroutine design_simple_update_design
266 subroutine design_simple_map_backward(this, sensitivity)
268 type(vector_t),
intent(in) :: sensitivity
270 end subroutine design_simple_map_backward
272 subroutine design_simple_write(this, idx)
274 integer,
intent(in) :: idx
277 integer :: npts, nx, ny, nz
278 character(len=100) :: filename
279 integer :: ios, funit
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.)
298 write(filename,
'(A,I4.4,A)')
'design_', idx,
'.vtk'
300 open(newunit=funit, file=filename, status=
'replace', action=
'write', iostat=ios)
302 call neko_error(
'Error opening file ' // filename)
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
313 write(funit,
'(A,1X,I0,1X,A)')
'POINTS', npts,
'float'
315 write(funit,
'(3(F20.12,1X))') this%x_coord%x(i), this%y_coord%x(i), &
320 write(funit,
'(A,1X,I0)')
'POINT_DATA', npts
321 write(funit,
'(A)')
'SCALARS density float 1'
322 write(funit,
'(A)')
'LOOKUP_TABLE default'
324 write(funit,
'(F20.12)') this%values%x(i)
328 end subroutine design_simple_write
330end 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.