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
63 module procedure set_optimization_ic_int
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
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
98 call json_get(params,
"type", type)
100 select case (trim(type))
102 call json_get(params,
'value', ic_value)
103 call set_optimization_ic_uniform(fld, ic_value)
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)
110 call set_optimization_ic_point_zone(fld, ic_value, read_str, zone_value)
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")
118 call set_optimization_ic_fld(fld, fname, interpolate, tol, mesh_fname)
121 call neko_error(
'Invalid initial condition')
125 call set_optimization_ic_common(fld, coef, gs)
127 end subroutine set_optimization_ic_int
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
142 if (neko_bcknd_device .eq. 1)
then
143 call device_memcpy(fld%x, fld%x_d, n, &
144 host_to_device, sync = .false.)
148 call gs%op(fld%x, n, gs_op_add)
150 if (neko_bcknd_device .eq. 1)
then
151 call device_col2(fld%x_d, coef%mult_d, n)
153 call col2(fld%x, coef%mult, n)
156 end subroutine set_optimization_ic_common
162 subroutine set_optimization_ic_uniform(fld, ic_value)
163 type(field_t),
intent(inout) :: fld
164 real(kind=rp),
intent(in) :: ic_value
166 character(len=LOG_SIZE) :: log_buf
168 call neko_log%message(
"Type : uniform")
169 write (log_buf,
'(A,ES12.6)')
"Value: ", ic_value
170 call neko_log%message(log_buf)
174 if (neko_bcknd_device .eq. 1)
then
175 call cfill(fld%x, ic_value, n)
178 end subroutine set_optimization_ic_uniform
187 subroutine set_optimization_ic_point_zone(fld, base_value, zone_name, &
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
195 character(len=LOG_SIZE) :: log_buf
196 class(point_zone_t),
pointer :: zone
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)
206 size = fld%dof%size()
207 zone => neko_point_zone_registry%get_point_zone(trim(zone_name))
209 call set_optimization_ic_uniform(fld, base_value)
210 call cfill_mask(fld%x, zone_value,
size, zone%mask, zone%size)
212 end subroutine set_optimization_ic_point_zone
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
236 character(len=LOG_SIZE) :: log_buf
237 integer :: sample_idx, sample_mesh_idx
238 type(fld_file_data_t) :: fld_data
240 logical :: mesh_mismatch
243 type(global_interpolation_t) :: global_interp
247 type(space_t) :: prev_Xh
248 type(interpolator_t) :: space_interp
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)
257 sample_idx = extract_fld_file_index(file_name, -1)
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'")
265 call filename_chsuffix(file_name, trim(file_name),
'fld')
268 f = file_t(trim(file_name))
270 if (interpolate)
then
273 if (mesh_file_name .eq.
"none")
then
274 mesh_file_name = trim(file_name)
275 sample_mesh_idx = sample_idx
279 sample_mesh_idx = extract_fld_file_index(mesh_file_name, -1)
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'")
286 write (log_buf,
'(A,ES12.6)')
"Tolerance : ", tolerance
287 call neko_log%message(log_buf)
288 write (log_buf,
'(A,A)')
"Mesh file : ", &
290 call neko_log%message(log_buf)
295 if (sample_mesh_idx .ne. sample_idx)
then
296 call f%set_counter(sample_mesh_idx)
297 call f%read(fld_data)
303 call f%set_counter(sample_idx)
304 call f%read(fld_data)
312 mesh_mismatch = (fld_data%glb_nelv .ne. fld%msh%glb_nelv .or. &
313 fld_data%gdim .ne. fld%msh%gdim)
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.")
325 if (interpolate)
then
327 select type (ft => f%file_type)
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.")
341 global_interp = fld_data%generate_interpolator(fld%dof, fld%msh, &
347 call global_interp%evaluate(fld%x, fld_data%u%x)
348 call global_interp%free()
353 call prev_xh%init(gll, fld_data%lx, fld_data%ly, fld_data%lz)
354 call space_interp%init(fld%Xh, prev_xh)
357 call space_interp%map_host(fld%x, fld_data%t%x, fld_data%nelv, fld%Xh)
359 call space_interp%free()
365 end subroutine set_optimization_ic_fld
Optimization initial condition.