Neko  0.8.1
A portable framework for high-order spectral element flow simulations
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 !
34 module usr_scalar
35  use num_types
36  use coefs
37  use dirichlet
38  use device
40  use utils
41  implicit none
42  private
43 
45  type, public, extends(dirichlet_t) :: usr_scalar_t
46  type(coef_t), pointer :: c => null()
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_coef => usr_scalar_set_coef
54  procedure, pass(this) :: set_eval => usr_scalar_set_eval
55  procedure, pass(this) :: apply_vector_dev => usr_scalar_apply_vector_dev
56  procedure, pass(this) :: apply_scalar_dev => usr_scalar_apply_scalar_dev
57  end type usr_scalar_t
58 
59  abstract interface
60 
77  subroutine usr_scalar_bc_eval(s, x, y, z, nx, ny, nz, &
78  ix, iy, iz, ie, t, tstep)
79  import rp
80  real(kind=rp), intent(inout) :: s
81  real(kind=rp), intent(in) :: x
82  real(kind=rp), intent(in) :: y
83  real(kind=rp), intent(in) :: z
84  real(kind=rp), intent(in) :: nx
85  real(kind=rp), intent(in) :: ny
86  real(kind=rp), intent(in) :: nz
87  integer, intent(in) :: ix
88  integer, intent(in) :: iy
89  integer, intent(in) :: iz
90  integer, intent(in) :: ie
91  real(kind=rp), intent(in) :: t
92  integer, intent(in) :: tstep
93  end subroutine usr_scalar_bc_eval
94  end interface
95 
96 
97  public :: usr_scalar_bc_eval
98 
99 contains
100 
101  subroutine usr_inflow_free(this)
102  type(usr_scalar_t), intent(inout) :: this
103 
104  if (c_associated(this%usr_x_d)) then
105  call device_free(this%usr_x_d)
106  end if
107 
108  end subroutine usr_inflow_free
109 
115  subroutine usr_scalar_apply_scalar(this, x, n, t, tstep)
116  class(usr_scalar_t), intent(inout) :: this
117  integer, intent(in) :: n
118  real(kind=rp), intent(inout), dimension(n) :: x
119  real(kind=rp), intent(in), optional :: t
120  integer, intent(in), optional :: tstep
121  integer :: i, m, k, idx(4), facet, tstep_
122  real(kind=rp) :: t_
123 
124  if (present(t)) then
125  t_ = t
126  else
127  t_ = 0.0_rp
128  end if
129 
130  if (present(tstep)) then
131  tstep_ = tstep
132  else
133  tstep_ = 1
134  end if
135 
136  associate(xc => this%c%dof%x, yc => this%c%dof%y, zc => this%c%dof%z, &
137  nx => this%c%nx, ny => this%c%ny, nz => this%c%nz, &
138  lx => this%c%Xh%lx)
139  m = this%msk(0)
140  do i = 1, m
141  k = this%msk(i)
142  facet = this%facet(i)
143  idx = nonlinear_index(k, lx, lx, lx)
144  select case(facet)
145  case(1,2)
146  call this%eval(x(k), &
147  xc(idx(1), idx(2), idx(3), idx(4)), &
148  yc(idx(1), idx(2), idx(3), idx(4)), &
149  zc(idx(1), idx(2), idx(3), idx(4)), &
150  nx(idx(2), idx(3), facet, idx(4)), &
151  ny(idx(2), idx(3), facet, idx(4)), &
152  nz(idx(2), idx(3), facet, idx(4)), &
153  idx(1), idx(2), idx(3), idx(4), &
154  t_, tstep_)
155  case(3,4)
156  call this%eval(x(k), &
157  xc(idx(1), idx(2), idx(3), idx(4)), &
158  yc(idx(1), idx(2), idx(3), idx(4)), &
159  zc(idx(1), idx(2), idx(3), idx(4)), &
160  nx(idx(1), idx(3), facet, idx(4)), &
161  ny(idx(1), idx(3), facet, idx(4)), &
162  nz(idx(1), idx(3), facet, idx(4)), &
163  idx(1), idx(2), idx(3), idx(4), &
164  t_, tstep_)
165  case(5,6)
166  call this%eval(x(k), &
167  xc(idx(1), idx(2), idx(3), idx(4)), &
168  yc(idx(1), idx(2), idx(3), idx(4)), &
169  zc(idx(1), idx(2), idx(3), idx(4)), &
170  nx(idx(1), idx(2), facet, idx(4)), &
171  ny(idx(1), idx(2), facet, idx(4)), &
172  nz(idx(1), idx(2), facet, idx(4)), &
173  idx(1), idx(2), idx(3), idx(4), &
174  t_, tstep_)
175  end select
176  end do
177  end associate
178  end subroutine usr_scalar_apply_scalar
179 
185  subroutine usr_scalar_apply_scalar_dev(this, x_d, t, tstep)
186  class(usr_scalar_t), intent(inout), target :: this
187  type(c_ptr) :: x_d
188  real(kind=rp), intent(in), optional :: t
189  integer, intent(in), optional :: tstep
190  integer :: i, m, k, idx(4), facet, tstep_
191  real(kind=rp) :: t_
192  integer(c_size_t) :: s
193  real(kind=rp), allocatable :: x(:)
194  m = this%msk(0)
195 
196  if (present(t)) then
197  t_ = t
198  else
199  t_ = 0.0_rp
200  end if
201 
202  if (present(tstep)) then
203  tstep_ = tstep
204  else
205  tstep_ = 1
206  end if
207 
208  associate(xc => this%c%dof%x, yc => this%c%dof%y, zc => this%c%dof%z, &
209  nx => this%c%nx, ny => this%c%ny, nz => this%c%nz, &
210  lx => this%c%Xh%lx, usr_x_d => this%usr_x_d)
211 
212 
213  ! Pretabulate values during first call to apply
214  if (.not. c_associated(usr_x_d)) then
215  allocate(x(m)) ! Temp arrays
216  s = m*rp
217 
218  call device_alloc(this%usr_x_d, s)
219 
220  do i = 1, m
221  k = this%msk(i)
222  facet = this%facet(i)
223  idx = nonlinear_index(k, lx, lx, lx)
224  select case(facet)
225  case(1,2)
226  call this%eval(x(i), &
227  xc(idx(1), idx(2), idx(3), idx(4)), &
228  yc(idx(1), idx(2), idx(3), idx(4)), &
229  zc(idx(1), idx(2), idx(3), idx(4)), &
230  nx(idx(2), idx(3), facet, idx(4)), &
231  ny(idx(2), idx(3), facet, idx(4)), &
232  nz(idx(2), idx(3), facet, idx(4)), &
233  idx(1), idx(2), idx(3), idx(4), &
234  t_, tstep_)
235  case(3,4)
236  call this%eval(x(i), &
237  xc(idx(1), idx(2), idx(3), idx(4)), &
238  yc(idx(1), idx(2), idx(3), idx(4)), &
239  zc(idx(1), idx(2), idx(3), idx(4)), &
240  nx(idx(1), idx(3), facet, idx(4)), &
241  ny(idx(1), idx(3), facet, idx(4)), &
242  nz(idx(1), idx(3), facet, idx(4)), &
243  idx(1), idx(2), idx(3), idx(4), &
244  t_, tstep_)
245  case(5,6)
246  call this%eval(x(i), &
247  xc(idx(1), idx(2), idx(3), idx(4)), &
248  yc(idx(1), idx(2), idx(3), idx(4)), &
249  zc(idx(1), idx(2), idx(3), idx(4)), &
250  nx(idx(1), idx(2), facet, idx(4)), &
251  ny(idx(1), idx(2), facet, idx(4)), &
252  nz(idx(1), idx(2), facet, idx(4)), &
253  idx(1), idx(2), idx(3), idx(4), &
254  t_, tstep_)
255  end select
256  end do
257 
258  call device_memcpy(x, this%usr_x_d, m, host_to_device, sync=.true.)
259 
260  deallocate(x)
261  end if
262 
263  call device_inhom_dirichlet_apply_scalar(this%msk_d, x_d, &
264  this%usr_x_d, m)
265  end associate
266 
267 
268  end subroutine usr_scalar_apply_scalar_dev
269 
271  subroutine usr_scalar_apply_vector(this, x, y, z, n, t, tstep)
272  class(usr_scalar_t), intent(inout) :: this
273  integer, intent(in) :: n
274  real(kind=rp), intent(inout), dimension(n) :: x
275  real(kind=rp), intent(inout), dimension(n) :: y
276  real(kind=rp), intent(inout), dimension(n) :: z
277  real(kind=rp), intent(in), optional :: t
278  integer, intent(in), optional :: tstep
279  integer :: i, m, k, idx(4), facet
280  end subroutine usr_scalar_apply_vector
281 
283  subroutine usr_scalar_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
284  class(usr_scalar_t), intent(inout), target :: this
285  type(c_ptr) :: x_d
286  type(c_ptr) :: y_d
287  type(c_ptr) :: z_d
288  real(kind=rp), intent(in), optional :: t
289  integer, intent(in), optional :: tstep
290  integer :: i, m, k, idx(4), facet
291  integer(c_size_t) :: s
292  real(kind=rp), allocatable :: x(:)
293  real(kind=rp), allocatable :: y(:)
294  real(kind=rp), allocatable :: z(:)
295  end subroutine usr_scalar_apply_vector_dev
296 
298  subroutine usr_scalar_set_coef(this, c)
299  class(usr_scalar_t), intent(inout) :: this
300  type(coef_t), target, intent(inout) :: c
301  this%c => c
302  end subroutine usr_scalar_set_coef
303 
306  subroutine usr_scalar_set_eval(this, user_scalar_bc)
307  class(usr_scalar_t), intent(inout) :: this
308  procedure(usr_scalar_bc_eval) :: user_scalar_bc
309  this%eval => user_scalar_bc
310  end subroutine usr_scalar_set_eval
311 
313  subroutine usr_scalar_validate(this)
314  class(usr_scalar_t), intent(inout) :: this
315  logical :: valid
316 
317  valid = .true. ! Assert it's going to be ok...
318  if (.not. associated(this%c)) then
319  call neko_warning('Missing coefficients')
320  valid = .false.
321  end if
322 
323  if (.not. associated(this%eval)) then
324  call neko_warning('Missing eval function')
325  valid = .false.
326  end if
327 
328  if (.not. valid) then
329  call neko_error('Invalid user defined scalar condition')
330  end if
331 
332  end subroutine usr_scalar_validate
333 
334 end 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...
Definition: usr_scalar.f90:77
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 a dirichlet boundary condition.
Definition: dirichlet.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines dirichlet conditions for scalars.
Definition: usr_scalar.f90:34
subroutine usr_scalar_apply_scalar_dev(this, x_d, t, tstep)
Scalar apply (device version) Just imitating inflow for now, but we should look this over Applies bou...
Definition: usr_scalar.f90:186
subroutine usr_scalar_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
No-op vector apply (device version)
Definition: usr_scalar.f90:284
subroutine usr_scalar_apply_vector(this, x, y, z, n, t, tstep)
No-op vector apply.
Definition: usr_scalar.f90:272
subroutine usr_scalar_set_eval(this, user_scalar_bc)
Assign user provided eval function.
Definition: usr_scalar.f90:307
subroutine usr_scalar_apply_scalar(this, x, n, t, tstep)
Scalar apply Just imitating inflow for now, but we should look this over Applies boundary conditions ...
Definition: usr_scalar.f90:116
subroutine usr_scalar_set_coef(this, c)
Assign coefficients (facet normals etc)
Definition: usr_scalar.f90:299
subroutine usr_inflow_free(this)
Definition: usr_scalar.f90:102
subroutine usr_scalar_validate(this)
Validate user scalar condition.
Definition: usr_scalar.f90:314
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
Generic Dirichlet boundary condition on .
Definition: dirichlet.f90:44
User defined dirichlet condition for scalars.
Definition: usr_scalar.f90:45