49 use json_module,
only : json_file
75 type(
field_t),
intent(inout) :: u
76 type(
field_t),
intent(inout) :: v
77 type(
field_t),
intent(inout) :: w
78 type(
field_t),
intent(inout) :: p
79 type(
coef_t),
intent(in) :: coef
80 type(
gs_t),
intent(inout) :: gs
81 character(len=*) :: type
82 type(json_file),
intent(inout) :: params
83 real(kind=
rp) :: delta, tol
84 real(kind=
rp),
allocatable :: uinf(:)
85 real(kind=
rp),
allocatable :: zone_value(:)
86 character(len=:),
allocatable :: read_str
87 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
88 logical :: interpolate
94 if (trim(type) .eq.
'uniform')
then
102 else if (trim(type) .eq.
'blasius')
then
104 call json_get(params,
'delta', delta)
105 call json_get(params,
'approximation', read_str)
106 call json_get(params,
'freestream_velocity', uinf)
113 else if (trim(type) .eq.
'point_zone')
then
115 call json_get(params,
'base_value', uinf)
116 call json_get(params,
'zone_name', read_str)
117 call json_get(params,
'zone_value', zone_value)
124 else if (trim(type) .eq.
'field')
then
126 call json_get(params,
'file_name', read_str)
127 fname = trim(read_str)
132 mesh_fname = trim(read_str)
146 type(
field_t),
target,
intent(inout) :: u
147 type(
field_t),
target,
intent(inout) :: v
148 type(
field_t),
target,
intent(inout) :: w
149 type(
field_t),
target,
intent(inout) :: p
150 type(
coef_t),
intent(in) :: coef
151 type(
gs_t),
intent(inout) :: gs
153 character(len=*),
intent(in) :: scheme_name
161 call fields%assign_to_field(1, u)
162 call fields%assign_to_field(2, v)
163 call fields%assign_to_field(3, w)
164 call fields%assign_to_field(4, p)
166 call user_proc(scheme_name, fields)
175 user_proc, scheme_name)
176 type(
field_t),
target,
intent(inout) :: rho
177 type(
field_t),
target,
intent(inout) :: u
178 type(
field_t),
target,
intent(inout) :: v
179 type(
field_t),
target,
intent(inout) :: w
180 type(
field_t),
target,
intent(inout) :: p
181 type(
coef_t),
intent(in) :: coef
182 type(
gs_t),
intent(inout) :: gs
184 character(len=*),
intent(in) :: scheme_name
189 call neko_log%message(
"Type: user (compressible flows)")
192 call fields%assign_to_field(1, rho)
193 call fields%assign_to_field(2, u)
194 call fields%assign_to_field(3, v)
195 call fields%assign_to_field(4, w)
196 call fields%assign_to_field(5, p)
197 call user_proc(scheme_name, fields)
210 call gs%op(p%x, p%dof%size(), gs_op_add)
211 call gs%op(rho%x, rho%dof%size(), gs_op_add)
214 call device_col2(rho%x_d, coef%mult_d, rho%dof%size())
217 call col2(rho%x, coef%mult, rho%dof%size())
218 call col2(p%x, coef%mult, p%dof%size())
224 type(
field_t),
intent(inout) :: u
225 type(
field_t),
intent(inout) :: v
226 type(
field_t),
intent(inout) :: w
227 type(
field_t),
intent(inout) :: p
228 type(
coef_t),
intent(in) :: coef
229 type(
gs_t),
intent(inout) :: gs
245 call gs%op(u%x, u%dof%size(), gs_op_add)
246 call gs%op(v%x, v%dof%size(), gs_op_add)
247 call gs%op(w%x, w%dof%size(), gs_op_add)
255 call col2(u%x, coef%mult, u%dof%size())
256 call col2(v%x, coef%mult, v%dof%size())
257 call col2(w%x, coef%mult, w%dof%size())
264 type(
field_t),
intent(inout) :: u
265 type(
field_t),
intent(inout) :: v
266 type(
field_t),
intent(inout) :: w
267 real(kind=
rp),
intent(in) :: uinf(3)
269 character(len=LOG_SIZE) :: log_buf
271 call neko_log%message(
"Type : uniform")
272 write (log_buf,
'(A, 3(ES12.6, A))')
"Value: [", (uinf(i),
", ", i=1, 2), &
281 call cfill(u%x, uinf(1), n)
282 call cfill(v%x, uinf(2), n)
283 call cfill(w%x, uinf(3), n)
291 type(
field_t),
intent(inout) :: u
292 type(
field_t),
intent(inout) :: v
293 type(
field_t),
intent(inout) :: w
294 real(kind=
rp),
intent(in) :: delta
295 real(kind=
rp),
intent(in) :: uinf(3)
296 character(len=*),
intent(in) :: type
299 character(len=LOG_SIZE) :: log_buf
301 call neko_log%message(
"Type : blasius")
302 write (log_buf,
'(A,ES12.6)')
"delta : ", delta
304 call neko_log%message(
"Approximation : " // trim(type))
305 write (log_buf,
'(A,"[",2(ES12.6,","),ES12.6,"]")')
"Value : ", &
306 uinf(1), uinf(2), uinf(3)
309 select case (trim(type))
323 call neko_error(
'Invalid Blasius approximation')
326 if ((uinf(1) .gt. 0.0_rp) .and.
abscmp(uinf(2), 0.0_rp) &
327 .and.
abscmp(uinf(3), 0.0_rp))
then
328 do i = 1, u%dof%size()
329 u%x(i,1,1,1) = bla(u%dof%z(i,1,1,1), delta, uinf(1))
330 v%x(i,1,1,1) = 0.0_rp
331 w%x(i,1,1,1) = 0.0_rp
333 else if (
abscmp(uinf(1), 0.0_rp) .and. (uinf(2) .gt. 0.0_rp) &
334 .and.
abscmp(uinf(3), 0.0_rp))
then
335 do i = 1, u%dof%size()
336 u%x(i,1,1,1) = 0.0_rp
337 v%x(i,1,1,1) = bla(u%dof%x(i,1,1,1), delta, uinf(2))
338 w%x(i,1,1,1) = 0.0_rp
340 else if (
abscmp(uinf(1), 0.0_rp) .and.
abscmp(uinf(2), 0.0_rp) &
341 .and. (uinf(3) .gt. 0.0_rp))
then
342 do i = 1, u%dof%size()
343 u%x(i,1,1,1) = 0.0_rp
344 v%x(i,1,1,1) = 0.0_rp
345 w%x(i,1,1,1) = bla(u%dof%y(i,1,1,1), delta, uinf(3))
361 type(
field_t),
intent(inout) :: u
362 type(
field_t),
intent(inout) :: v
363 type(
field_t),
intent(inout) :: w
364 real(kind=
rp),
intent(in),
dimension(3) :: base_value
365 character(len=*),
intent(in) :: zone_name
366 real(kind=
rp),
intent(in) :: zone_value(:)
367 character(len=LOG_SIZE) :: log_buf
373 call neko_log%message(
"Type : point_zone")
374 write (log_buf,
'(A,ES12.6)')
"Base value : ", base_value
376 call neko_log%message(
"Zone name : " // trim(zone_name))
377 write (log_buf,
'(A,"[",2(ES12.6,","),ES12.6," ]")')
"Value : ", &
378 zone_value(1), zone_value(2), zone_value(3)
386 call cfill_mask(u%x, zone_value(1),
size, zone%mask%get(), zone%size)
387 call cfill_mask(v%x, zone_value(2),
size, zone%mask%get(), zone%size)
388 call cfill_mask(w%x, zone_value(3),
size, zone%mask%get(), zone%size)
410 interpolate, tolerance, mesh_file_name)
411 type(
field_t),
intent(inout) :: u
412 type(
field_t),
intent(inout) :: v
413 type(
field_t),
intent(inout) :: w
414 type(
field_t),
intent(inout) :: p
415 character(len=*),
intent(in) :: file_name
416 logical,
intent(in) :: interpolate
417 real(kind=
rp),
intent(in) :: tolerance
418 character(len=*),
intent(inout) :: mesh_file_name
420 character(len=LOG_SIZE) :: log_buf
421 integer :: sample_idx, sample_mesh_idx
422 integer :: last_index
425 logical :: mesh_mismatch
436 call neko_log%message(
"Type : field")
437 call neko_log%message(
"File name : " // trim(file_name))
438 write (log_buf,
'(A,L1)')
"Interpolation : ", interpolate
440 if (interpolate)
then
446 if (sample_idx .eq. -1) &
447 call neko_error(
"Invalid file name for the initial condition. The&
448 & file format must be e.g. 'mean0.f00001'")
454 call f%init(trim(file_name))
456 if (interpolate)
then
459 if (mesh_file_name .eq.
"none")
then
460 mesh_file_name = trim(file_name)
461 sample_mesh_idx = sample_idx
467 if (sample_mesh_idx .eq. -1)
then
468 call neko_error(
"Invalid file name for the initial condition. &
469 &The file format must be e.g. 'mean0.f00001'")
472 write (log_buf,
'(A,ES12.6)')
"Tolerance : ", tolerance
474 write (log_buf,
'(A,A)')
"Mesh file : ", &
481 if (sample_mesh_idx .ne. sample_idx)
then
482 call f%set_counter(sample_mesh_idx)
483 call f%read(fld_data)
489 call f%set_counter(sample_idx)
490 call f%read(fld_data)
498 mesh_mismatch = (fld_data%glb_nelv .ne. u%msh%glb_nelv .or. &
499 fld_data%gdim .ne. u%msh%gdim)
501 if (mesh_mismatch .and. .not. interpolate)
then
502 call neko_error(
"The fld file must match the current mesh! &
503 &Use 'interpolate': 'true' to enable interpolation.")
504 else if (.not. mesh_mismatch .and. interpolate)
then
505 call neko_log%warning(
"You have activated interpolation but you might &
506 &still be using the same mesh.")
510 if (interpolate)
then
513 select type (ft => f%file_type)
515 if (.not. ft%dp_precision)
then
516 call neko_warning(
"The coordinates read from the field file are &
517 &in single precision.")
518 call neko_log%message(
"It is recommended to use a mesh in double &
519 &precision for better interpolation results.")
520 call neko_log%message(
"If the interpolation does not work, you&
521 &can try to increase the tolerance.")
528 call device_memcpy(fld_data%x%x, fld_data%x%x_d, fld_data%x%size(),&
530 call device_memcpy(fld_data%y%x, fld_data%y%x_d, fld_data%y%size(),&
532 call device_memcpy(fld_data%z%x, fld_data%z%x_d, fld_data%z%size(),&
537 call fld_data%generate_interpolator(global_interp, u%dof, u%msh, &
545 call global_interp%evaluate(u%x(:,1,1,1), fld_data%u%x, .false.)
546 call global_interp%evaluate(v%x(:,1,1,1), fld_data%v%x, .false.)
547 call global_interp%evaluate(w%x(:,1,1,1), fld_data%w%x, .false.)
548 call global_interp%evaluate(p%x(:,1,1,1), fld_data%p%x, .false.)
554 call global_interp%free
559 call prev_xh%init(
gll, fld_data%lx, fld_data%ly, fld_data%lz)
560 call space_interp%init(u%Xh, prev_xh)
563 call space_interp%map_host(u%x, fld_data%u%x, fld_data%nelv, u%Xh)
564 call space_interp%map_host(v%x, fld_data%v%x, fld_data%nelv, u%Xh)
565 call space_interp%map_host(w%x, fld_data%w%x, fld_data%nelv, u%Xh)
566 call space_interp%map_host(p%x, fld_data%p%x, fld_data%nelv, u%Xh)
568 call space_interp%free
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 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, public set_flow_ic_fld(u, v, w, p, file_name, interpolate, tolerance, mesh_file_name)
Set the initial condition of the flow based on a field. @detail The fields are read from an fld file....
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.
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.
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.