38  use, 
intrinsic :: iso_c_binding, only : c_ptr
 
   43  use json_module, 
only : json_file
 
   72     procedure, pass(this) :: init_from_components => &
 
   74     procedure, pass(this) :: set_stress_scalar => &
 
   76     procedure, pass(this) :: set_stress_array => &
 
   79     generic :: set_stress => set_stress_scalar, set_stress_array
 
 
   91    integer, 
intent(in) :: n
 
   92    real(kind=
rp), 
intent(inout), 
dimension(n) :: x
 
   94    logical, 
intent(in), 
optional :: strong
 
   95    integer :: i, m, k, facet
 
   99    call neko_error(
"The shear stress bc is not applicable to scalar fields.")
 
 
  107    integer, 
intent(in) :: n
 
  108    real(kind=
rp), 
intent(inout), 
dimension(n) :: x
 
  109    real(kind=
rp), 
intent(inout), 
dimension(n) :: y
 
  110    real(kind=
rp), 
intent(inout), 
dimension(n) :: z
 
  112    logical, 
intent(in), 
optional :: strong
 
  115    if (
present(strong)) 
then 
  122       call this%symmetry%apply_vector(x, y, z, n, strong = .true.)
 
  124       call this%neumann_x%apply_scalar(x, n, strong = .false.)
 
  125       call this%neumann_y%apply_scalar(y, n, strong = .false.)
 
  126       call this%neumann_z%apply_scalar(z, n, strong = .false.)
 
 
  135    type(c_ptr), 
intent(inout) :: x_d
 
  137    logical, 
intent(in), 
optional :: strong
 
  138    type(c_ptr), 
intent(inout) :: strm
 
  140    call neko_error(
"The shear stress bc is not applicable to scalar fields.")
 
 
  149    type(c_ptr), 
intent(inout) :: x_d
 
  150    type(c_ptr), 
intent(inout) :: y_d
 
  151    type(c_ptr), 
intent(inout) :: z_d
 
  153    logical, 
intent(in), 
optional :: strong
 
  154    type(c_ptr), 
intent(inout) :: strm
 
  157    if (
present(strong)) 
then 
  164       call this%symmetry%apply_vector_dev(x_d, y_d, z_d, strong = .true., &
 
  167       call this%neumann_x%apply_scalar_dev(x_d, strong = .false., strm = strm)
 
  168       call this%neumann_y%apply_scalar_dev(y_d, strong = .false., strm = strm)
 
  169       call this%neumann_z%apply_scalar_dev(z_d, strong = .false., strm = strm)
 
 
  179    type(
coef_t), 
target, 
intent(in) :: coef
 
  180    type(json_file), 
intent(inout) ::json
 
  181    real(kind=
rp), 
allocatable :: value(:)
 
  185    if (
size(
value) .ne. 3) 
then 
  186       call neko_error (
"The shear stress vector provided for the shear stress & 
  187       & boundary condition should have 3 components.")
 
  190    call this%init_from_components(coef, 
value)
 
 
  198    type(
coef_t), 
intent(in) :: coef
 
  199    real(kind=
rp), 
intent(in) :: value(3)
 
  201    call this%init_base(coef)
 
  202    this%strong = .false.
 
  204    call this%symmetry%free()
 
  205    call this%symmetry%init_from_components(this%coef)
 
  207    call this%neumann_x%free()
 
  208    call this%neumann_y%free()
 
  209    call this%neumann_z%free()
 
  211    call this%neumann_x%init_from_components(this%coef, value(1))
 
  212    call this%neumann_y%init_from_components(this%coef, value(2))
 
  213    call this%neumann_z%init_from_components(this%coef, value(3))
 
 
  219    logical, 
optional, 
intent(in) :: only_facets
 
  220    logical :: only_facets_
 
  222    if (
present(only_facets)) 
then 
  223       only_facets_ = only_facets
 
  225       only_facets_ = .false.
 
  228    call this%finalize_base(only_facets_)
 
  230    call this%symmetry%mark_facets(this%marked_facet)
 
  231    call this%symmetry%finalize()
 
  234    call this%neumann_x%mark_facets(this%marked_facet)
 
  235    call this%neumann_y%mark_facets(this%marked_facet)
 
  236    call this%neumann_z%mark_facets(this%marked_facet)
 
  238    call this%neumann_x%finalize(only_facets_)
 
  239    call this%neumann_y%finalize(only_facets_)
 
  240    call this%neumann_z%finalize(only_facets_)
 
 
  247    real(kind=
rp), 
intent(in) :: tau_x
 
  248    real(kind=
rp), 
intent(in) :: tau_y
 
  249    real(kind=
rp), 
intent(in) :: tau_z
 
  252    call this%neumann_x%set_flux(tau_x)
 
  253    call this%neumann_y%set_flux(tau_y)
 
  254    call this%neumann_z%set_flux(tau_z)
 
 
  269    call this%neumann_x%set_flux(tau_x)
 
  270    call this%neumann_y%set_flux(tau_y)
 
  271    call this%neumann_z%set_flux(tau_z)
 
 
  279    call this%symmetry%free
 
  281    call this%neumann_x%free
 
  282    call this%neumann_y%free
 
  283    call this%neumann_z%free
 
 
Retrieves a parameter by name or throws an error.
 
Defines a boundary condition.
 
Utilities for retrieving parameters from the case files.
 
Defines a Neumann boundary condition.
 
integer, parameter, public rp
Global precision used in computations.
 
Defines a shear stress boundary condition for a vector field. Maintainer: Timofey Mukha.
 
subroutine shear_stress_init_from_components(this, coef, value)
Constructor from components.
 
subroutine shear_stress_free(this)
Destructor.
 
subroutine shear_stress_apply_scalar(this, x, n, time, strong)
Apply shear stress for a scalar field x.
 
subroutine shear_stress_set_stress_array(this, tau_x, tau_y, tau_z)
Set the shear stress components.
 
subroutine shear_stress_apply_scalar_dev(this, x_d, time, strong, strm)
Boundary condition apply for a generic shear_stress condition to a vector x (device version)
 
subroutine shear_stress_set_stress_scalar(this, tau_x, tau_y, tau_z)
Set the value of the shear stress vector using 3 scalars.
 
subroutine shear_stress_finalize(this, only_facets)
 
subroutine shear_stress_init(this, coef, json)
Constructor.
 
subroutine shear_stress_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
Boundary condition apply for a generic shear_stress condition to vectors x, y and z (device version)
 
subroutine shear_stress_apply_vector(this, x, y, z, n, time, strong)
Boundary condition apply for a generic shear_stress condition to vectors x, y and z.
 
Mixed Dirichlet-Neumann axis aligned symmetry plane.
 
Module with things related to the simulation time.
 
Base type for a boundary condition.
 
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
 
A Neumann boundary condition for scalar fields. Sets the flux of the field to the chosen value.
 
A shear stress boundary condition.
 
Mixed Dirichlet-Neumann symmetry plane condition.
 
A struct that contains all info about the time, expand as needed.