Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
field_neumann.f90
Go to the documentation of this file.
1! Copyright (c) 2026, 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!
35 use num_types, only : rp
36 use coefs, only : coef_t
37 use bc, only : bc_t
38 use field, only : field_t
39 use field_list, only : field_list_t
40 use vector, only : vector_t
42 use json_module, only : json_file
43 use json_utils, only : json_get
44 use math, only : masked_gather_copy_0
48 use, intrinsic :: iso_c_binding, only : c_ptr
49 use time_state, only : time_state_t
50 implicit none
51 private
52
59 type, public, extends(bc_t) :: field_neumann_t
61 type(field_t) :: field_bc
65 type(vector_t) :: flux
67 procedure(field_neumann_update), nopass, pointer :: update => null()
68 contains
70 procedure, pass(this) :: init => field_neumann_init
72 procedure, pass(this) :: init_from_components => field_neumann_init_from_components
74 procedure, pass(this) :: free => field_neumann_free
76 procedure, pass(this) :: finalize => field_neumann_finalize
78 procedure, pass(this) :: apply_scalar => field_neumann_apply_scalar
80 procedure, pass(this) :: apply_vector => field_neumann_apply_vector
82 procedure, pass(this) :: apply_vector_dev => field_neumann_apply_vector_dev
84 procedure, pass(this) :: apply_scalar_dev => field_neumann_apply_scalar_dev
86 procedure, pass(this), private :: gather_flux => field_neumann_gather_flux
87
88 end type field_neumann_t
89
91 abstract interface
92 subroutine field_neumann_update(fields, bc, time)
94 type(field_list_t), intent(inout) :: fields
95 type(field_neumann_t), intent(in) :: bc
96 type(time_state_t), intent(in) :: time
97 end subroutine field_neumann_update
98 end interface
99
100 public :: field_neumann_update
101
102contains
106 subroutine field_neumann_init(this, coef, json)
107 class(field_neumann_t), intent(inout), target :: this
108 type(coef_t), target, intent(in) :: coef
109 type(json_file), intent(inout) :: json
110 character(len=:), allocatable :: field_name
111
112 call json_get(json, "field_name", field_name)
113 call this%init_from_components(coef, field_name)
114 if (allocated(field_name)) then
115 deallocate(field_name)
116 end if
117
118 end subroutine field_neumann_init
119
123 subroutine field_neumann_init_from_components(this, coef, field_name)
124 class(field_neumann_t), intent(inout), target :: this
125 type(coef_t), intent(in) :: coef
126 character(len=*), intent(in) :: field_name
127
128 call this%init_base(coef)
129 this%strong = .false.
130
131 call this%field_bc%init(this%dof, field_name)
132 call this%field_list%init(1)
133 call this%field_list%assign_to_field(1, this%field_bc)
134
136
138 subroutine field_neumann_free(this)
139 class(field_neumann_t), target, intent(inout) :: this
140
141 call this%free_base
142 call this%field_bc%free()
143 call this%field_list%free()
144 call this%flux%free()
145
146 if (associated(this%update)) then
147 this%update => null()
148 end if
149
150 end subroutine field_neumann_free
151
154 class(field_neumann_t), intent(inout) :: this
155
156 if (this%msk(0) .gt. 0) then
157 if (neko_bcknd_device .eq. 1) then
158 call device_masked_gather_copy_0(this%flux%x_d, this%field_bc%x_d, &
159 this%msk_d, this%field_bc%dof%size(), this%msk(0))
160 else
161 call masked_gather_copy_0(this%flux%x, this%field_bc%x, &
162 this%msk, this%field_bc%dof%size(), this%msk(0))
163 end if
164 end if
165
166 end subroutine field_neumann_gather_flux
167
172 subroutine field_neumann_apply_scalar(this, x, n, time, strong)
173 class(field_neumann_t), intent(inout) :: this
174 integer, intent(in) :: n
175 real(kind=rp), intent(inout), dimension(n) :: x
176 type(time_state_t), intent(in), optional :: time
177 logical, intent(in), optional :: strong
178 integer :: i, m, k, facet
179 integer :: idx(4)
180 logical :: strong_
181
182 if (present(strong)) then
183 strong_ = strong
184 else
185 strong_ = .true.
186 end if
187
188 if (.not. strong_) then
189
190 if (.not. this%updated) then
191 call this%update(this%field_list, this, time)
192 call this%gather_flux()
193 this%updated = .true.
194 end if
195
196 m = this%msk(0)
197 !$omp parallel do private(k, facet, idx)
198 do i = 1, m
199 k = this%msk(i)
200 facet = this%facet(i)
201 idx = nonlinear_index(k, this%coef%Xh%lx, this%coef%Xh%lx, &
202 this%coef%Xh%lx)
203 select case (facet)
204 case (1,2)
205 x(k) = x(k) + &
206 this%flux%x(i) * &
207 this%coef%area(idx(2), idx(3), facet, idx(4))
208 case (3,4)
209 x(k) = x(k) + &
210 this%flux%x(i) * &
211 this%coef%area(idx(1), idx(3), facet, idx(4))
212 case (5,6)
213 x(k) = x(k) + &
214 this%flux%x(i) * &
215 this%coef%area(idx(1), idx(2), facet, idx(4))
216 end select
217 end do
218 !$omp end parallel do
219 end if
220
221 end subroutine field_neumann_apply_scalar
222
227 subroutine field_neumann_apply_scalar_dev(this, x_d, time, strong, strm)
228 class(field_neumann_t), intent(inout), target :: this
229 type(c_ptr), intent(inout) :: x_d
230 type(time_state_t), intent(in), optional :: time
231 logical, intent(in), optional :: strong
232 type(c_ptr), intent(inout) :: strm
233 logical :: strong_
234
235 if (present(strong)) then
236 strong_ = strong
237 else
238 strong_ = .true.
239 end if
240
241 if (.not. strong_) then
242 if (.not. this%updated) then
243 call this%update(this%field_list, this, time)
244 call this%gather_flux()
245 this%updated = .true.
246 end if
247
248 if (this%msk(0) .gt. 0) then
249 call device_neumann_apply_scalar(this%msk_d, this%facet_d, x_d, &
250 this%flux%x_d, this%coef%area_d, this%coef%Xh%lx, &
251 size(this%msk), strm)
252 end if
253 end if
254
255 end subroutine field_neumann_apply_scalar_dev
256
258 subroutine field_neumann_apply_vector(this, x, y, z, n, time, strong)
259 class(field_neumann_t), intent(inout) :: this
260 integer, intent(in) :: n
261 real(kind=rp), intent(inout), dimension(n) :: x
262 real(kind=rp), intent(inout), dimension(n) :: y
263 real(kind=rp), intent(inout), dimension(n) :: z
264 type(time_state_t), intent(in), optional :: time
265 logical, intent(in), optional :: strong
266
267 call neko_error("field_neumann cannot apply vector BCs.")
268
269 end subroutine field_neumann_apply_vector
270
272 subroutine field_neumann_apply_vector_dev(this, x_d, y_d, z_d, time, &
273 strong, strm)
274 class(field_neumann_t), intent(inout), target :: this
275 type(c_ptr), intent(inout) :: x_d
276 type(c_ptr), intent(inout) :: y_d
277 type(c_ptr), intent(inout) :: z_d
278 type(time_state_t), intent(in), optional :: time
279 logical, intent(in), optional :: strong
280 type(c_ptr), intent(inout) :: strm
281
282 call neko_error("field_neumann cannot apply vector BCs.")
283
284 end subroutine field_neumann_apply_vector_dev
285
287 subroutine field_neumann_finalize(this, only_facets)
288 class(field_neumann_t), target, intent(inout) :: this
289 logical, optional, intent(in) :: only_facets
290
291 if (present(only_facets)) then
292 if (.not. only_facets) then
293 call neko_error("For field_neumann_t, only_facets has to be true.")
294 end if
295 end if
296
297 call this%finalize_base(.true.)
298 call this%flux%init(this%msk(0))
299
300 end subroutine field_neumann_finalize
301
302end module field_neumann
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
Abstract interface defining a neumann condition on a list of fields.
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_masked_gather_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
subroutine, public device_neumann_apply_scalar(msk, facet, x, flux, area, lx, m, strm)
Defines user neumann condition for a scalar field.
subroutine field_neumann_gather_flux(this)
Gather field-defined values into compact boundary flux storage.
subroutine field_neumann_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
(No-op) Apply vector (device).
subroutine field_neumann_free(this)
Destructor.
subroutine field_neumann_finalize(this, only_facets)
Finalize.
subroutine field_neumann_apply_scalar_dev(this, x_d, time, strong, strm)
Apply scalar (device).
subroutine field_neumann_apply_vector(this, x, y, z, n, time, strong)
(No-op) Apply vector.
subroutine field_neumann_init(this, coef, json)
Constructor.
subroutine field_neumann_apply_scalar(this, x, n, time, strong)
Apply scalar by adding weak neumann contribution.
subroutine field_neumann_init_from_components(this, coef, field_name)
Constructor from components.
Defines a field.
Definition field.f90:34
Utilities for retrieving parameters from the case files.
Definition math.f90:60
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:313
Build configurations.
integer, parameter neko_bcknd_device
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
field_list_t, To be able to group fields together
User defined neumann condition, for which the user can work with an entire field. The type stores a s...
A struct that contains all info about the time, expand as needed.