50 use json_module,
only : json_file
78 type(
field_t),
intent(inout) :: u
79 type(
field_t),
intent(inout) :: v
80 type(
field_t),
intent(inout) :: w
81 type(
field_t),
intent(inout) :: p
82 type(
coef_t),
intent(in) :: coef
83 type(
gs_t),
intent(inout) :: gs
84 character(len=*) :: type
85 type(json_file),
intent(inout) :: params
86 real(kind=
rp) :: delta
87 real(kind=
rp),
allocatable :: uinf(:)
88 real(kind=
rp),
allocatable :: zone_value(:)
89 character(len=:),
allocatable :: read_str
95 if (trim(type) .eq.
'uniform')
then
103 else if (trim(type) .eq.
'blasius')
then
106 call json_get(params,
'approximation', read_str)
114 else if (trim(type) .eq.
'point_zone')
then
117 call json_get(params,
'zone_name', read_str)
125 else if (trim(type) .eq.
'field')
then
128 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
129 logical :: interpolate
130 type(json_file) :: interp_subdict
132 call json_get(params,
'file_name', read_str)
133 fname = trim(read_str)
139 mesh_fname = trim(read_str)
144 mesh_fname, interp_subdict)
157 type(field_t),
target,
intent(inout) :: u
158 type(field_t),
target,
intent(inout) :: v
159 type(field_t),
target,
intent(inout) :: w
160 type(field_t),
target,
intent(inout) :: p
161 type(coef_t),
intent(in) :: coef
162 type(gs_t),
intent(inout) :: gs
163 procedure(user_initial_conditions_intf) :: user_proc
164 character(len=*),
intent(in) :: scheme_name
166 type(field_list_t) :: fields
169 call neko_log%message(
"Type: user")
172 call fields%assign_to_field(1, u)
173 call fields%assign_to_field(2, v)
174 call fields%assign_to_field(3, w)
175 call fields%assign_to_field(4, p)
177 call user_proc(scheme_name, fields)
186 user_proc, scheme_name)
187 type(field_t),
target,
intent(inout) :: rho
188 type(field_t),
target,
intent(inout) :: u
189 type(field_t),
target,
intent(inout) :: v
190 type(field_t),
target,
intent(inout) :: w
191 type(field_t),
target,
intent(inout) :: p
192 type(coef_t),
intent(in) :: coef
193 type(gs_t),
intent(inout) :: gs
194 procedure(user_initial_conditions_intf) :: user_proc
195 character(len=*),
intent(in) :: scheme_name
197 type(field_list_t) :: fields
200 call neko_log%message(
"Type: user (compressible flows)")
203 call fields%assign_to_field(1, rho)
204 call fields%assign_to_field(2, u)
205 call fields%assign_to_field(3, v)
206 call fields%assign_to_field(4, w)
207 call fields%assign_to_field(5, p)
208 call user_proc(scheme_name, fields)
214 if (neko_bcknd_device .eq. 1)
then
215 call device_memcpy(p%x, p%x_d, n, host_to_device, sync = .false.)
216 call device_memcpy(rho%x, rho%x_d, n, host_to_device, sync = .false.)
221 call gs%op(p%x, p%dof%size(), gs_op_add)
222 call gs%op(rho%x, rho%dof%size(), gs_op_add)
224 if (neko_bcknd_device .eq. 1)
then
225 call device_col2(rho%x_d, coef%mult_d, rho%dof%size())
226 call device_col2(p%x_d, coef%mult_d, p%dof%size())
228 call col2(rho%x, coef%mult, rho%dof%size())
229 call col2(p%x, coef%mult, p%dof%size())
235 type(field_t),
intent(inout) :: u
236 type(field_t),
intent(inout) :: v
237 type(field_t),
intent(inout) :: w
238 type(field_t),
intent(inout) :: p
239 type(coef_t),
intent(in) :: coef
240 type(gs_t),
intent(inout) :: gs
245 if (neko_bcknd_device .eq. 1)
then
246 call u%copy_from(host_to_device, sync = .false.)
247 call v%copy_from(host_to_device, sync = .false.)
248 call w%copy_from(host_to_device, sync = .false.)
251 call p%copy_from(host_to_device, sync = .true.)
255 call rotate_cyc(u%x, v%x, w%x, 1, coef)
256 call gs%op(u%x, u%dof%size(), gs_op_add)
257 call gs%op(v%x, v%dof%size(), gs_op_add)
258 call gs%op(w%x, w%dof%size(), gs_op_add)
259 call rotate_cyc(u%x, v%x, w%x, 0, coef)
261 if (neko_bcknd_device .eq. 1)
then
262 call device_col2(u%x_d, coef%mult_d, u%dof%size())
263 call device_col2(v%x_d, coef%mult_d, v%dof%size())
264 call device_col2(w%x_d, coef%mult_d, w%dof%size())
266 call col2(u%x, coef%mult, u%dof%size())
267 call col2(v%x, coef%mult, v%dof%size())
268 call col2(w%x, coef%mult, w%dof%size())
275 type(field_t),
intent(inout) :: u
276 type(field_t),
intent(inout) :: v
277 type(field_t),
intent(inout) :: w
278 real(kind=rp),
intent(in) :: uinf(3)
280 character(len=LOG_SIZE) :: log_buf
282 call neko_log%message(
"Type : uniform")
283 write (log_buf,
'(A, 3(ES12.6, A))')
"Value: [", &
284 (uinf(i),
", ", i = 1, 2), uinf(3),
"]"
285 call neko_log%message(log_buf)
291 if (neko_bcknd_device .eq. 1)
then
292 call cfill(u%x, uinf(1), n)
293 call cfill(v%x, uinf(2), n)
294 call cfill(w%x, uinf(3), n)
302 type(field_t),
intent(inout) :: u
303 type(field_t),
intent(inout) :: v
304 type(field_t),
intent(inout) :: w
305 real(kind=rp),
intent(in) :: delta
306 real(kind=rp),
intent(in) :: uinf(3)
307 character(len=*),
intent(in) :: type
308 procedure(blasius_profile),
pointer :: bla => null()
310 character(len=LOG_SIZE) :: log_buf
312 call neko_log%message(
"Type : blasius")
313 write (log_buf,
'(A,ES12.6)')
"delta : ", delta
314 call neko_log%message(log_buf)
315 call neko_log%message(
"Approximation : " // trim(type))
316 write (log_buf,
'(A,"[",2(ES12.6,","),ES12.6,"]")')
"Value : ", &
317 uinf(1), uinf(2), uinf(3)
318 call neko_log%message(log_buf)
320 select case (trim(type))
322 bla => blasius_linear
324 bla => blasius_quadratic
328 bla => blasius_quartic
334 call neko_error(
'Invalid Blasius approximation')
337 if ((uinf(1) .gt. 0.0_rp) .and. abscmp(uinf(2), 0.0_rp) &
338 .and. abscmp(uinf(3), 0.0_rp))
then
339 do i = 1, u%dof%size()
340 u%x(i,1,1,1) = bla(u%dof%z(i,1,1,1), delta, uinf(1))
341 v%x(i,1,1,1) = 0.0_rp
342 w%x(i,1,1,1) = 0.0_rp
344 else if (abscmp(uinf(1), 0.0_rp) .and. (uinf(2) .gt. 0.0_rp) &
345 .and. abscmp(uinf(3), 0.0_rp))
then
346 do i = 1, u%dof%size()
347 u%x(i,1,1,1) = 0.0_rp
348 v%x(i,1,1,1) = bla(u%dof%x(i,1,1,1), delta, uinf(2))
349 w%x(i,1,1,1) = 0.0_rp
351 else if (abscmp(uinf(1), 0.0_rp) .and. abscmp(uinf(2), 0.0_rp) &
352 .and. (uinf(3) .gt. 0.0_rp))
then
353 do i = 1, u%dof%size()
354 u%x(i,1,1,1) = 0.0_rp
355 v%x(i,1,1,1) = 0.0_rp
356 w%x(i,1,1,1) = bla(u%dof%y(i,1,1,1), delta, uinf(3))
372 type(field_t),
intent(inout) :: u
373 type(field_t),
intent(inout) :: v
374 type(field_t),
intent(inout) :: w
375 real(kind=rp),
intent(in),
dimension(3) :: base_value
376 character(len=*),
intent(in) :: zone_name
377 real(kind=rp),
intent(in) :: zone_value(:)
378 character(len=LOG_SIZE) :: log_buf
381 class(point_zone_t),
pointer :: zone
384 call neko_log%message(
"Type : point_zone")
385 write (log_buf,
'(A,ES12.6)')
"Base value : ", base_value
386 call neko_log%message(log_buf)
387 call neko_log%message(
"Zone name : " // trim(zone_name))
388 write (log_buf,
'(A,"[",2(ES12.6,","),ES12.6," ]")')
"Value : ", &
389 zone_value(1), zone_value(2), zone_value(3)
390 call neko_log%message(log_buf)
395 zone => neko_point_zone_registry%get_point_zone(trim(zone_name))
397 call cfill_mask(u%x, zone_value(1),
size, zone%mask%get(), zone%size)
398 call cfill_mask(v%x, zone_value(2),
size, zone%mask%get(), zone%size)
399 call cfill_mask(w%x, zone_value(3),
size, zone%mask%get(), zone%size)
419 interpolate, mesh_file_name, global_interp_subdict)
420 type(field_t),
target,
intent(inout) :: u
421 type(field_t),
target,
intent(inout) :: v
422 type(field_t),
target,
intent(inout) :: w
423 type(field_t),
target,
intent(inout) :: p
424 character(len=*),
intent(inout) :: file_name
425 logical,
intent(in) :: interpolate
426 character(len=*),
intent(inout) :: mesh_file_name
427 type(json_file),
intent(inout) :: global_interp_subdict
429 type(field_t),
pointer :: us, vs, ws, ps
436 call import_fields(file_name, global_interp_subdict, mesh_file_name, &
437 u = us, v = vs, w = ws, p = ps, &
438 interpolate = interpolate)
440 nullify(us, vs, ws, ps)
444 call u%copy_from(device_to_host, sync = .false.)
445 call v%copy_from(device_to_host, sync = .false.)
446 call w%copy_from(device_to_host, sync = .false.)
447 call p%copy_from(device_to_host, sync = .true.)
Copy data between host and device (or device and device)
Synchronize a device or stream.
Abstract interface for computing a Blasius flow profile.
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...
subroutine set_flow_ic_usr(u, v, w, p, coef, gs, user_proc, scheme_name)
Set intial flow condition (user defined)
subroutine, public set_flow_ic_fld(u, v, w, p, file_name, interpolate, mesh_file_name, global_interp_subdict)
Set the initial condition of the flow based on a field. @detail The fields are read from an fld file....
subroutine set_flow_ic_point_zone(u, v, w, base_value, zone_name, zone_value)
Set the initial condition of the flow based on a point zone.
subroutine set_flow_ic_int(u, v, w, p, coef, gs, type, params)
Set initial flow condition (builtin)
subroutine set_compressible_flow_ic_usr(rho, u, v, w, p, coef, gs, user_proc, scheme_name)
Set intial flow condition (user defined) for compressible flows.
subroutine set_flow_ic_uniform(u, v, w, uinf)
Uniform initial condition.
subroutine set_flow_ic_common(u, v, w, p, coef, gs)
subroutine set_flow_ic_blasius(u, v, w, delta, uinf, type)
Set a Blasius profile as initial condition.
real(kind=rp) function, public blasius_quadratic(y, delta, u)
Quadratic approximate Blasius Profile .
real(kind=rp) function, public blasius_quartic(y, delta, u)
Quartic approximate Blasius Profile .
real(kind=rp) function, public blasius_sin(y, delta, u)
Sinusoidal approximate Blasius Profile .
real(kind=rp) function, public blasius_cubic(y, delta, u)
Cubic approximate Blasius Profile .
real(kind=rp) function, public blasius_tanh(y, delta, u)
Hyperbolic tangent approximate Blasius Profile from O. Savas (2012) where is the 99 percent thickne...
real(kind=rp) function, public blasius_linear(y, delta, u)
Linear approximate Blasius profile .
Implements global_interpolation given a dofmap.
Importation of fields from fld files.
Routines to interpolate between different spaces.
Utilities for retrieving parameters from the case files.
subroutine, public json_get_subdict_or_empty(json, key, output)
Extract a sub-object from a json object and returns an empty object if the key is missing.
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.
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.
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.