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
148
152 subroutine neumann_init_from_components_single(this, coef, flux)
153 class(neumann_t), intent(inout), target :: this
154 type(coef_t), intent(in) :: coef
155 real(kind=rp), intent(in) :: flux
156
157 call this%init_base(coef)
158 allocate(this%init_flux_(1))
159 this%init_flux_(1) = flux
160
162
165 subroutine neumann_apply_scalar(this, x, n, time, strong)
166 class(neumann_t), intent(inout) :: this
167 integer, intent(in) :: n
168 real(kind=rp), intent(inout), dimension(n) :: x
169 type(time_state_t), intent(in), optional :: time
170 logical, intent(in), optional :: strong
171 integer :: i, m, k, facet
172 ! Store non-linear index
173 integer :: idx(4)
174 logical :: strong_
175
176 if (present(strong)) then
177 strong_ = strong
178 else
179 strong_ = .true.
180 end if
181
182 m = this%msk(0)
183 if (.not. strong_) then
184 do i = 1, m
185 k = this%msk(i)
186 facet = this%facet(i)
187 idx = nonlinear_index(k, this%coef%Xh%lx, this%coef%Xh%lx,&
188 this%coef%Xh%lx)
189 select case (facet)
190 case (1,2)
191 x(k) = x(k) + &
192 this%flux(1)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
193 case (3,4)
194 x(k) = x(k) + &
195 this%flux(1)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
196 case (5,6)
197 x(k) = x(k) + &
198 this%flux(1)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
199 end select
200 end do
201 end if
202 end subroutine neumann_apply_scalar
203
206 subroutine neumann_apply_vector(this, x, y, z, n, time, strong)
207 class(neumann_t), intent(inout) :: this
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
212 type(time_state_t), intent(in), optional :: time
213 logical, intent(in), optional :: strong
214 integer :: i, m, k, facet
215 ! Store non-linear index
216 integer :: idx(4)
217 logical :: strong_
218
219 if (present(strong)) then
220 strong_ = strong
221 else
222 strong_ = .true.
223 end if
224
225 m = this%msk(0)
226 if (.not. strong_) then
227 do i = 1, m
228 k = this%msk(i)
229 facet = this%facet(i)
230 idx = nonlinear_index(k, this%coef%Xh%lx, this%coef%Xh%lx,&
231 this%coef%Xh%lx)
232 select case (facet)
233 case (1,2)
234 x(k) = x(k) + &
235 this%flux(1)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
236 y(k) = y(k) + &
237 this%flux(2)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
238 z(k) = z(k) + &
239 this%flux(3)%x(i)*this%coef%area(idx(2), idx(3), facet, idx(4))
240 case (3,4)
241 x(k) = x(k) + &
242 this%flux(1)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
243 y(k) = y(k) + &
244 this%flux(2)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
245 z(k) = z(k) + &
246 this%flux(3)%x(i)*this%coef%area(idx(1), idx(3), facet, idx(4))
247 case (5,6)
248 x(k) = x(k) + &
249 this%flux(1)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
250 y(k) = y(k) + &
251 this%flux(2)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
252 z(k) = z(k) + &
253 this%flux(3)%x(i)*this%coef%area(idx(1), idx(2), facet, idx(4))
254 end select
255 end do
256 end if
257 end subroutine neumann_apply_vector
258
261 subroutine neumann_apply_scalar_dev(this, x_d, time, strong, strm)
262 class(neumann_t), intent(inout), target :: this
263 type(c_ptr), intent(inout) :: x_d
264 type(time_state_t), intent(in), optional :: time
265 logical, intent(in), optional :: strong
266 type(c_ptr), intent(inout) :: strm
267 logical :: strong_
268
269 if (present(strong)) then
270 strong_ = strong
271 else
272 strong_ = .true.
273 end if
274
275 if (.not. this%uniform_0 .and. this%msk(0) .gt. 0 .and. &
276 .not. strong_) then
277 call device_neumann_apply_scalar(this%msk_d, this%facet_d, x_d, &
278 this%flux(1)%x_d, this%coef%area_d, this%coef%Xh%lx, &
279 size(this%msk), strm)
280 end if
281 end subroutine neumann_apply_scalar_dev
282
285 subroutine neumann_apply_vector_dev(this, x_d, y_d, z_d, &
286 time, strong, 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
291 type(time_state_t), intent(in), optional :: time
292 logical, intent(in), optional :: strong
293 type(c_ptr), intent(inout) :: strm
294 logical :: strong_
295
296 if (present(strong)) then
297 strong_ = strong
298 else
299 strong_ = .true.
300 end if
301
302 if (.not. this%uniform_0 .and. this%msk(0) .gt. 0 .and. &
303 .not. strong_) then
304 call device_neumann_apply_vector(this%msk_d, this%facet_d, &
305 x_d, y_d, z_d, &
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)
309 end if
310
311 end subroutine neumann_apply_vector_dev
312
314 subroutine neumann_free(this)
315 class(neumann_t), target, intent(inout) :: this
316
317 call this%free_base()
318
319 end subroutine neumann_free
320
322 subroutine neumann_finalize(this, only_facets)
323 class(neumann_t), target, intent(inout) :: this
324 logical, optional, intent(in) :: only_facets
325 integer :: i, j
326
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.")
330 end if
331 end if
332
333 call this%finalize_base(.true.)
334
335 ! Allocate flux vectors and assign to initial constant values
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)
339 end do
340
341 this%uniform_0 = .true.
342
343 do i = 1, 3
344 this%uniform_0 = abscmp(this%init_flux_(i), 0.0_rp) .and. this%uniform_0
345 end do
346 end subroutine neumann_finalize
347
351 subroutine neumann_set_flux_scalar(this, flux, comp)
352 class(neumann_t), intent(inout) :: this
353 real(kind=rp), intent(in) :: flux
354 integer, intent(in) :: comp
355
356 if (size(this%flux) .lt. comp) then
357 call neko_error("Component index out of bounds in " // &
358 "neumann_set_flux_scalar")
359 end if
360
361 this%flux(comp) = flux
362 ! If we were uniform zero before, and this comp is set to zero, we are still
363 ! uniform zero
364 this%uniform_0 = abscmp(flux, 0.0_rp) .and. this%uniform_0
365
366 end subroutine neumann_set_flux_scalar
367
371 subroutine neumann_set_flux_array(this, flux, comp)
372 class(neumann_t), intent(inout) :: this
373 type(vector_t), intent(in) :: flux
374 integer, intent(in) :: comp
375 integer :: i
376
377 if (size(this%flux) .lt. comp) then
378 call neko_error("Component index out of bounds in " // &
379 "neuman_set_flux_array")
380 end if
381
382 this%flux(comp) = flux
383
384 ! Once a flux is set explicitly, we no longer assume it is uniform zero.
385 this%uniform_0 = .false.
386
387 end subroutine neumann_set_flux_array
388end 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:352
subroutine neumann_set_flux_array(this, flux, comp)
Set a flux component using a vector_t of values.
Definition neumann.f90:372
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:287
subroutine neumann_finalize(this, only_facets)
Finalize by setting the flux.
Definition neumann.f90:323
subroutine neumann_init_from_components_single(this, coef, flux)
Constructor from components, using an signle flux.
Definition neumann.f90:153
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:207
subroutine neumann_apply_scalar(this, x, n, time, strong)
Boundary condition apply for a generic Neumann condition to a vector x.
Definition neumann.f90:166
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:262
subroutine neumann_free(this)
Destructor.
Definition neumann.f90:315
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.