Neko  0.8.1
A portable framework for high-order spectral element flow simulations
usr_inflow.f90
Go to the documentation of this file.
1 ! Copyright (c) 2020-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 !
34 module usr_inflow
35  use num_types
36  use coefs
37  use inflow
38  use device
40  use utils
41  implicit none
42  private
43 
45  type, public, extends(inflow_t) :: usr_inflow_t
46  type(coef_t), pointer :: c => null()
47  procedure(usr_inflow_eval), nopass, pointer :: eval => null()
48  type(c_ptr), private :: usr_x_d = c_null_ptr
49  type(c_ptr), private :: usr_y_d = c_null_ptr
50  type(c_ptr), private :: usr_z_d = c_null_ptr
51  contains
52  procedure, pass(this) :: apply_scalar => usr_inflow_apply_scalar
53  procedure, pass(this) :: apply_vector => usr_inflow_apply_vector
54  procedure, pass(this) :: validate => usr_inflow_validate
55  procedure, pass(this) :: set_coef => usr_inflow_set_coef
56  procedure, pass(this) :: set_eval => usr_inflow_set_eval
57  procedure, pass(this) :: apply_vector_dev => usr_inflow_apply_vector_dev
58  procedure, pass(this) :: apply_scalar_dev => usr_inflow_apply_scalar_dev
59  end type usr_inflow_t
60 
61  abstract interface
62 
79  subroutine usr_inflow_eval(u, v, w, x, y, z, nx, ny, nz, &
80  ix, iy, iz, ie, t, tstep)
81  import rp
82  real(kind=rp), intent(inout) :: u
83  real(kind=rp), intent(inout) :: v
84  real(kind=rp), intent(inout) :: w
85  real(kind=rp), intent(in) :: x
86  real(kind=rp), intent(in) :: y
87  real(kind=rp), intent(in) :: z
88  real(kind=rp), intent(in) :: nx
89  real(kind=rp), intent(in) :: ny
90  real(kind=rp), intent(in) :: nz
91  integer, intent(in) :: ix
92  integer, intent(in) :: iy
93  integer, intent(in) :: iz
94  integer, intent(in) :: ie
95  real(kind=rp), intent(in) :: t
96  integer, intent(in) :: tstep
97  end subroutine usr_inflow_eval
98  end interface
99 
100  public :: usr_inflow_eval
101 
102 contains
103 
104  subroutine usr_inflow_free(this)
105  type(usr_inflow_t), intent(inout) :: this
106 
107  if (c_associated(this%usr_x_d)) then
108  call device_free(this%usr_x_d)
109  end if
110 
111  if (c_associated(this%usr_y_d)) then
112  call device_free(this%usr_y_d)
113  end if
114 
115  if (c_associated(this%usr_z_d)) then
116  call device_free(this%usr_z_d)
117  end if
118 
119  end subroutine usr_inflow_free
120 
122  subroutine usr_inflow_apply_scalar(this, x, n, t, tstep)
123  class(usr_inflow_t), intent(inout) :: this
124  integer, intent(in) :: n
125  real(kind=rp), intent(inout), dimension(n) :: x
126  real(kind=rp), intent(in), optional :: t
127  integer, intent(in), optional :: tstep
128  end subroutine usr_inflow_apply_scalar
129 
131  subroutine usr_inflow_apply_scalar_dev(this, x_d, t, tstep)
132  class(usr_inflow_t), intent(inout), target :: this
133  type(c_ptr) :: x_d
134  real(kind=rp), intent(in), optional :: t
135  integer, intent(in), optional :: tstep
136  end subroutine usr_inflow_apply_scalar_dev
137 
139  subroutine usr_inflow_apply_vector(this, x, y, z, n, t, tstep)
140  class(usr_inflow_t), intent(inout) :: this
141  integer, intent(in) :: n
142  real(kind=rp), intent(inout), dimension(n) :: x
143  real(kind=rp), intent(inout), dimension(n) :: y
144  real(kind=rp), intent(inout), dimension(n) :: z
145  real(kind=rp), intent(in), optional :: t
146  integer, intent(in), optional :: tstep
147  integer :: i, m, k, idx(4), facet, tstep_
148  real(kind=rp) :: t_
149 
150  if (present(t)) then
151  t_ = t
152  else
153  t_ = 0.0_rp
154  end if
155 
156  if (present(tstep)) then
157  tstep_ = tstep
158  else
159  tstep_ = 1
160  end if
161 
162  associate(xc => this%c%dof%x, yc => this%c%dof%y, zc => this%c%dof%z, &
163  nx => this%c%nx, ny => this%c%ny, nz => this%c%nz, &
164  lx => this%c%Xh%lx)
165  m = this%msk(0)
166  do i = 1, m
167  k = this%msk(i)
168  facet = this%facet(i)
169  idx = nonlinear_index(k, lx, lx, lx)
170  select case(facet)
171  case(1,2)
172  call this%eval(x(k), y(k), z(k), &
173  xc(idx(1), idx(2), idx(3), idx(4)), &
174  yc(idx(1), idx(2), idx(3), idx(4)), &
175  zc(idx(1), idx(2), idx(3), idx(4)), &
176  nx(idx(2), idx(3), facet, idx(4)), &
177  ny(idx(2), idx(3), facet, idx(4)), &
178  nz(idx(2), idx(3), facet, idx(4)), &
179  idx(1), idx(2), idx(3), idx(4), &
180  t_, tstep_)
181  case(3,4)
182  call this%eval(x(k), y(k), z(k), &
183  xc(idx(1), idx(2), idx(3), idx(4)), &
184  yc(idx(1), idx(2), idx(3), idx(4)), &
185  zc(idx(1), idx(2), idx(3), idx(4)), &
186  nx(idx(1), idx(3), facet, idx(4)), &
187  ny(idx(1), idx(3), facet, idx(4)), &
188  nz(idx(1), idx(3), facet, idx(4)), &
189  idx(1), idx(2), idx(3), idx(4), &
190  t_, tstep_)
191  case(5,6)
192  call this%eval(x(k), y(k), z(k), &
193  xc(idx(1), idx(2), idx(3), idx(4)), &
194  yc(idx(1), idx(2), idx(3), idx(4)), &
195  zc(idx(1), idx(2), idx(3), idx(4)), &
196  nx(idx(1), idx(2), facet, idx(4)), &
197  ny(idx(1), idx(2), facet, idx(4)), &
198  nz(idx(1), idx(2), facet, idx(4)), &
199  idx(1), idx(2), idx(3), idx(4), &
200  t_, tstep_)
201  end select
202  end do
203  end associate
204 
205  end subroutine usr_inflow_apply_vector
206 
207  subroutine usr_inflow_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
208  class(usr_inflow_t), intent(inout), target :: this
209  type(c_ptr) :: x_d
210  type(c_ptr) :: y_d
211  type(c_ptr) :: z_d
212  real(kind=rp), intent(in), optional :: t
213  integer, intent(in), optional :: tstep
214  integer :: i, m, k, idx(4), facet, tstep_
215  integer(c_size_t) :: s
216  real(kind=rp) :: t_
217  real(kind=rp), allocatable :: x(:)
218  real(kind=rp), allocatable :: y(:)
219  real(kind=rp), allocatable :: z(:)
220 
221  if (present(t)) then
222  t_ = t
223  else
224  t_ = 0.0_rp
225  end if
226 
227  if (present(tstep)) then
228  tstep_ = tstep
229  else
230  tstep_ = 1
231  end if
232 
233  associate(xc => this%c%dof%x, yc => this%c%dof%y, zc => this%c%dof%z, &
234  nx => this%c%nx, ny => this%c%ny, nz => this%c%nz, &
235  lx => this%c%Xh%lx, usr_x_d => this%usr_x_d, usr_y_d => this%usr_y_d, &
236  usr_z_d => this%usr_z_d)
237 
238  m = this%msk(0)
239 
240 
241  ! Pretabulate values during first call to apply
242  if (.not. c_associated(usr_x_d)) then
243  allocate(x(m), y(m), z(m)) ! Temp arrays
244 
245  s = m*rp
246 
247  call device_alloc(usr_x_d, s)
248  call device_alloc(usr_y_d, s)
249  call device_alloc(usr_z_d, s)
250 
251  associate(xc => this%c%dof%x, yc => this%c%dof%y, zc => this%c%dof%z, &
252  nx => this%c%nx, ny => this%c%ny, nz => this%c%nz, &
253  lx => this%c%Xh%lx)
254  do i = 1, m
255  k = this%msk(i)
256  facet = this%facet(i)
257  idx = nonlinear_index(k, lx, lx, lx)
258  select case(facet)
259  case(1,2)
260  call this%eval(x(i), y(i), z(i), &
261  xc(idx(1), idx(2), idx(3), idx(4)), &
262  yc(idx(1), idx(2), idx(3), idx(4)), &
263  zc(idx(1), idx(2), idx(3), idx(4)), &
264  nx(idx(2), idx(3), facet, idx(4)), &
265  ny(idx(2), idx(3), facet, idx(4)), &
266  nz(idx(2), idx(3), facet, idx(4)), &
267  idx(1), idx(2), idx(3), idx(4), &
268  t_, tstep_)
269  case(3,4)
270  call this%eval(x(i), y(i), z(i), &
271  xc(idx(1), idx(2), idx(3), idx(4)), &
272  yc(idx(1), idx(2), idx(3), idx(4)), &
273  zc(idx(1), idx(2), idx(3), idx(4)), &
274  nx(idx(1), idx(3), facet, idx(4)), &
275  ny(idx(1), idx(3), facet, idx(4)), &
276  nz(idx(1), idx(3), facet, idx(4)), &
277  idx(1), idx(2), idx(3), idx(4), &
278  t_, tstep_)
279  case(5,6)
280  call this%eval(x(i), y(i), z(i), &
281  xc(idx(1), idx(2), idx(3), idx(4)), &
282  yc(idx(1), idx(2), idx(3), idx(4)), &
283  zc(idx(1), idx(2), idx(3), idx(4)), &
284  nx(idx(1), idx(2), facet, idx(4)), &
285  ny(idx(1), idx(2), facet, idx(4)), &
286  nz(idx(1), idx(2), facet, idx(4)), &
287  idx(1), idx(2), idx(3), idx(4), &
288  t_, tstep_)
289  end select
290  end do
291  end associate
292 
293  call device_memcpy(x, usr_x_d, m, host_to_device, sync=.false.)
294  call device_memcpy(y, usr_y_d, m, host_to_device, sync=.false.)
295  call device_memcpy(z, usr_z_d, m, host_to_device, sync=.true.)
296 
297  deallocate(x, y, z)
298  end if
299 
300  call device_inhom_dirichlet_apply_vector(this%msk_d, x_d, y_d, z_d, &
301  usr_x_d, usr_y_d, usr_z_d, m)
302 
303  end associate
304 
305  end subroutine usr_inflow_apply_vector_dev
306 
308  subroutine usr_inflow_set_coef(this, c)
309  class(usr_inflow_t), intent(inout) :: this
310  type(coef_t), target, intent(inout) :: c
311  this%c => c
312  end subroutine usr_inflow_set_coef
313 
316  subroutine usr_inflow_set_eval(this, usr_eval)
317  class(usr_inflow_t), intent(inout) :: this
318  procedure(usr_inflow_eval) :: usr_eval
319  this%eval => usr_eval
320  end subroutine usr_inflow_set_eval
321 
323  subroutine usr_inflow_validate(this)
324  class(usr_inflow_t), intent(inout) :: this
325  logical :: valid
326 
327  valid = .true. ! Assert it's going to be ok...
328  if (.not. associated(this%c)) then
329  call neko_warning('Missing coefficients')
330  valid = .false.
331  end if
332 
333  if (.not. associated(this%eval)) then
334  call neko_warning('Missing eval function')
335  valid = .false.
336  end if
337 
338  if (.not. valid) then
339  call neko_error('Invalid user defined inflow condition')
340  end if
341 
342  end subroutine usr_inflow_validate
343 
344 end module usr_inflow
__device__ void nonlinear_index(const int idx, const int lx, int *index)
Abstract interface defining a user defined inflow condition (pointwise)
Definition: usr_inflow.f90:79
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:172
Defines inflow dirichlet conditions.
Definition: inflow.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines inflow dirichlet conditions.
Definition: usr_inflow.f90:34
subroutine usr_inflow_apply_scalar(this, x, n, t, tstep)
No-op scalar apply.
Definition: usr_inflow.f90:123
subroutine usr_inflow_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
Definition: usr_inflow.f90:208
subroutine usr_inflow_apply_vector(this, x, y, z, n, t, tstep)
Apply user defined inflow conditions (vector valued)
Definition: usr_inflow.f90:140
subroutine usr_inflow_apply_scalar_dev(this, x_d, t, tstep)
No-op scalar apply (device version)
Definition: usr_inflow.f90:132
subroutine usr_inflow_validate(this)
Validate user inflow condition.
Definition: usr_inflow.f90:324
subroutine usr_inflow_free(this)
Definition: usr_inflow.f90:105
subroutine usr_inflow_set_eval(this, usr_eval)
Assign user provided eval function.
Definition: usr_inflow.f90:317
subroutine usr_inflow_set_coef(this, c)
Assign coefficients (facet normals etc)
Definition: usr_inflow.f90:309
Utilities.
Definition: utils.f90:35
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:54
Dirichlet condition for inlet (vector valued)
Definition: inflow.f90:43
User defined dirichlet condition for inlet (vector valued)
Definition: usr_inflow.f90:45