Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
usr_scalar.f90
Go to the documentation of this file.
1! Copyright (c) 2023, 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 device
41 use json_module, only : json_file
42 implicit none
43 private
44
46 type, public, extends(bc_t) :: usr_scalar_t
47 procedure(usr_scalar_bc_eval), nopass, pointer :: eval => null()
48 type(c_ptr), private :: usr_x_d = c_null_ptr
49 contains
50 procedure, pass(this) :: apply_scalar => usr_scalar_apply_scalar
51 procedure, pass(this) :: apply_vector => usr_scalar_apply_vector
52 procedure, pass(this) :: validate => usr_scalar_validate
53 procedure, pass(this) :: set_eval => usr_scalar_set_eval
54 procedure, pass(this) :: apply_vector_dev => usr_scalar_apply_vector_dev
55 procedure, pass(this) :: apply_scalar_dev => usr_scalar_apply_scalar_dev
57 procedure, pass(this) :: init => usr_scalar_init
59 procedure, pass(this) :: free => usr_scalar_free
61 procedure, pass(this) :: finalize => usr_scalar_finalize
62 end type usr_scalar_t
63
64 abstract interface
65
82 subroutine usr_scalar_bc_eval(s, x, y, z, nx, ny, nz, &
83 ix, iy, iz, ie, t, tstep)
84 import rp
85 real(kind=rp), intent(inout) :: s
86 real(kind=rp), intent(in) :: x
87 real(kind=rp), intent(in) :: y
88 real(kind=rp), intent(in) :: z
89 real(kind=rp), intent(in) :: nx
90 real(kind=rp), intent(in) :: ny
91 real(kind=rp), intent(in) :: nz
92 integer, intent(in) :: ix
93 integer, intent(in) :: iy
94 integer, intent(in) :: iz
95 integer, intent(in) :: ie
96 real(kind=rp), intent(in) :: t
97 integer, intent(in) :: tstep
98 end subroutine usr_scalar_bc_eval
99 end interface
100
101
102 public :: usr_scalar_bc_eval
103
104contains
105
109 subroutine usr_scalar_init(this, coef, json)
110 class(usr_scalar_t), intent(inout), target :: this
111 type(coef_t), intent(in) :: coef
112 type(json_file), intent(inout) :: json
113
114 call this%init_base(coef)
115 end subroutine usr_scalar_init
116
117 subroutine usr_scalar_free(this)
118 class(usr_scalar_t), target, intent(inout) :: this
119
120 call this%free_base
121
122 if (c_associated(this%usr_x_d)) then
123 call device_free(this%usr_x_d)
124 end if
125
126 end subroutine usr_scalar_free
127
133 subroutine usr_scalar_apply_scalar(this, x, n, t, tstep, strong)
134 class(usr_scalar_t), intent(inout) :: this
135 integer, intent(in) :: n
136 real(kind=rp), intent(inout), dimension(n) :: x
137 real(kind=rp), intent(in), optional :: t
138 integer, intent(in), optional :: tstep
139 logical, intent(in), optional :: strong
140 integer :: i, m, k, idx(4), facet, tstep_
141 real(kind=rp) :: t_
142 logical :: strong_ = .true.
143
144 if (present(strong)) strong_ = strong
145
146 if (present(t)) then
147 t_ = t
148 else
149 t_ = 0.0_rp
150 end if
151
152 if (present(tstep)) then
153 tstep_ = tstep
154 else
155 tstep_ = 1
156 end if
157
158 associate(xc => this%coef%dof%x, yc => this%coef%dof%y, &
159 zc => this%coef%dof%z, nx => this%coef%nx, ny => this%coef%ny, &
160 nz => this%coef%nz, lx => this%coef%Xh%lx)
161 m = this%msk(0)
162 if (strong_) then
163 do i = 1, m
164 k = this%msk(i)
165 facet = this%facet(i)
166 idx = nonlinear_index(k, lx, lx, lx)
167 select case (facet)
168 case (1, 2)
169 call this%eval(x(k), &
170 xc(idx(1), idx(2), idx(3), idx(4)), &
171 yc(idx(1), idx(2), idx(3), idx(4)), &
172 zc(idx(1), idx(2), idx(3), idx(4)), &
173 nx(idx(2), idx(3), facet, idx(4)), &
174 ny(idx(2), idx(3), facet, idx(4)), &
175 nz(idx(2), idx(3), facet, idx(4)), &
176 idx(1), idx(2), idx(3), idx(4), &
177 t_, tstep_)
178 case (3, 4)
179 call this%eval(x(k), &
180 xc(idx(1), idx(2), idx(3), idx(4)), &
181 yc(idx(1), idx(2), idx(3), idx(4)), &
182 zc(idx(1), idx(2), idx(3), idx(4)), &
183 nx(idx(1), idx(3), facet, idx(4)), &
184 ny(idx(1), idx(3), facet, idx(4)), &
185 nz(idx(1), idx(3), facet, idx(4)), &
186 idx(1), idx(2), idx(3), idx(4), &
187 t_, tstep_)
188 case (5, 6)
189 call this%eval(x(k), &
190 xc(idx(1), idx(2), idx(3), idx(4)), &
191 yc(idx(1), idx(2), idx(3), idx(4)), &
192 zc(idx(1), idx(2), idx(3), idx(4)), &
193 nx(idx(1), idx(2), facet, idx(4)), &
194 ny(idx(1), idx(2), facet, idx(4)), &
195 nz(idx(1), idx(2), facet, idx(4)), &
196 idx(1), idx(2), idx(3), idx(4), &
197 t_, tstep_)
198 end select
199 end do
200 end if
201 end associate
202 end subroutine usr_scalar_apply_scalar
203
209 subroutine usr_scalar_apply_scalar_dev(this, x_d, t, tstep, strong)
210 class(usr_scalar_t), intent(inout), target :: this
211 type(c_ptr) :: x_d
212 real(kind=rp), intent(in), optional :: t
213 integer, intent(in), optional :: tstep
214 logical, intent(in), optional :: strong
215 integer :: i, m, k, idx(4), facet, tstep_
216 real(kind=rp) :: t_
217 integer(c_size_t) :: s
218 real(kind=rp), allocatable :: x(:)
219 logical :: strong_ = .true.
220
221 m = this%msk(0)
222
223 if (present(strong)) strong_ = strong
224
225 if (present(t)) then
226 t_ = t
227 else
228 t_ = 0.0_rp
229 end if
230
231 if (present(tstep)) then
232 tstep_ = tstep
233 else
234 tstep_ = 1
235 end if
236
237 associate(xc => this%coef%dof%x, yc => this%coef%dof%y, &
238 zc => this%coef%dof%z, nx => this%coef%nx, ny => this%coef%ny, &
239 nz => this%coef%nz, lx => this%coef%Xh%lx, &
240 usr_x_d => this%usr_x_d)
241
242
243 ! Pretabulate values during first call to apply
244 if (.not. c_associated(usr_x_d) .and. strong_) then
245 allocate(x(m)) ! Temp arrays
246 s = m*rp
247
248 call device_alloc(this%usr_x_d, s)
249
250 do i = 1, m
251 k = this%msk(i)
252 facet = this%facet(i)
253 idx = nonlinear_index(k, lx, lx, lx)
254 select case(facet)
255 case (1,2)
256 call this%eval(x(i), &
257 xc(idx(1), idx(2), idx(3), idx(4)), &
258 yc(idx(1), idx(2), idx(3), idx(4)), &
259 zc(idx(1), idx(2), idx(3), idx(4)), &
260 nx(idx(2), idx(3), facet, idx(4)), &
261 ny(idx(2), idx(3), facet, idx(4)), &
262 nz(idx(2), idx(3), facet, idx(4)), &
263 idx(1), idx(2), idx(3), idx(4), &
264 t_, tstep_)
265 case (3,4)
266 call this%eval(x(i), &
267 xc(idx(1), idx(2), idx(3), idx(4)), &
268 yc(idx(1), idx(2), idx(3), idx(4)), &
269 zc(idx(1), idx(2), idx(3), idx(4)), &
270 nx(idx(1), idx(3), facet, idx(4)), &
271 ny(idx(1), idx(3), facet, idx(4)), &
272 nz(idx(1), idx(3), facet, idx(4)), &
273 idx(1), idx(2), idx(3), idx(4), &
274 t_, tstep_)
275 case (5,6)
276 call this%eval(x(i), &
277 xc(idx(1), idx(2), idx(3), idx(4)), &
278 yc(idx(1), idx(2), idx(3), idx(4)), &
279 zc(idx(1), idx(2), idx(3), idx(4)), &
280 nx(idx(1), idx(2), facet, idx(4)), &
281 ny(idx(1), idx(2), facet, idx(4)), &
282 nz(idx(1), idx(2), facet, idx(4)), &
283 idx(1), idx(2), idx(3), idx(4), &
284 t_, tstep_)
285 end select
286 end do
287
288 call device_memcpy(x, this%usr_x_d, m, host_to_device, sync = .true.)
289
290 deallocate(x)
291 end if
292
293 if (strong_) then
294 call device_inhom_dirichlet_apply_scalar(this%msk_d, x_d, &
295 this%usr_x_d, m)
296 end if
297 end associate
298
299
300 end subroutine usr_scalar_apply_scalar_dev
301
303 subroutine usr_scalar_apply_vector(this, x, y, z, n, t, tstep, strong)
304 class(usr_scalar_t), intent(inout) :: this
305 integer, intent(in) :: n
306 real(kind=rp), intent(inout), dimension(n) :: x
307 real(kind=rp), intent(inout), dimension(n) :: y
308 real(kind=rp), intent(inout), dimension(n) :: z
309 real(kind=rp), intent(in), optional :: t
310 integer, intent(in), optional :: tstep
311 logical, intent(in), optional :: strong
312
313 end subroutine usr_scalar_apply_vector
314
316 subroutine usr_scalar_apply_vector_dev(this, x_d, y_d, z_d, t, tstep, strong)
317 class(usr_scalar_t), intent(inout), target :: this
318 type(c_ptr) :: x_d
319 type(c_ptr) :: y_d
320 type(c_ptr) :: z_d
321 real(kind=rp), intent(in), optional :: t
322 integer, intent(in), optional :: tstep
323 logical, intent(in), optional :: strong
324
325 end subroutine usr_scalar_apply_vector_dev
326
329 subroutine usr_scalar_set_eval(this, user_scalar_bc)
330 class(usr_scalar_t), intent(inout) :: this
331 procedure(usr_scalar_bc_eval) :: user_scalar_bc
332 this%eval => user_scalar_bc
333 end subroutine usr_scalar_set_eval
334
336 subroutine usr_scalar_validate(this)
337 class(usr_scalar_t), intent(inout) :: this
338 logical :: valid
339
340 valid = .true. ! Assert it's going to be ok...
341 if (.not. associated(this%coef)) then
342 call neko_warning('Missing coefficients')
343 valid = .false.
344 end if
345
346 if (.not. associated(this%eval)) then
347 call neko_warning('Missing eval function')
348 valid = .false.
349 end if
350
351 if (.not. valid) then
352 call neko_error('Invalid user defined scalar condition')
353 end if
354
355 end subroutine usr_scalar_validate
356
358 subroutine usr_scalar_finalize(this)
359 class(usr_scalar_t), target, intent(inout) :: this
360
361 call this%finalize_base()
362 end subroutine usr_scalar_finalize
363
364end module usr_scalar
__device__ void nonlinear_index(const int idx, const int lx, int *index)
Abstract interface defining a user defined scalar boundary condition (pointwise) Just imitating inflo...
Defines a boundary condition.
Definition bc.f90:34
Coefficients.
Definition coef.f90:34
Device abstraction, common interface for various accelerators.
Definition device.F90:34
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:200
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines dirichlet conditions for scalars.
subroutine usr_scalar_apply_vector_dev(this, x_d, y_d, z_d, t, tstep, strong)
No-op vector apply (device version)
subroutine usr_scalar_init(this, coef, json)
Constructor.
subroutine usr_scalar_set_eval(this, user_scalar_bc)
Assign user provided eval function.
subroutine usr_scalar_apply_scalar(this, x, n, t, tstep, strong)
Scalar apply Just imitating inflow for now, but we should look this over Applies boundary conditions ...
subroutine usr_scalar_apply_scalar_dev(this, x_d, t, tstep, strong)
Scalar apply (device version) Just imitating inflow for now, but we should look this over Applies bou...
subroutine usr_scalar_validate(this)
Validate user scalar condition.
subroutine usr_scalar_free(this)
subroutine usr_scalar_apply_vector(this, x, y, z, n, t, tstep, strong)
No-op vector apply.
subroutine usr_scalar_finalize(this)
Finalize.
Utilities.
Definition utils.f90:35
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:266
Base type for a boundary condition.
Definition bc.f90:57
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
User defined dirichlet condition for scalars.