49  use json_module, 
only : json_file
 
   74    type(
field_t), 
intent(inout) :: u
 
   75    type(
field_t), 
intent(inout) :: v
 
   76    type(
field_t), 
intent(inout) :: w
 
   77    type(
field_t), 
intent(inout) :: p
 
   78    type(
coef_t), 
intent(in) :: coef
 
   79    type(
gs_t), 
intent(inout) :: gs
 
   80    character(len=*) :: type
 
   81    type(json_file), 
intent(inout) :: params
 
   82    real(kind=
rp) :: delta, tol
 
   83    real(kind=
rp), 
allocatable :: uinf(:)
 
   84    real(kind=
rp), 
allocatable :: zone_value(:)
 
   85    character(len=:), 
allocatable :: read_str
 
   86    character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
 
   87    logical :: interpolate
 
   93    if (trim(type) .eq. 
'uniform') 
then 
  101    else if (trim(type) .eq. 
'blasius') 
then 
  103       call json_get(params, 
'delta', delta)
 
  104       call json_get(params, 
'approximation', read_str)
 
  105       call json_get(params, 
'freestream_velocity', uinf)
 
  112    else if (trim(type) .eq. 
'point_zone') 
then 
  114       call json_get(params, 
'base_value', uinf)
 
  115       call json_get(params, 
'zone_name', read_str)
 
  116       call json_get(params, 
'zone_value', zone_value)
 
  123    else if (trim(type) .eq. 
'field') 
then 
  125       call json_get(params, 
'file_name', read_str)
 
  126       fname = trim(read_str)
 
  131       mesh_fname = trim(read_str)
 
 
  145    type(
field_t), 
target, 
intent(inout) :: u
 
  146    type(
field_t), 
target, 
intent(inout) :: v
 
  147    type(
field_t), 
target, 
intent(inout) :: w
 
  148    type(
field_t), 
target, 
intent(inout) :: p
 
  149    type(
coef_t), 
intent(in) :: coef
 
  150    type(
gs_t), 
intent(inout) :: gs
 
  152    character(len=*), 
intent(in) :: scheme_name
 
  160    call fields%assign_to_field(1, u)
 
  161    call fields%assign_to_field(2, v)
 
  162    call fields%assign_to_field(3, w)
 
  163    call fields%assign_to_field(4, p)
 
  165    call user_proc(scheme_name, fields)
 
 
  174       user_proc, scheme_name)
 
  175    type(
field_t), 
target, 
intent(inout) :: rho
 
  176    type(
field_t), 
target, 
intent(inout) :: u
 
  177    type(
field_t), 
target, 
intent(inout) :: v
 
  178    type(
field_t), 
target, 
intent(inout) :: w
 
  179    type(
field_t), 
target, 
intent(inout) :: p
 
  180    type(
coef_t), 
intent(in) :: coef
 
  181    type(
gs_t), 
intent(inout) :: gs
 
  183    character(len=*), 
intent(in) :: scheme_name
 
  188    call neko_log%message(
"Type: user (compressible flows)")
 
  191    call fields%assign_to_field(1, rho)
 
  192    call fields%assign_to_field(2, u)
 
  193    call fields%assign_to_field(3, v)
 
  194    call fields%assign_to_field(4, w)
 
  195    call fields%assign_to_field(5, p)
 
  196    call user_proc(scheme_name, fields)
 
  209    call gs%op(p%x, p%dof%size(), gs_op_add)
 
  210    call gs%op(rho%x, rho%dof%size(), gs_op_add)
 
  213       call device_col2(rho%x_d, coef%mult_d, rho%dof%size())
 
  216       call col2(rho%x, coef%mult, rho%dof%size())
 
  217       call col2(p%x, coef%mult, p%dof%size())
 
 
  223    type(
field_t), 
intent(inout) :: u
 
  224    type(
field_t), 
intent(inout) :: v
 
  225    type(
field_t), 
intent(inout) :: w
 
  226    type(
field_t), 
intent(inout) :: p
 
  227    type(
coef_t), 
intent(in) :: coef
 
  228    type(
gs_t), 
intent(inout) :: gs
 
  243    call gs%op(u%x, u%dof%size(), gs_op_add)
 
  244    call gs%op(v%x, v%dof%size(), gs_op_add)
 
  245    call gs%op(w%x, w%dof%size(), gs_op_add)
 
  252       call col2(u%x, coef%mult, u%dof%size())
 
  253       call col2(v%x, coef%mult, v%dof%size())
 
  254       call col2(w%x, coef%mult, w%dof%size())
 
 
  261    type(
field_t), 
intent(inout) :: u
 
  262    type(
field_t), 
intent(inout) :: v
 
  263    type(
field_t), 
intent(inout) :: w
 
  264    real(kind=
rp), 
intent(in) :: uinf(3)
 
  266    character(len=LOG_SIZE) :: log_buf
 
  268    call neko_log%message(
"Type : uniform")
 
  269    write (log_buf, 
'(A, 3(ES12.6, A))') 
"Value: [", (uinf(i), 
", ", i=1, 2), &
 
  278       call cfill(u%x, uinf(1), n)
 
  279       call cfill(v%x, uinf(2), n)
 
  280       call cfill(w%x, uinf(3), n)
 
 
  288    type(
field_t), 
intent(inout) :: u
 
  289    type(
field_t), 
intent(inout) :: v
 
  290    type(
field_t), 
intent(inout) :: w
 
  291    real(kind=
rp), 
intent(in) :: delta
 
  292    real(kind=
rp), 
intent(in) :: uinf(3)
 
  293    character(len=*), 
