38  use, 
intrinsic :: iso_c_binding, only : c_ptr
 
   65     procedure, pass(this) :: init_shear_stress => &
 
   67     procedure, pass(this) :: set_stress_scalar => &
 
   69     procedure, pass(this) :: set_stress_array => &
 
   71     generic :: set_stress => set_stress_scalar, set_stress_array
 
 
   80    integer, 
intent(in) :: n
 
   81    real(kind=
rp), 
intent(inout),  
dimension(n) :: x
 
   82    real(kind=
rp), 
intent(in), 
optional :: t
 
   83    integer, 
intent(in), 
optional :: tstep
 
   84    integer :: i, m, k, facet
 
   88    call neko_error(
"The shear stress bc is not applicable to scalar fields.")
 
 
   96    integer, 
intent(in) :: n
 
   97    real(kind=
rp), 
intent(inout),  
dimension(n) :: x
 
   98    real(kind=
rp), 
intent(inout),  
dimension(n) :: y
 
   99    real(kind=
rp), 
intent(inout),  
dimension(n) :: z
 
  100    real(kind=
rp), 
intent(in), 
optional :: t
 
  101    integer, 
intent(in), 
optional :: tstep
 
  103    call this%neumann_x%apply_scalar(x, n, t, tstep)
 
  104    call this%neumann_y%apply_scalar(y, n, t, tstep)
 
  105    call this%neumann_z%apply_scalar(z, n, t, tstep)
 
 
  114    real(kind=
rp), 
intent(in), 
optional :: t
 
  115    integer, 
intent(in), 
optional :: tstep
 
  117    call neko_error(
"shear_stress bc not implemented on the device")
 
 
  128    real(kind=
rp), 
intent(in), 
optional :: t
 
  129    integer, 
intent(in), 
optional :: tstep
 
  131    call neko_error(
"shear_stress bc not implemented on the device")
 
 
  140    type(
coef_t), 
target, 
intent(in) :: coef
 
  144    call this%symmetry%free()
 
  145    call this%symmetry%init_base(coef)
 
  146    call this%symmetry%mark_facets(this%marked_facet)
 
  147    call this%symmetry%finalize()
 
  148    call this%symmetry%init(coef)
 
  150    call this%neumann_x%init_base(coef)
 
  151    call this%neumann_y%init_base(coef)
 
  152    call this%neumann_z%init_base(coef)
 
  154    call this%neumann_x%mark_facets(this%marked_facet)
 
  155    call this%neumann_y%mark_facets(this%marked_facet)
 
  156    call this%neumann_z%mark_facets(this%marked_facet)
 
  158    call this%neumann_x%finalize_neumann(0.0_rp)
 
  159    call this%neumann_y%finalize_neumann(0.0_rp)
 
  160    call this%neumann_z%finalize_neumann(0.0_rp)
 
 
  168    real(kind=
rp), 
intent(in) :: tau_x
 
  169    real(kind=
rp), 
intent(in) :: tau_y
 
  170    real(kind=
rp), 
intent(in) :: tau_z
 
  173    call this%neumann_x%set_flux(tau_x)
 
  174    call this%neumann_y%set_flux(tau_y)
 
  175    call this%neumann_z%set_flux(tau_z)
 
 
  183    real(kind=
rp), 
intent(in) :: tau_x(this%msk(0))
 
  184    real(kind=
rp), 
intent(in) :: tau_y(this%msk(0))
 
  185    real(kind=
rp), 
intent(in) :: tau_z(this%msk(0))
 
  187    call this%neumann_x%set_flux(tau_x)
 
  188    call this%neumann_y%set_flux(tau_y)
 
  189    call this%neumann_z%set_flux(tau_z)
 
 
  197    call this%symmetry%free
 
  199    call this%neumann_x%free
 
  200    call this%neumann_y%free
 
  201    call this%neumann_z%free
 
 
Defines a boundary condition.
 
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_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
Boundary condition apply for a generic shear_stress condition to vectors x, y and z (device version)
 
subroutine shear_stress_free(this)
Destructor.
 
subroutine shear_stress_init_shear_stress(this, coef)
Additional constructor that should be run after the finalization of the bc. Similar to the symmetry c...
 
subroutine shear_stress_set_stress_array(this, tau_x, tau_y, tau_z)
Set the shear stress components.
 
subroutine shear_stress_apply_scalar(this, x, n, t, tstep)
Apply shear stress for a scalar field x.
 
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_apply_vector(this, x, y, z, n, t, tstep)
Boundary condition apply for a generic shear_stress condition to vectors x, y and z.
 
subroutine shear_stress_apply_scalar_dev(this, x_d, t, tstep)
Boundary condition apply for a generic shear_stress condition to a vector x (device version)
 
Mixed Dirichlet-Neumann axis aligned symmetry plane.
 
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.