Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
set_optimization_ic.f90
1! Copyright (c) 2021, 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!
35 use gather_scatter, only : gs_t, gs_op_add
36 use neko_config, only : neko_bcknd_device
37 use num_types, only : rp
38 use device_math, only : device_col2
39 use device, only : device_memcpy, host_to_device
40 use field, only : field_t
41 use utils, only : neko_error, filename_chsuffix, filename_suffix, &
42 neko_warning, neko_fname_len, extract_fld_file_index
43 use coefs, only : coef_t
44 use math, only : col2, cfill, cfill_mask
45 use user_intf, only : useric_scalar
46 use json_module, only : json_file
47 use json_utils, only: json_get, json_get_or_default
48 use point_zone, only: point_zone_t
49 use point_zone_registry, only: neko_point_zone_registry
50 use field_registry, only: neko_field_registry
51 use logger, only: neko_log, log_size
52 use fld_file_data, only: fld_file_data_t
53 use fld_file, only: fld_file_t
54 use checkpoint, only: chkp_t
55 use file, only: file_t
56 use global_interpolation, only: global_interpolation_t
57 use interpolation, only: interpolator_t
58 use space, only: space_t, gll
59 implicit none
60 private
61
63 module procedure set_optimization_ic_int
64 end interface set_optimization_ic
65
66 public :: set_optimization_ic
67
68contains
69
84 subroutine set_optimization_ic_int(fld, coef, gs, params)
85 type(field_t), intent(inout) :: fld
86 type(coef_t), intent(in) :: coef
87 type(gs_t), intent(inout) :: gs
88 type(json_file), intent(inout) :: params
89
90 ! Variables for retrieving JSON parameters
91 character(len=:), allocatable :: type
92 real(kind=rp) :: ic_value
93 character(len=:), allocatable :: read_str
94 character(len=:), allocatable :: fname, mesh_fname
95 real(kind=rp) :: zone_value, tol
96 logical :: interpolate
97
98 call json_get(params, "type", type)
99
100 select case (trim(type))
101 case ('uniform')
102 call json_get(params, 'value', ic_value)
103 call set_optimization_ic_uniform(fld, ic_value)
104
105 case ('point_zone')
106 call json_get(params, 'base_value', ic_value)
107 call json_get(params, 'zone_name', read_str)
108 call json_get(params, 'zone_value', zone_value)
109
110 call set_optimization_ic_point_zone(fld, ic_value, read_str, zone_value)
111
112 case ('field')
113 call json_get(params, 'file_name', fname)
114 call json_get_or_default(params, 'interpolate', interpolate, .false.)
115 call json_get_or_default(params, 'tolerance', tol, 0.000001_rp)
116 call json_get_or_default(params, 'mesh_file_name', mesh_fname, "none")
117
118 call set_optimization_ic_fld(fld, fname, interpolate, tol, mesh_fname)
119
120 case default
121 call neko_error('Invalid initial condition')
122
123 end select
124
125 call set_optimization_ic_common(fld, coef, gs)
126
127 end subroutine set_optimization_ic_int
128
135 subroutine set_optimization_ic_common(fld, coef, gs)
136 type(field_t), intent(inout) :: fld
137 type(coef_t), intent(in) :: coef
138 type(gs_t), intent(inout) :: gs
139 integer :: n
140
141 n = fld%dof%size()
142 if (neko_bcknd_device .eq. 1) then
143 call device_memcpy(fld%x, fld%x_d, n, &
144 host_to_device, sync = .false.)
145 end if
146
147 ! Ensure continuity across elements for initial conditions
148 call gs%op(fld%x, n, gs_op_add)
149
150 if (neko_bcknd_device .eq. 1) then
151 call device_col2(fld%x_d, coef%mult_d, n)
152 else
153 call col2(fld%x, coef%mult, n)
154 end if
155
156 end subroutine set_optimization_ic_common
157
162 subroutine set_optimization_ic_uniform(fld, ic_value)
163 type(field_t), intent(inout) :: fld
164 real(kind=rp), intent(in) :: ic_value
165 integer :: n
166 character(len=LOG_SIZE) :: log_buf
167
168 call neko_log%message("Type : uniform")
169 write (log_buf, '(A,ES12.6)') "Value: ", ic_value
170 call neko_log%message(log_buf)
171
172 fld = ic_value
173 n = fld%dof%size()
174 if (neko_bcknd_device .eq. 1) then
175 call cfill(fld%x, ic_value, n)
176 end if
177
178 end subroutine set_optimization_ic_uniform
179
187 subroutine set_optimization_ic_point_zone(fld, base_value, zone_name, &
188 zone_value)
189 type(field_t), intent(inout) :: fld
190 real(kind=rp), intent(in) :: base_value
191 character(len=*), intent(in) :: zone_name
192 real(kind=rp), intent(in) :: zone_value
193
194 ! Internal variables
195 character(len=LOG_SIZE) :: log_buf
196 class(point_zone_t), pointer :: zone
197 integer :: size
198
199 call neko_log%message("Type : point_zone")
200 write (log_buf, '(A,ES12.6)') "Base value: ", base_value
201 call neko_log%message(log_buf)
202 call neko_log%message("Zone name : " // trim(zone_name))
203 write (log_buf, '(A,ES12.6)') "Zone value: ", zone_value
204 call neko_log%message(log_buf)
205
206 size = fld%dof%size()
207 zone => neko_point_zone_registry%get_point_zone(trim(zone_name))
208
209 call set_optimization_ic_uniform(fld, base_value)
210 call cfill_mask(fld%x, zone_value, size, zone%mask, zone%size)
211
212 end subroutine set_optimization_ic_point_zone
213
228 subroutine set_optimization_ic_fld(fld, file_name, &
229 interpolate, tolerance, mesh_file_name)
230 type(field_t), intent(inout) :: fld
231 character(len=*), intent(in) :: file_name
232 logical, intent(in) :: interpolate
233 real(kind=rp), intent(in) :: tolerance
234 character(len=*), intent(inout) :: mesh_file_name
235
236 character(len=LOG_SIZE) :: log_buf
237 integer :: sample_idx, sample_mesh_idx
238 type(fld_file_data_t) :: fld_data
239 type(file_t) :: f
240 logical :: mesh_mismatch
241
242 ! ---- For the mesh to mesh interpolation
243 type(global_interpolation_t) :: global_interp
244 ! -----
245
246 ! ---- For space to space interpolation
247 type(space_t) :: prev_Xh
248 type(interpolator_t) :: space_interp
249 ! ----
250
251 call neko_log%message("Type : field")
252 call neko_log%message("File name : " // trim(file_name))
253 write (log_buf, '(A,L1)') "Interpolation: ", interpolate
254 call neko_log%message(log_buf)
255
256 ! Extract sample index from the file name
257 sample_idx = extract_fld_file_index(file_name, -1)
258
259 if (sample_idx .eq. -1) then
260 call neko_error("Invalid file name for the initial condition. The " // &
261 "file format must be e.g. 'mean0.f00001'")
262 end if
263
264 ! Change from "field0.f000*" to "field0.fld" for the fld reader
265 call filename_chsuffix(file_name, trim(file_name), 'fld')
266
267 call fld_data%init()
268 f = file_t(trim(file_name))
269
270 if (interpolate) then
271
272 ! If no mesh file is specified, use the default file name
273 if (mesh_file_name .eq. "none") then
274 mesh_file_name = trim(file_name)
275 sample_mesh_idx = sample_idx
276 else
277
278 ! Extract sample index from the mesh file name
279 sample_mesh_idx = extract_fld_file_index(mesh_file_name, -1)
280
281 if (sample_mesh_idx .eq. -1) then
282 call neko_error("Invalid file name for the initial condition." // &
283 " The file format must be e.g. 'mean0.f00001'")
284 end if
285
286 write (log_buf, '(A,ES12.6)') "Tolerance : ", tolerance
287 call neko_log%message(log_buf)
288 write (log_buf, '(A,A)') "Mesh file : ", &
289 trim(mesh_file_name)
290 call neko_log%message(log_buf)
291
292 end if ! if mesh_file_name .eq. none
293
294 ! Read the mesh coordinates if they are not in our fld file
295 if (sample_mesh_idx .ne. sample_idx) then
296 call f%set_counter(sample_mesh_idx)
297 call f%read(fld_data)
298 end if
299
300 end if
301
302 ! Read the field file containing (u,v,w,p)
303 call f%set_counter(sample_idx)
304 call f%read(fld_data)
305
306 !
307 ! Check that the data in the fld file matches the current case.
308 ! Note that this is a safeguard and there are corner cases where
309 ! two different meshes have the same dimension and same # of elements
310 ! but this should be enough to cover obvious cases.
311 !
312 mesh_mismatch = (fld_data%glb_nelv .ne. fld%msh%glb_nelv .or. &
313 fld_data%gdim .ne. fld%msh%gdim)
314
315 if (mesh_mismatch .and. .not. interpolate) then
316 call neko_error("The fld file must match the current mesh! " // &
317 "Use 'interpolate': 'true' to enable interpolation.")
318 else if (.not. mesh_mismatch .and. interpolate) then
319 call neko_log%warning("You have activated interpolation but you " // &
320 "might still be using the same mesh.")
321 end if
322
323
324 ! Mesh interpolation if specified
325 if (interpolate) then
326 ! Issue a warning if the mesh is in single precision
327 select type (ft => f%file_type)
328 type is (fld_file_t)
329 if (.not. ft%dp_precision) then
330 call neko_warning("The coordinates read from the field file " // &
331 "are in single precision.")
332 call neko_log%message("It is recommended to use a mesh in " // &
333 "double precision for better interpolation results.")
334 call neko_log%message("If the interpolation does not work, " // &
335 "you can try to increase the tolerance.")
336 end if
337 class default
338 end select
339
340 ! Generates an interpolator object and performs the point search
341 global_interp = fld_data%generate_interpolator(fld%dof, fld%msh, &
342 tolerance)
343
344 ! Evaluate design
345 ! NOTE. Is is currently the convention in the design0.f* files that
346 ! the unfiltered design should be stored in the u component.
347 call global_interp%evaluate(fld%x, fld_data%u%x)
348 call global_interp%free()
349
350 else ! No interpolation, just potentially from different spaces
351
352 ! Build a space_t object from the data in the fld file
353 call prev_xh%init(gll, fld_data%lx, fld_data%ly, fld_data%lz)
354 call space_interp%init(fld%Xh, prev_xh)
355
356 ! Do the space-to-space interpolation
357 call space_interp%map_host(fld%x, fld_data%t%x, fld_data%nelv, fld%Xh)
358
359 call space_interp%free()
360
361 end if
362
363 call fld_data%free()
364
365 end subroutine set_optimization_ic_fld
366
367end module optimization_ic
Optimization initial condition.