37  use, 
intrinsic :: iso_c_binding, only : c_ptr
 
   49     real(kind=
rp), 
allocatable, 
private :: flux_(:)
 
   63     generic :: set_flux => set_flux_scalar, set_flux_array
 
 
   74    integer, 
intent(in) :: n
 
   75    real(kind=
rp), 
intent(inout),  
dimension(n) :: x
 
   76    real(kind=
rp), 
intent(in), 
optional :: t
 
   77    integer, 
intent(in), 
optional :: tstep
 
   78    integer :: i, m, k, facet
 
   90          x(k) = x(k) + this%flux_(i)*this%coef%area(idx(2), idx(3), facet, &
 
   93          x(k) = x(k) + this%flux_(i)*this%coef%area(idx(1), idx(3), facet, &
 
   96          x(k) = x(k) + this%flux_(i)*this%coef%area(idx(1), idx(2), facet, &
 
 
  106    integer, 
intent(in) :: n
 
  107    real(kind=
rp), 
intent(inout),  
dimension(n) :: x
 
  108    real(kind=
rp), 
intent(inout),  
dimension(n) :: y
 
  109    real(kind=
rp), 
intent(inout),  
dimension(n) :: z
 
  110    real(kind=
rp), 
intent(in), 
optional :: t
 
  111    integer, 
intent(in), 
optional :: tstep
 
  113    call neko_error(
"Neumann bc not implemented for vectors")
 
 
  120    class(
neumann_t), 
intent(inout), 
target :: this
 
  122    real(kind=
rp), 
intent(in), 
optional :: t
 
  123    integer, 
intent(in), 
optional :: tstep
 
  125    call neko_error(
"Neumann bc not implemented on the device")
 
 
  132    class(
neumann_t), 
intent(inout), 
target :: this
 
  136    real(kind=
rp), 
intent(in), 
optional :: t
 
  137    integer, 
intent(in), 
optional :: tstep
 
  139    call neko_error(
"Neumann bc not implemented on the device")
 
 
  145    class(
neumann_t), 
target, 
intent(inout) :: this
 
  147    call this%free_base()
 
 
  154    real(kind=
rp) :: flux(this%msk(0))
 
 
  162    real(kind=
rp), 
intent(in) :: flux
 
 
  171    real(kind=
rp), 
intent(in) :: flux(this%msk(0))
 
  173    call copy(this%flux_, flux, this%msk(0))
 
 
  180    real(kind=
rp), 
intent(in) :: flux
 
  183    allocate(this%flux_(this%msk(0)))
 
  185    call cfill(this%flux_, flux, this%msk(0))
 
 
__device__ void nonlinear_index(const int idx, const int lx, int *index)
 
Defines a boundary condition.
 
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
 
subroutine, public copy(a, b, n)
Copy a vector .
 
Defines a Neumann boundary condition.
 
subroutine neumann_set_flux_array(this, flux)
Set the flux using an array of values.
 
subroutine neumann_apply_vector(this, x, y, z, n, t, tstep)
Boundary condition apply for a generic Neumann condition to vectors x, y and z.
 
subroutine neumann_apply_scalar(this, x, n, t, tstep)
Boundary condition apply for a generic Neumann condition to a vector x.
 
subroutine neumann_apply_scalar_dev(this, x_d, t, tstep)
Boundary condition apply for a generic Neumann condition to a vector x (device version)
 
subroutine neumann_set_flux_scalar(this, flux)
Set the flux using a scalar.
 
pure real(kind=rp) function, dimension(this%msk(0)) neumann_flux(this)
Get the flux.
 
subroutine neumann_finalize_neumann(this, flux)
Finalize by setting the flux.
 
subroutine neumann_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
Boundary condition apply for a generic Neumann condition to vectors x, y and z (device version)
 
subroutine neumann_free(this)
Destructor.
 
integer, parameter, public rp
Global precision used in computations.
 
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.