Neko 1.99.2
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
41 use json_utils, only : json_get
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(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 do i = 1, m
186 k = this%msk(i)
187 facet = this%facet(i)
188 idx = nonlinear_index(k, this%coef%Xh%lx, this%coef%Xh%lx,&
189 this%coef%Xh%lx)
190 select case (facet)
191 case (1,2)
192 x(k) = x(k) + &
193 this%flux(1)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
194 case (3,4)
195 x(k) = x(k) + &
196 this%flux(1)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
197 case (5,6)
198 x(k) = x(k) + &
199 this%flux(1)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
200 end select
201 end do
202 end if
203 end subroutine neumann_apply_scalar
204
207 subroutine neumann_apply_vector(this, x, y, z, n, time, strong)
208 class(neumann_t), intent(inout) :: this
209 integer, intent(in) :: n
210 real(kind=rp), intent(inout), dimension(n) :: x
211 real(kind=rp), intent(inout), dimension(n) :: y
212 real(kind=rp), intent(inout), dimension(n) :: z
213 type(time_state_t), intent(in), optional :: time
214 logical, intent(in), optional :: strong
215 integer :: i, m, k, facet
216 ! Store non-linear index
217 integer :: idx(4)
218 logical :: strong_
219
220 if (present(strong)) then
221 strong_ = strong
222 else
223 strong_ = .true.
224 end if
225
226 m = this%msk(0)
227 if (.not. strong_) then
228 do i = 1, m
229 k = this%msk(i)
230 facet = this%facet(i)
231 idx = nonlinear_index(k, this%coef%Xh%lx, this%coef%Xh%lx,&
232 this%coef%Xh%lx)
233 select case (facet)
234 case (1,2)
235 x(k) = x(k) + &
236 this%flux(1)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
237 y(k) = y(k) + &
238 this%flux(2)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
239 z(k) = z(k) + &
240 this%flux(3)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
241 case (3,4)
242 x(k) = x(k) + &
243 this%flux(1)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
244 y(k) = y(k) + &
245 this%flux(2)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
246 z(k) = z(k) + &
247 this%flux(3)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
248 case (5,6)
249 x(k) = x(k) + &
250 this%flux(1)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
251 y(k) = y(k) + &
252 this%flux(2)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
253 z(k) = z(k) + &
254 this%flux(3)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
255 end select
256 end do
257 end if
258 end subroutine neumann_apply_vector
259
262 subroutine neumann_apply_scalar_dev(this, x_d, time, strong, strm)
263 class(neumann_t), intent(inout), target :: this
264 type(c_ptr), intent(inout) :: x_d
265 type(time_state_t), intent(in), optional :: time
266 logical, intent(in), optional :: strong
267 type(c_ptr), intent(inout) :: strm
268 logical :: strong_
269
270 if (present(strong)) then
271 strong_ = strong
272 else
273 strong_ = .true.
274 end if
275
276 if (.not. this%uniform_0 .and. this%msk(0) .gt. 0 .and. &
277 .not. strong_) then
278 call device_neumann_apply_scalar(this%msk_d, this%facet_d, x_d, &
279 this%flux(1)%x_d, this%coef%area_d, this%coef%Xh%lx, &
280 size(this%msk), strm)
281 end if
282 end subroutine neumann_apply_scalar_dev
283
286 subroutine neumann_apply_vector_dev(this, x_d, y_d, z_d, &
287 time, strong, strm)
288 class(neumann_t), intent(inout), target :: this
289 type(c_ptr), intent(inout) :: x_d
290 type(c_ptr), intent(inout) :: y_d
291 type(c_ptr), intent(inout) :: z_d
292 type(time_state_t), intent(in), optional :: time
293 logical, intent(in), optional :: strong
294 type(c_ptr), intent(inout) :: strm
295 logical :: strong_
296
297 if (present(strong)) then
298 strong_ = strong
299 else
300 strong_ = .true.
301 end if
302
303 if (.not. this%uniform_0 .and. this%msk(0) .gt. 0 .and. &
304 .not. strong_) then
305 call device_neumann_apply_vector(this%msk_d, this%facet_d, &
306 x_d, y_d, z_d, &
307 this%flux(1)%x_d, this%flux(2)%x_d, this%flux(3)%x_d, &
308 this%coef%area_d, this%coef%Xh%lx, &
309 size(this%msk), strm)
310 end if
311
312 end subroutine neumann_apply_vector_dev
313
315 subroutine neumann_free(this)
316 class(neumann_t), target, intent(inout) :: this
317
318 call this%free_base()
319
320 end subroutine neumann_free
321
323 subroutine neumann_finalize(this, only_facets)
324 class(neumann_t), target, intent(inout) :: this
325 logical, optional, intent(in) :: only_facets
326 integer :: i, j
327
328 if (present(only_facets)) then
329 if (only_facets .eqv. .false.) then
330 call neko_error("For neumann_t, only_facets has to be true.")
331 end if
332 end if
333
334 call this%finalize_base(.true.)
335
336 ! Allocate flux vectors and assign to initial constant values
337 do i = 1,size(this%init_flux_)
338 call this%flux(i)%init(this%msk(0))
339 this%flux(i) = this%init_flux_(i)
340 end do
341
342 this%uniform_0 = .true.
343
344 do i = 1, 3
345 this%uniform_0 = abscmp(this%init_flux_(i), 0.0_rp) .and. this%uniform_0
346 end do
347 end subroutine neumann_finalize
348
352 subroutine neumann_set_flux_scalar(this, flux, comp)
353 class(neumann_t), intent(inout) :: this
354 real(kind=rp), intent(in) :: flux
355 integer, intent(in) :: comp
356
357 if (size(this%flux) .lt. comp) then
358 call neko_error("Component index out of bounds in " // &
359 "neumann_set_flux_scalar")
360 end if
361
362 this%flux(comp) = flux
363 ! If we were uniform zero before, and this comp is set to zero, we are still
364 ! uniform zero
365 this%uniform_0 = abscmp(flux, 0.0_rp) .and. this%uniform_0
366
367 end subroutine neumann_set_flux_scalar
368
372 subroutine neumann_set_flux_array(this, flux, comp)
373 class(neumann_t), intent(inout) :: this
374 type(vector_t), intent(in) :: flux
375 integer, intent(in) :: comp
376 integer :: i
377
378 if (size(this%flux) .lt. comp) then
379 call neko_error("Component index out of bounds in " // &
380 "neuman_set_flux_array")
381 end if
382
383 this%flux(comp) = flux
384
385 ! Once a flux is set explicitly, we no longer assume it is uniform zero.
386 this%uniform_0 = .false.
387
388 end subroutine neumann_set_flux_array
389end 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
Retrieves a parameter by name or throws an error.
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:487
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:249
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:353
subroutine neumann_set_flux_array(this, flux, comp)
Set a flux component using a vector_t of values.
Definition neumann.f90:373
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:288
subroutine neumann_finalize(this, only_facets)
Finalize by setting the flux.
Definition neumann.f90:324
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:208
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:263
subroutine neumann_free(this)
Destructor.
Definition neumann.f90:316
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:61
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
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.