37  use, 
intrinsic :: iso_c_binding, only : c_ptr
 
   43  use json_module, 
only : json_file
 
   63     procedure, pass(this) :: apply_scalar_dev => &
 
   65     procedure, pass(this) :: apply_vector_dev => &
 
 
   74    integer, 
intent(in) :: n
 
   75    real(kind=
rp), 
intent(inout), 
dimension(n) :: x
 
   77    logical, 
intent(in), 
optional :: strong
 
   79    call neko_error(
"The wall model bc is not applicable to scalar fields.")
 
 
   92    integer, 
intent(in) :: n
 
   93    real(kind=
rp), 
intent(inout), 
dimension(n) :: x
 
   94    real(kind=
rp), 
intent(inout), 
dimension(n) :: y
 
   95    real(kind=
rp), 
intent(inout), 
dimension(n) :: z
 
   97    logical, 
intent(in), 
optional :: strong
 
   98    integer :: i, m, k, fid
 
   99    real(kind=
rp) :: magtau
 
  102    if (
present(strong)) 
then 
  108    if (
present(time) .eqv. .false.) 
then 
  109       call neko_error(
"wall_model_bc_apply_vector: time is required.")
 
  112    if (.not. strong_) 
then 
  114       call this%wall_model%compute(time%t, time%tstep)
 
  117       call this%wall_model%compute_mag_field()
 
  121       call this%set_stress(this%wall_model%tau_x, this%wall_model%tau_y, &
 
  122            this%wall_model%tau_z)
 
  126    call this%shear_stress_t%apply_vector(x, y, z, n, time, strong_)
 
 
  134    type(c_ptr), 
intent(inout) :: x_d
 
  136    logical, 
intent(in), 
optional :: strong
 
  137    type(c_ptr), 
intent(inout) :: strm
 
  139    call neko_error(
"The wall model bc is not applicable to scalar fields.")
 
 
  148    type(c_ptr), 
intent(inout) :: x_d
 
  149    type(c_ptr), 
intent(inout) :: y_d
 
  150    type(c_ptr), 
intent(inout) :: z_d
 
  152    logical, 
intent(in), 
optional :: strong
 
  154    type(c_ptr), 
intent(inout) :: strm
 
  156    if (
present(strong)) 
then 
  162    if (
present(time) .eqv. .false.) 
then 
  163       call neko_error(
"wall_model_bc_apply_vector_dev: time is required.")
 
  166    if (.not. strong_) 
then 
  168       call this%wall_model%compute(time%t, time%tstep)
 
  171       call this%wall_model%compute_mag_field()
 
  175       call this%set_stress(this%wall_model%tau_x, this%wall_model%tau_y, &
 
  176            this%wall_model%tau_z)
 
  180    call this%shear_stress_t%apply_vector_dev(x_d, y_d, z_d, &
 
 
  190    type(
coef_t), 
target, 
intent(in) :: coef
 
  191    type(json_file), 
intent(inout) :: json
 
  192    real(kind=
rp) :: value(3) = [0, 0, 0]
 
  193    character(len=:), 
allocatable :: type_name
 
  196    call this%shear_stress_t%init_from_components(coef, 
value)
 
  198    call json_get(json, 
"model", type_name)
 
  199    call wall_model_allocator(this%wall_model, type_name)
 
  201    call this%wall_model%partial_init(coef, json)
 
 
  208    call this%shear_stress_t%free()
 
  209    call this%wall_model%free()
 
 
  216    logical, 
optional, 
intent(in) :: only_facets
 
  218    if (
present(only_facets)) 
then 
  219       if (only_facets .eqv. .false.) 
then 
  220          call neko_error(
"For wall_model_bc_t, only_facets has to be true.")
 
  224    call this%shear_stress_t%finalize(.true.)
 
  225    call this%wall_model%finalize(this%msk, this%facet)
 
 
Retrieves a parameter by name or throws an error.
 
Utilities for retrieving parameters from the case files.
 
integer, parameter, public rp
Global precision used in computations.
 
Defines a shear stress boundary condition for a vector field. Maintainer: Timofey Mukha.
 
Module with things related to the simulation time.
 
Defines the wall_model_bc_t type. Maintainer: Timofey Mukha.
 
subroutine wall_model_bc_init(this, coef, json)
Constructor.
 
subroutine wall_model_bc_apply_scalar(this, x, n, time, strong)
Apply shear stress for a scalar field x.
 
subroutine wall_model_bc_apply_scalar_dev(this, x_d, time, strong, strm)
Boundary condition apply for a generic wall_model_bc condition to a vector x (device version)
 
subroutine wall_model_bc_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
Boundary condition apply for a generic wall_model_bc condition to vectors x, y and z (device version)
 
subroutine wall_model_bc_apply_vector(this, x, y, z, n, time, strong)
Apply the boundary condition to the right-hand side.
 
subroutine wall_model_bc_free(this)
Destructor.
 
subroutine wall_model_bc_finalize(this, only_facets)
Finalize by building mask arrays and init'ing the wall model.
 
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
 
A shear stress boundary condition.
 
A struct that contains all info about the time, expand as needed.
 
Base abstract type for wall-stress models for wall-modelled LES.
 
A shear stress boundary condition, computing the stress values using a wall model.