intent(in) :: type
 
  296    character(len=LOG_SIZE) :: log_buf
 
  298    call neko_log%message(
"Type         : blasius")
 
  299    write (log_buf, 
'(A,ES12.6)') 
"delta        : ", delta
 
  301    call neko_log%message(
"Approximation : " // trim(type))
 
  302    write (log_buf, 
'(A,"[",2(ES12.6,","),ES12.6,"]")') 
"Value         : ", &
 
  303         uinf(1), uinf(2), uinf(3)
 
  306    select case (trim(type))
 
  320       call neko_error(
'Invalid Blasius approximation')
 
  323    if ((uinf(1) .gt. 0.0_rp) .and. 
abscmp(uinf(2), 0.0_rp) &
 
  324         .and. 
abscmp(uinf(3), 0.0_rp)) 
then 
  325       do i = 1, u%dof%size()
 
  326          u%x(i,1,1,1) = bla(u%dof%z(i,1,1,1), delta, uinf(1))
 
  327          v%x(i,1,1,1) = 0.0_rp
 
  328          w%x(i,1,1,1) = 0.0_rp
 
  330    else if (
abscmp(uinf(1), 0.0_rp) .and. (uinf(2) .gt. 0.0_rp) &
 
  331         .and. 
abscmp(uinf(3), 0.0_rp)) 
then 
  332       do i = 1, u%dof%size()
 
  333          u%x(i,1,1,1) = 0.0_rp
 
  334          v%x(i,1,1,1) = bla(u%dof%x(i,1,1,1), delta, uinf(2))
 
  335          w%x(i,1,1,1) = 0.0_rp
 
  337    else if (
abscmp(uinf(1), 0.0_rp) .and. 
abscmp(uinf(2), 0.0_rp) &
 
  338         .and. (uinf(3) .gt. 0.0_rp)) 
then 
  339       do i = 1, u%dof%size()
 
  340          u%x(i,1,1,1) = 0.0_rp
 
  341          v%x(i,1,1,1) = 0.0_rp
 
  342          w%x(i,1,1,1) = bla(u%dof%y(i,1,1,1), delta, uinf(3))
 
 
  358    type(
field_t), 
intent(inout) :: u
 
  359    type(
field_t), 
intent(inout) :: v
 
  360    type(
field_t), 
intent(inout) :: w
 
  361    real(kind=
rp), 
intent(in), 
dimension(3) :: base_value
 
  362    character(len=*), 
intent(in) :: zone_name
 
  363    real(kind=
rp), 
intent(in) :: zone_value(:)
 
  364    character(len=LOG_SIZE) :: log_buf
 
  370    call neko_log%message(
"Type       : point_zone")
 
  371    write (log_buf, 
'(A,ES12.6)') 
"Base value : ", base_value
 
  373    call neko_log%message(
"Zone name : " // trim(zone_name))
 
  374    write (log_buf, 
'(A,"[",2(ES12.6,","),ES12.6," ]")') 
"Value      : ", &
 
  375         zone_value(1), zone_value(2), zone_value(3)
 
  383    call cfill_mask(u%x, zone_value(1), 
size, zone%mask%get(), zone%size)
 
  384    call cfill_mask(v%x, zone_value(2), 
size, zone%mask%get(), zone%size)
 
  385    call cfill_mask(w%x, zone_value(3), 
size, zone%mask%get(), zone%size)
 
 
  407       interpolate, tolerance, mesh_file_name)
 
  408    type(
field_t), 
intent(inout) :: u
 
  409    type(
field_t), 
intent(inout) :: v
 
  410    type(
field_t), 
intent(inout) :: w
 
  411    type(
field_t), 
intent(inout) :: p
 
  412    character(len=*), 
intent(in) :: file_name
 
  413    logical, 
intent(in) :: interpolate
 
  414    real(kind=
rp), 
intent(in) :: tolerance
 
  415    character(len=*), 
intent(inout) :: mesh_file_name
 
  417    character(len=LOG_SIZE) :: log_buf
 
  418    integer :: sample_idx, sample_mesh_idx
 
  419    integer :: last_index
 
  422    logical :: mesh_mismatch
 
  433    call neko_log%message(
"Type          : field")
 
  434    call neko_log%message(
"File name     : " // trim(file_name))
 
  435    write (log_buf, 
'(A,L1)') 
"Interpolation : ", interpolate
 
  437    if (interpolate) 
then 
  443    if (sample_idx .eq. -1) &
 
  444         call neko_error(
"Invalid file name for the initial condition. The& 
  445    & file format must be e.g. 'mean0.f00001'")
 
  451    call f%init(trim(file_name))
 
  453    if (interpolate) 
then 
  456       if (mesh_file_name .eq. 
"none") 
then 
  457          mesh_file_name = trim(file_name)
 
  458          sample_mesh_idx = sample_idx
 
  464          if (sample_mesh_idx .eq. -1) 
then 
  465             call neko_error(
"Invalid file name for the initial condition. & 
  466             &The file format must be e.g. 'mean0.f00001'")
 
  469          write (log_buf, 
'(A,ES12.6)') 
"Tolerance     : ", tolerance
 
  471          write (log_buf, 
'(A,A)') 
"Mesh file     : ", &
 
  478       if (sample_mesh_idx .ne. sample_idx) 
then 
  479          call f%set_counter(sample_mesh_idx)
 
  480          call f%read(fld_data)
 
  486    call f%set_counter(sample_idx)
 
  487    call f%read(fld_data)
 
  495    mesh_mismatch = (fld_data%glb_nelv .ne. u%msh%glb_nelv .or. &
 
  496         fld_data%gdim .ne. u%msh%gdim)
 
  498    if (mesh_mismatch .and. .not. interpolate) 
then 
  499       call neko_error(
"The fld file must match the current mesh! & 
  500       &Use 'interpolate': 'true' to enable interpolation.")
 
  501    else if (.not. mesh_mismatch .and. interpolate) 
then 
  502       call neko_log%warning(
"You have activated interpolation but you might & 
  503       &still be using the same mesh.")
 
  507    if (interpolate) 
then 
  510       select type (ft => f%file_type)
 
  512          if (.not. ft%dp_precision) 
then 
  513             call neko_warning(
"The coordinates read from the field file are & 
  514             &in single precision.")
 
  515             call neko_log%message(
"It is recommended to use a mesh in double & 
  516             &precision for better interpolation results.")
 
  517             call neko_log%message(
"If the interpolation does not work, you& 
  518             &can try to increase the tolerance.")
 
  525          call device_memcpy(fld_data%x%x, fld_data%x%x_d, fld_data%x%size(),&
 
  527          call device_memcpy(fld_data%y%x, fld_data%y%x_d, fld_data%y%size(),&
 
  529          call device_memcpy(fld_data%z%x, fld_data%z%x_d, fld_data%z%size(),&
 
  534       call fld_data%generate_interpolator(global_interp, u%dof, u%msh, &
 
  542       call global_interp%evaluate(u%x(:,1,1,1), fld_data%u%x, .false.)
 
  543       call global_interp%evaluate(v%x(:,1,1,1), fld_data%v%x, .false.)
 
  544       call global_interp%evaluate(w%x(:,1,1,1), fld_data%w%x, .false.)
 
  545       call global_interp%evaluate(p%x(:,1,1,1), fld_data%p%x, .false.)
 
  551       call global_interp%free
 
  556       call prev_xh%init(
gll, fld_data%lx, fld_data%ly, fld_data%lz)
 
  557       call space_interp%init(u%Xh, prev_xh)
 
  560       call space_interp%map_host(u%x, fld_data%u%x, fld_data%nelv, u%Xh)
 
  561       call space_interp%map_host(v%x, fld_data%v%x, fld_data%nelv, u%Xh)
 
  562       call space_interp%map_host(w%x, fld_data%w%x, fld_data%nelv, u%Xh)
 
  563       call space_interp%map_host(p%x, fld_data%p%x, fld_data%nelv, u%Xh)
 
  565       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_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_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.
 
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.