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" // &
153 class(
neumann_t),
intent(inout),
target :: this
154 type(
coef_t),
intent(in) :: coef
155 real(kind=
rp),
intent(in) :: flux
157 call this%init_base(coef)
158 allocate(this%init_flux_(1))
159 this%init_flux_(1) = flux
167 integer,
intent(in) :: n
168 real(kind=
rp),
intent(inout),
dimension(n) :: x
170 logical,
intent(in),
optional :: strong
171 integer :: i, m, k, facet
176 if (
present(strong))
then
183 if (.not. strong_)
then
186 facet = this%facet(i)
192 this%flux(1)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
195 this%flux(1)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
198 this%flux(1)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
208 integer,
intent(in) :: n
209 real(kind=
rp),
intent(inout),
dimension(n) :: x
210 real(kind=
rp),
intent(inout),
dimension(n) :: y
211 real(kind=
rp),
intent(inout),
dimension(n) :: z
213 logical,
intent(in),
optional :: strong
214 integer :: i, m, k, facet
219 if (
present(strong))
then
226 if (.not. strong_)
then
229 facet = this%facet(i)
235 this%flux(1)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
237 this%flux(2)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
239 this%flux(3)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
242 this%flux(1)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
244 this%flux(2)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
246 this%flux(3)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
249 this%flux(1)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
251 this%flux(2)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
253 this%flux(3)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
262 class(
neumann_t),
intent(inout),
target :: this
263 type(c_ptr),
intent(inout) :: x_d
265 logical,
intent(in),
optional :: strong
266 type(c_ptr),
intent(inout) :: strm
269 if (
present(strong))
then
275 if (.not. this%uniform_0 .and. this%msk(0) .gt. 0 .and. &
278 this%flux(1)%x_d, this%coef%area_d, this%coef%Xh%lx, &
279 size(this%msk), strm)
287 class(
neumann_t),
intent(inout),
target :: this
288 type(c_ptr),
intent(inout) :: x_d
289 type(c_ptr),
intent(inout) :: y_d
290 type(c_ptr),
intent(inout) :: z_d
292 logical,
intent(in),
optional :: strong
293 type(c_ptr),
intent(inout) :: strm
296 if (
present(strong))
then
302 if (.not. this%uniform_0 .and. this%msk(0) .gt. 0 .and. &
306 this%flux(1)%x_d, this%flux(2)%x_d, this%flux(3)%x_d, &
307 this%coef%area_d, this%coef%Xh%lx, &
308 size(this%msk), strm)
315 class(
neumann_t),
target,
intent(inout) :: this
317 call this%free_base()
323 class(
neumann_t),
target,
intent(inout) :: this
324 logical,
optional,
intent(in) :: only_facets
327 if (
present(only_facets))
then
328 if (only_facets .eqv. .false.)
then
329 call neko_error(
"For neumann_t, only_facets has to be true.")
333 call this%finalize_base(.true.)
336 do i = 1,
size(this%init_flux_)
337 call this%flux(i)%init(this%msk(0))
338 this%flux(i) = this%init_flux_(i)
341 this%uniform_0 = .true.
344 this%uniform_0 =
abscmp(this%init_flux_(i), 0.0_rp) .and. this%uniform_0
353 real(kind=
rp),
intent(in) :: flux
354 integer,
intent(in) :: comp
356 if (
size(this%flux) .lt. comp)
then
357 call neko_error(
"Component index out of bounds in " // &
358 "neumann_set_flux_scalar")
361 this%flux(comp) = flux
364 this%uniform_0 =
abscmp(flux, 0.0_rp) .and. this%uniform_0
374 integer,
intent(in) :: comp
377 if (
size(this%flux) .lt. comp)
then
378 call neko_error(
"Component index out of bounds in " // &
379 "neuman_set_flux_array")
382 this%flux(comp) = flux
385 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)
Retrieves a parameter by name or throws an error.
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.