46 use json_module,
only : json_file
83 type(
field_t),
intent(inout) :: s
84 type(
coef_t),
intent(in) :: coef
85 type(
gs_t),
intent(inout) :: gs
86 character(len=*) :: type
87 type(json_file),
intent(inout) :: params
88 integer,
intent(in) :: i
91 real(kind=
rp) :: ic_value
92 character(len=:),
allocatable :: read_str
93 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
94 real(kind=
rp) :: zone_value, tol
95 logical :: interpolate
97 if (trim(type) .eq.
'uniform')
then
99 call json_get(params,
'value', ic_value)
102 else if (trim(type) .eq.
'point_zone')
then
104 call json_get(params,
'base_value', ic_value)
105 call json_get(params,
'zone_name', read_str)
106 call json_get(params,
'zone_value', zone_value)
110 else if (trim(type) .eq.
'field')
then
112 call json_get(params,
'file_name', read_str)
113 fname = trim(read_str)
119 mesh_fname = trim(read_str)
139 character(len=*),
intent(in) :: scheme_name
140 type(
field_t),
target,
intent(inout) :: s
141 type(
coef_t),
intent(in) :: coef
142 type(
gs_t),
intent(inout) :: gs
149 call fields%assign_to_field(1, s)
151 call user_proc(scheme_name, fields)
163 type(
field_t),
intent(inout) :: s
164 type(
coef_t),
intent(in) :: coef
165 type(
gs_t),
intent(inout) :: gs
174 call gs%op(s%x, n, gs_op_add)
179 call col2(s%x, coef%mult, n)
189 type(
field_t),
intent(inout) :: s
190 real(kind=
rp),
intent(in) :: ic_value
192 character(len=LOG_SIZE) :: log_buf
194 call neko_log%message(
"Type : uniform")
195 write (log_buf,
'(A,ES12.6)')
"Value: ", ic_value
201 call cfill(s%x, ic_value, n)
214 type(
field_t),
intent(inout) :: s
215 real(kind=
rp),
intent(in) :: base_value
216 character(len=*),
intent(in) :: zone_name
217 real(kind=
rp),
intent(in) :: zone_value
220 character(len=LOG_SIZE) :: log_buf
224 call neko_log%message(
"Type : point_zone")
225 write (log_buf,
'(A,ES12.6)')
"Base value: ", base_value
227 call neko_log%message(
"Zone name : " // trim(zone_name))
228 write (log_buf,
'(A,ES12.6)')
"Zone value: ", zone_value
235 call cfill_mask(s%x, zone_value,
size, zone%mask%get(), zone%size)
255 interpolate, tolerance, mesh_file_name, i)
256 type(
field_t),
intent(inout) :: s
257 character(len=*),
intent(in) :: file_name
258 logical,
intent(in) :: interpolate
259 real(kind=
rp),
intent(in) :: tolerance
260 character(len=*),
intent(inout) :: mesh_file_name
261 integer,
intent(in) :: i
263 character(len=LOG_SIZE) :: log_buf
264 integer :: sample_idx, sample_mesh_idx
265 integer :: last_index
268 logical :: mesh_mismatch
279 call neko_log%message(
"Type : field")
280 call neko_log%message(
"File name : " // trim(file_name))
281 write (log_buf,
'(A,L1)')
"Interpolation : ", interpolate
287 if (sample_idx .eq. -1)
then
288 call neko_error(
"Invalid file name for the initial condition. The " // &
289 "file format must be e.g. 'mean0.f00001'")
296 call f%init(trim(file_name))
298 if (interpolate)
then
301 if (mesh_file_name .eq.
"none")
then
302 mesh_file_name = trim(file_name)
303 sample_mesh_idx = sample_idx
309 if (sample_mesh_idx .eq. -1)
then
310 call neko_error(
"Invalid file name for the initial condition." // &
311 " The file format must be e.g. 'mean0.f00001'")
314 write (log_buf,
'(A,ES12.6)')
"Tolerance : ", tolerance
316 write (log_buf,
'(A,A)')
"Mesh file : ", &
323 if (sample_mesh_idx .ne. sample_idx)
then
324 call f%set_counter(sample_mesh_idx)
325 call f%read(fld_data)
331 call f%set_counter(sample_idx)
332 call f%read(fld_data)
340 mesh_mismatch = (fld_data%glb_nelv .ne. s%msh%glb_nelv .or. &
341 fld_data%gdim .ne. s%msh%gdim)
343 if (mesh_mismatch .and. .not. interpolate)
then
344 call neko_error(
"The fld file must match the current mesh! " // &
345 "Use 'interpolate': 'true' to enable interpolation.")
346 else if (.not. mesh_mismatch .and. interpolate)
then
347 call neko_log%warning(
"You have activated interpolation but you " // &
348 "might still be using the same mesh.")
353 if (interpolate)
then
355 select type (ft => f%file_type)
357 if (.not. ft%dp_precision)
then
358 call neko_warning(
"The coordinates read from the field file " // &
359 "are in single precision.")
360 call neko_log%message(
"It is recommended to use a mesh in " // &
361 "double precision for better interpolation results.")
362 call neko_log%message(
"If the interpolation does not work, " // &
363 "you can try to increase the tolerance.")
375 call fld_data%generate_interpolator(global_interp, s%dof, s%msh, &
382 call global_interp%evaluate(s%x, fld_data%s(i)%x, .false.)
384 call global_interp%evaluate(s%x, fld_data%t%x, .false.)
387 call global_interp%free
395 call prev_xh%init(
gll, fld_data%lx, fld_data%ly, fld_data%lz)
396 call space_interp%init(s%Xh, prev_xh)
401 call space_interp%map_host(s%x, fld_data%s(i)%x, fld_data%nelv, s%Xh)
403 call space_interp%map_host(s%x, fld_data%t%x, fld_data%nelv, s%Xh)
406 call space_interp%free
Copy data between host and device (or device and device)
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Abstract interface for user defined initial conditions.
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
integer, parameter, public device_to_host
Module for file I/O operations.
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
Implements global_interpolation given a dofmap.
Routines to interpolate between different spaces.
Utilities for retrieving parameters from the case files.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
subroutine, public col2(a, b, n)
Vector multiplication .
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
Scalar initial condition.
subroutine set_scalar_ic_uniform(s, ic_value)
Uniform initial condition.
subroutine set_scalar_ic_point_zone(s, base_value, zone_name, zone_value)
Point zone initial condition.
subroutine set_scalar_ic_common(s, coef, gs)
Set scalar initial condition (common)
subroutine set_scalar_ic_int(s, coef, gs, type, params, i)
Set scalar initial condition (builtin)
subroutine set_scalar_ic_fld(s, file_name, interpolate, tolerance, mesh_file_name, i)
Set the initial condition of the scalar based on a field. @detail The field is read from an fld file....
subroutine set_scalar_ic_usr(scheme_name, s, coef, gs, user_proc)
Set scalar intial condition (user defined)
Defines a function space.
integer, parameter, public gll
Interfaces for user interaction with NEKO.
integer function, public extract_fld_file_index(fld_filename, default_index)
Extracts the index of a field file. For example, "myfield.f00045" will return 45. If the suffix of th...
integer, parameter, public neko_fname_len
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
subroutine, public filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
subroutine, public filename_suffix(fname, suffix)
Extract a filename's suffix.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
field_list_t, To be able to group fields together
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Interface for NEKTON fld files.
Implements global interpolation for arbitrary points in the domain.
Interpolation between two space::space_t.
Base abstract type for point zones.
The function space for the SEM solution fields.