Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
neumann.f90
Go to the documentation of this file.
1! Copyright (c) 2024-2025, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
34module neumann
35 use num_types, only : rp
36 use bc, only : bc_t
37 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
39 use coefs, only : coef_t
40 use json_module, only : json_file
42 use math, only : cfill, copy, abscmp
43 use vector, only : vector_t
49 use time_state, only : time_state_t
50 implicit none
51 private
52
60 type, public, extends(bc_t) :: neumann_t
63 type(vector_t), allocatable :: flux(:)
68 real(kind=rp), allocatable, private :: init_flux_(:)
71 logical :: uniform_0 = .false.
72 contains
73 procedure, pass(this) :: apply_scalar => neumann_apply_scalar
74 procedure, pass(this) :: apply_vector => neumann_apply_vector
75 procedure, pass(this) :: apply_scalar_dev => neumann_apply_scalar_dev
76 procedure, pass(this) :: apply_vector_dev => neumann_apply_vector_dev
78 procedure, pass(this) :: init => neumann_init
80 procedure, pass(this) :: neumann_init_from_components_array
82 procedure, pass(this) :: neumann_init_from_components_single
83 generic :: init_from_components => &
87 procedure, pass(this) :: set_flux_scalar => neumann_set_flux_scalar
89 procedure, pass(this) :: set_flux_array => neumann_set_flux_array
91 generic :: set_flux => set_flux_scalar, set_flux_array
93 procedure, pass(this) :: free => neumann_free
95 procedure, pass(this) :: finalize => neumann_finalize
96 end type neumann_t
97
98contains
99
103 subroutine neumann_init(this, coef, json)
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
108 logical :: found
109
110 call this%init_base(coef)
111 this%strong = .false.
112
113 ! Try to read array from json
114 call json%get("flux", this%init_flux_, found)
115
116 ! If we haven't found an array, try to read a single value
117 if (.not. found) then
118 call json_get_or_lookup(json, "flux", flux)
119 allocate(this%init_flux_(1))
120 this%init_flux_(1) = flux
121 end if
122
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" // &
126 " vector.")
127 end if
128
129 allocate(this%flux(size(this%init_flux_)))
130 end subroutine neumann_init
131
135 subroutine neumann_init_from_components_array(this, coef, flux)
136 class(neumann_t), intent(inout), target :: this
137 type(coef_t), intent(in) :: coef
138 real(kind=rp), intent(in) :: flux(3)
139
140 call this%init_base(coef)
141 this%init_flux_ = flux
142
143 if ((size(this%init_flux_) .ne. 3)) then
144 call neko_error("Neumann BC flux must be a scalar or a 3-component" // &
145 " vector.")
146 end if
147 allocate(this%flux(size(this%init_flux_)))
149
153 subroutine neumann_init_from_components_single(this, coef, flux)
154 class(neumann_t), intent(inout), target :: this
155 type(coef_t), intent(in) :: coef
156 real(kind=rp), intent(in) :: flux
157
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_)))
163
166 subroutine neumann_apply_scalar(this, x, n, time, strong)
167 class(neumann_t), intent(inout) :: this
168 integer, intent(in) :: n
169 real(kind=rp), intent(inout), dimension(n) :: x
170 type(time_state_t), intent(in), optional :: time
171 logical, intent(in), optional :: strong
172 integer :: i, m, k, facet
173 ! Store non-linear index
174 integer :: idx(4)
175 logical :: strong_
176
177 if (present(strong)) then
178 strong_ = strong
179 else
180 strong_ = .true.
181 end if
182
183 m = this%msk(0)
184 if (.not. strong_) then
185 !$omp parallel do private(k, facet, idx)
186 do i = 1,m
187 k = this%msk(i)
188 facet = this%facet(i)
189 idx = nonlinear_index(k, this%coef%Xh%lx, this%coef%Xh%lx, &
190 this%coef%Xh%lx)
191 select case (facet)
192 case (1,2)
193 x(k) = x(k) + &
194 this%flux(1)%x(i) * &
195 this%coef%area(idx(2), idx(3), facet, idx(4))
196 case (3,4)
197 x(k) = x(k) + &
198 this%flux(1)%x(i) * &
199 this%coef%area(idx(1), idx(3), facet, idx(4))
200 case (5,6)
201 x(k) = x(k) + &
202 this%flux(1)%x(i) * &
203 this%coef%area(idx(1), idx(2), facet, idx(4))
204 end select
205 end do
206 !$omp end parallel do
207 end if
208 end subroutine neumann_apply_scalar
209
212 subroutine neumann_apply_vector(this, x, y, z, n, time, strong)
213 class(neumann_t), intent(inout) :: this
214 integer, intent(in) :: n
215 real(kind=rp), intent(inout), dimension(n) :: x
216 real(kind=rp), intent(inout), dimension(n) :: y
217 real(kind=rp), intent(inout), dimension(n) :: z
218 type(time_state_t), intent(in), optional :: time
219 logical, intent(in), optional :: strong
220 integer :: i, m, k, facet
221 ! Store non-linear index
222 integer :: idx(4)
223 logical :: strong_
224
225 if (present(strong)) then
226 strong_ = strong
227 else
228 strong_ = .true.
229 end if
230
231 m = this%msk(0)
232 if (.not. strong_) then
233 !$omp parallel do private(k, facet, idx)
234 do i = 1, m
235 k = this%msk(i)
236 facet = this%facet(i)
237 idx = nonlinear_index(k, this%coef%Xh%lx, this%coef%Xh%lx, &
238 this%coef%Xh%lx)
239 select case (facet)
240 case (1,2)
241 x(k) = x(k) + &
242 this%flux(1)%x(i) * &
243 this%coef%area(idx(2), idx(3), facet, idx(4))
244 y(k) = y(k) + &
245 this%flux(2)%x(i) * &
246 this%coef%area(idx(2), idx(3), facet, idx(4))
247 z(k) = z(k) + &
248 this%flux(3)%x(i) * &
249 this%coef%area(idx(2), idx(3), facet, idx(4))
250 case (3,4)
251 x(k) = x(k) + &
252 this%flux(1)%x(i) * &
253 this%coef%area(idx(1), idx(3), facet, idx(4))
254 y(k) = y(k) + &
255 this%flux(2)%x(i) * &
256 this%coef%area(idx(1), idx(3), facet, idx(4))
257 z(k) = z(k) + &
258 this%flux(3)%x(i) * &
259 this%coef%area(idx(1), idx(3), facet, idx(4))
260 case (5,6)
261 x(k) = x(k) + &
262 this%flux(1)%x(i) * &
263 this%coef%area(idx(1), idx(2), facet, idx(4))
264 y(k) = y(k) + &
265 this%flux(2)%x(i) * &
266 this%coef%area(idx(1), idx(2), facet, idx(4))
267 z(k) = z(k) + &
268 this%flux(3)%x(i) * &
269 this%coef%area(idx(1), idx(2), facet, idx(4))
270 end select
271 end do
272 !$omp end parallel do
273 end if
274 end subroutine neumann_apply_vector
275
278 subroutine neumann_apply_scalar_dev(this, x_d, time, strong, strm)
279 class(neumann_t), intent(inout), target :: this
280 type(c_ptr), intent(inout) :: x_d
281 type(time_state_t), intent(in), optional :: time
282 logical, intent(in), optional :: strong
283 type(c_ptr), intent(inout) :: strm
284 logical :: strong_
285
286 if (present(strong)) then
287 strong_ = strong
288 else
289 strong_ = .true.
290 end if
291
292 if (.not. this%uniform_0 .and. this%msk(0) .gt. 0 .and. &
293 .not. strong_) then
294 call device_neumann_apply_scalar(this%msk_d, this%facet_d, x_d, &
295 this%flux(1)%x_d, this%coef%area_d, this%coef%Xh%lx, &
296 size(this%msk), strm)
297 end if
298 end subroutine neumann_apply_scalar_dev
299
302 subroutine neumann_apply_vector_dev(this, x_d, y_d, z_d, &
303 time, strong, strm)
304 class(neumann_t), intent(inout), target :: this
305 type(c_ptr), intent(inout) :: x_d
306 type(c_ptr), intent(inout) :: y_d
307 type(c_ptr), intent(inout) :: z_d
308 type(time_state_t), intent(in), optional :: time
309 logical, intent(in), optional :: strong
310 type(c_ptr), intent(inout) :: strm
311 logical :: strong_
312
313 if (present(strong)) then
314 strong_ = strong
315 else
316 strong_ = .true.
317 end if
318
319 if (.not. this%uniform_0 .and. this%msk(0) .gt. 0 .and. &
320 .not. strong_) then
321 call device_neumann_apply_vector(this%msk_d, this%facet_d, &
322 x_d, y_d, z_d, &
323 this%flux(1)%x_d, this%flux(2)%x_d, this%flux(3)%x_d, &
324 this%coef%area_d, this%coef%Xh%lx, &
325 size(this%msk), strm)
326 end if
327
328 end subroutine neumann_apply_vector_dev
329
331 subroutine neumann_free(this)
332 class(neumann_t), target, intent(inout) :: this
333
334 call this%free_base()
335
336 end subroutine neumann_free
337
339 subroutine neumann_finalize(this, only_facets)
340 class(neumann_t), target, intent(inout) :: this
341 logical, optional, intent(in) :: only_facets
342 integer :: i, j
343
344 if (present(only_facets)) then
345 if (.not. only_facets) then
346 call neko_error("For neumann_t, only_facets has to be true.")
347 end if
348 end if
349
350 call this%finalize_base(.true.)
351
352 ! Allocate flux vectors and assign to initial constant values
353 do i = 1, size(this%init_flux_)
354 call this%flux(i)%init(this%msk(0))
355 this%flux(i) = this%init_flux_(i)
356 end do
357
358 this%uniform_0 = .true.
359
360 do i = 1, 3
361 this%uniform_0 = abscmp(this%init_flux_(i), 0.0_rp) .and. this%uniform_0
362 end do
363 end subroutine neumann_finalize
364
368 subroutine neumann_set_flux_scalar(this, flux, comp)
369 class(neumann_t), intent(inout) :: this
370 real(kind=rp), intent(in) :: flux
371 integer, intent(in) :: comp
372
373 if (size(this%flux) .lt. comp) then
374 call neko_error("Component index out of bounds in " // &
375 "neumann_set_flux_scalar")
376 end if
377
378 this%flux(comp) = flux
379 ! If we were uniform zero before, and this comp is set to zero, we are still
380 ! uniform zero
381 this%uniform_0 = abscmp(flux, 0.0_rp) .and. this%uniform_0
382
383 end subroutine neumann_set_flux_scalar
384
388 subroutine neumann_set_flux_array(this, flux, comp)
389 class(neumann_t), intent(inout) :: this
390 type(vector_t), intent(in) :: flux
391 integer, intent(in) :: comp
392 integer :: i
393
394 if (size(this%flux) .lt. comp) then
395 call neko_error("Component index out of bounds in " // &
396 "neuman_set_flux_array")
397 end if
398
399 this%flux(comp) = flux
400
401 ! Once a flux is set explicitly, we no longer assume it is uniform zero.
402 this%uniform_0 = .false.
403
404 end subroutine neumann_set_flux_array
405end module neumann
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
Copy data between host and device (or device and device)
Definition device.F90:71
Defines a boundary condition.
Definition bc.f90:34
Coefficients.
Definition coef.f90:34
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.
Definition device.F90:34
integer, parameter, public device_to_host
Definition device.F90:47
Utilities for retrieving parameters from the case files.
Definition math.f90:60
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:488
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:250
Build configurations.
integer, parameter neko_bcknd_device
Defines a Neumann boundary condition.
Definition neumann.f90:34
subroutine neumann_init_from_components_array(this, coef, flux)
Constructor from components, using a flux array for vector components.
Definition neumann.f90:136
subroutine neumann_set_flux_scalar(this, flux, comp)
Set the flux using a scalar.
Definition neumann.f90:369
subroutine neumann_set_flux_array(this, flux, comp)
Set a flux component using a vector_t of values.
Definition neumann.f90:389
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)
Definition neumann.f90:304
subroutine neumann_finalize(this, only_facets)
Finalize by setting the flux.
Definition neumann.f90:340
subroutine neumann_init_from_components_single(this, coef, flux)
Constructor from components, using an signle flux.
Definition neumann.f90:154
subroutine neumann_init(this, coef, json)
Constructor.
Definition neumann.f90:104
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.
Definition neumann.f90:213
subroutine neumann_apply_scalar(this, x, n, time, strong)
Boundary condition apply for a generic Neumann condition to a vector x.
Definition neumann.f90:167
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)
Definition neumann.f90:279
subroutine neumann_free(this)
Destructor.
Definition neumann.f90:332
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
Defines a vector.
Definition vector.f90:34
Base type for a boundary condition.
Definition bc.f90:62
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:62
A Neumann boundary condition. Sets the flux of the field to the chosen values.
Definition neumann.f90:60
A struct that contains all info about the time, expand as needed.