37 use,
intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
40 use json_module,
only : json_file
68 real(kind=
rp),
allocatable,
private :: init_flux_(:)
71 logical :: uniform_0 = .false.
83 generic :: init_from_components => &
91 generic :: set_flux => set_flux_scalar, set_flux_array
104 class(
neumann_t),
intent(inout),
target :: this
105 type(
coef_t),
target,
intent(in) :: coef
106 type(json_file),
intent(inout) :: json
107 real(kind=
rp) :: flux
110 call this%init_base(coef)
111 this%strong = .false.
114 call json%get(
"flux", this%init_flux_, found)
117 if (.not. found)
then
119 allocate(this%init_flux_(1))
120 this%init_flux_(1) = flux
123 if ((
size(this%init_flux_) .ne. 1) &
124 .and. (
size(this%init_flux_) .ne. 3))
then
125 call neko_error(
"Neumann BC flux must be a scalar or a 3-component" // &
129 allocate(this%flux(
size(this%init_flux_)))
136 class(
neumann_t),
intent(inout),
target :: this
137 type(
coef_t),
intent(in) :: coef
138 real(kind=
rp),
intent(in) :: flux(3)
140 call this%init_base(coef)
141 this%init_flux_ = flux
143 if ((
size(this%init_flux_) .ne. 3))
then
144 call neko_error(
"Neumann BC flux must be a scalar or a 3-component" // &
147 allocate(this%flux(
size(this%init_flux_)))
154 class(
neumann_t),
intent(inout),
target :: this
155 type(
coef_t),
intent(in) :: coef
156 real(kind=
rp),
intent(in) :: flux
158 call this%init_base(coef)
159 allocate(this%init_flux_(1))
160 this%init_flux_(1) = flux
161 allocate(this%flux(
size(this%init_flux_)))
168 integer,
intent(in) :: n
169 real(kind=
rp),
intent(inout),
dimension(n) :: x
171 logical,
intent(in),
optional :: strong
172 integer :: i, m, k, facet
177 if (
present(strong))
then
184 if (.not. strong_)
then
187 facet = this%facet(i)
193 this%flux(1)%x(i) * &
194 this%coef%area(idx(2), idx(3), facet, idx(4))
197 this%flux(1)%x(i) * &
198 this%coef%area(idx(1), idx(3), facet, idx(4))
201 this%flux(1)%x(i) * &
202 this%coef%area(idx(1), idx(2), facet, idx(4))
212 integer,
intent(in) :: n
213 real(kind=
rp),
intent(inout),
dimension(n) :: x
214 real(kind=
rp),
intent(inout),
dimension(n) :: y
215 real(kind=
rp),
intent(inout),
dimension(n) :: z
217 logical,
intent(in),
optional :: strong
218 integer :: i, m, k, facet
223 if (
present(strong))
then
230 if (.not. strong_)
then
233 facet = this%facet(i)
239 this%flux(1)%x(i) * &
240 this%coef%area(idx(2), idx(3), facet, idx(4))
242 this%flux(2)%x(i) * &
243 this%coef%area(idx(2), idx(3), facet, idx(4))
245 this%flux(3)%x(i) * &
246 this%coef%area(idx(2), idx(3), facet, idx(4))
249 this%flux(1)%x(i) * &
250 this%coef%area(idx(1), idx(3), facet, idx(4))
252 this%flux(2)%x(i) * &
253 this%coef%area(idx(1), idx(3), facet, idx(4))
255 this%flux(3)%x(i) * &
256 this%coef%area(idx(1), idx(3), facet, idx(4))
259 this%flux(1)%x(i) * &
260 this%coef%area(idx(1), idx(2), facet, idx(4))
262 this%flux(2)%x(i) * &
263 this%coef%area(idx(1), idx(2), facet, idx(4))
265 this%flux(3)%x(i) * &
266 this%coef%area(idx(1), idx(2), facet, idx(4))
275 class(
neumann_t),
intent(inout),
target :: this
276 type(c_ptr),
intent(inout) :: x_d
278 logical,
intent(in),
optional :: strong
279 type(c_ptr),
intent(inout) :: strm
282 if (
present(strong))
then
288 if (.not. this%uniform_0 .and. this%msk(0) .gt. 0 .and. &
291 this%flux(1)%x_d, this%coef%area_d, this%coef%Xh%lx, &
292 size(this%msk), strm)
300 class(
neumann_t),
intent(inout),
target :: this
301 type(c_ptr),
intent(inout) :: x_d
302 type(c_ptr),
intent(inout) :: y_d
303 type(c_ptr),
intent(inout) :: z_d
305 logical,
intent(in),
optional :: strong
306 type(c_ptr),
intent(inout) :: strm
309 if (
present(strong))
then
315 if (.not. this%uniform_0 .and. this%msk(0) .gt. 0 .and. &
319 this%flux(1)%x_d, this%flux(2)%x_d, this%flux(3)%x_d, &
320 this%coef%area_d, this%coef%Xh%lx, &
321 size(this%msk), strm)
328 class(
neumann_t),
target,
intent(inout) :: this
330 call this%free_base()
336 class(
neumann_t),
target,
intent(inout) :: this
337 logical,
optional,
intent(in) :: only_facets
340 if (
present(only_facets))
then
341 if (.not. only_facets)
then
342 call neko_error(
"For neumann_t, only_facets has to be true.")
346 call this%finalize_base(.true.)
349 do i = 1,
size(this%init_flux_)
350 call this%flux(i)%init(this%msk(0))
351 this%flux(i) = this%init_flux_(i)
354 this%uniform_0 = .true.
357 this%uniform_0 =
abscmp(this%init_flux_(i), 0.0_rp) .and. this%uniform_0
366 real(kind=
rp),
intent(in) :: flux
367 integer,
intent(in) :: comp
369 if (
size(this%flux) .lt. comp)
then
370 call neko_error(
"Component index out of bounds in " // &
371 "neumann_set_flux_scalar")
374 this%flux(comp) = flux
377 this%uniform_0 =
abscmp(flux, 0.0_rp) .and. this%uniform_0
387 integer,
intent(in) :: comp
390 if (
size(this%flux) .lt. comp)
then
391 call neko_error(
"Component index out of bounds in " // &
392 "neuman_set_flux_array")
395 this%flux(comp) = flux
398 this%uniform_0 = .false.
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Copy data between host and device (or device and device)
Defines a boundary condition.
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
subroutine, public device_neumann_apply_scalar(msk, facet, x, flux, area, lx, m, strm)
subroutine, public device_neumann_apply_vector(msk, facet, x, y, z, flux_x, flux_y, flux_z, area, lx, m, strm)
Device abstraction, common interface for various accelerators.
integer, parameter, public device_to_host
Utilities for retrieving parameters from the case files.
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
subroutine, public copy(a, b, n)
Copy a vector .
integer, parameter neko_bcknd_device
Defines a Neumann boundary condition.
subroutine neumann_init_from_components_array(this, coef, flux)
Constructor from components, using a flux array for vector components.
subroutine neumann_set_flux_scalar(this, flux, comp)
Set the flux using a scalar.
subroutine neumann_set_flux_array(this, flux, comp)
Set a flux component using a vector_t of values.
subroutine neumann_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
Boundary condition apply for a generic Neumann condition to vectors x, y and z (device version)
subroutine neumann_finalize(this, only_facets)
Finalize by setting the flux.
subroutine neumann_init_from_components_single(this, coef, flux)
Constructor from components, using an signle flux.
subroutine neumann_init(this, coef, json)
Constructor.
subroutine neumann_apply_vector(this, x, y, z, n, time, strong)
Boundary condition apply for a generic Neumann condition to vectors x, y and z.
subroutine neumann_apply_scalar(this, x, n, time, strong)
Boundary condition apply for a generic Neumann condition to a vector x.
subroutine neumann_apply_scalar_dev(this, x_d, time, strong, strm)
Boundary condition apply for a generic Neumann condition to a vector x (device version)
subroutine neumann_free(this)
Destructor.
integer, parameter, public rp
Global precision used in computations.
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. Sets the flux of the field to the chosen values.
A struct that contains all info about the time, expand as needed.