Neko  0.9.0
A portable framework for high-order spectral element flow simulations
source_scalar.f90
Go to the documentation of this file.
1 ! Copyright (c) 2020-2021, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source_scalar and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source_scalar 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 neko_config
36  use num_types
37  use dofmap
38  use utils
39  use device
40  use device_math
41  use, intrinsic :: iso_c_binding
42  implicit none
43 
46  real(kind=rp), allocatable :: s(:,:,:,:)
47  type(dofmap_t), pointer :: dm
48  type(c_ptr) :: s_d = c_null_ptr
49  procedure(source_scalar_term), pass(f), pointer :: eval => null()
50  procedure(source_scalar_term_pw), nopass, pointer :: eval_pw => null()
51  end type source_scalar_t
52 
54  abstract interface
55  subroutine source_scalar_term(f, t)
56  import source_scalar_t
57  import rp
58  class(source_scalar_t), intent(inout) :: f
59  real(kind=rp), intent(in) :: t
60  end subroutine source_scalar_term
61  end interface
62 
64  abstract interface
65  subroutine source_scalar_term_pw(s, j, k, l, e, t)
66  import rp
67  real(kind=rp), intent(inout) :: s
68  integer, intent(in) :: j
69  integer, intent(in) :: k
70  integer, intent(in) :: l
71  integer, intent(in) :: e
72  real(kind=rp), intent(in) :: t
73  end subroutine source_scalar_term_pw
74  end interface
75 
76 contains
77 
79  subroutine source_scalar_init(f, dm)
80  type(source_scalar_t), intent(inout) :: f
81  type(dofmap_t), intent(inout), target :: dm
82 
83  call source_scalar_free(f)
84 
85  f%dm => dm
86 
87  allocate(f%s(dm%Xh%lx, dm%Xh%ly, dm%Xh%lz, dm%msh%nelv))
88 
89  f%s = 0d0
90 
91  if (neko_bcknd_device .eq. 1) then
92  call device_map(f%s, f%s_d, dm%size())
93  end if
94 
95  end subroutine source_scalar_init
96 
98  subroutine source_scalar_free(f)
99  type(source_scalar_t), intent(inout) :: f
100 
101  if (allocated(f%s)) then
102  deallocate(f%s)
103  end if
104 
105  nullify(f%dm)
106 
107  if (c_associated(f%s_d)) then
108  call device_free(f%s_d)
109  end if
110 
111  end subroutine source_scalar_free
112 
114  subroutine source_scalar_set_type(f, f_eval)
115  type(source_scalar_t), intent(inout) :: f
116  procedure(source_scalar_term) :: f_eval
117  f%eval => f_eval
118  end subroutine source_scalar_set_type
119 
121  subroutine source_scalar_set_pw_type(f, f_eval_pw)
122  type(source_scalar_t), intent(inout) :: f
123  procedure(source_scalar_term_pw) :: f_eval_pw
124  if (neko_bcknd_device .eq. 1) then
125  call neko_error('Pointwise source_scalar terms not supported on accelerators')
126  end if
127  f%eval => source_scalar_eval_pw
128  f%eval_pw => f_eval_pw
129  end subroutine source_scalar_set_pw_type
130 
133  subroutine source_scalar_eval_noforce(f, t)
134  class(source_scalar_t), intent(inout) :: f
135  real(kind=rp), intent(in) :: t
136  if (neko_bcknd_device .eq. 1) then
137  call device_rzero(f%s_d, f%dm%size())
138  else
139  f%s = 0d0
140  end if
141  end subroutine source_scalar_eval_noforce
142 
144  subroutine source_scalar_eval_pw(f, t)
145  class(source_scalar_t), intent(inout) :: f
146  real(kind=rp), intent(in) :: t
147  integer :: j, k, l, e
148  integer :: jj,kk,ll,ee
149 
150  do e = 1, f%dm%msh%nelv
151  ee = e
152  do l = 1, f%dm%Xh%lz
153  ll = l
154  do k = 1, f%dm%Xh%ly
155  kk = k
156  do j = 1, f%dm%Xh%lx
157  jj =j
158  call f%eval_pw(f%s(j,k,l,e), jj, kk, ll, ee, t)
159  end do
160  end do
161  end do
162  end do
163 
164  end subroutine source_scalar_eval_pw
165 
166 end module source_scalar
Map a Fortran array to a device (allocate and associate)
Definition: device.F90:57
Abstract interface defining how to compute a source_scalar term pointwise.
Abstract interface defining how to compute a source_scalar term.
subroutine, public device_rzero(a_d, n)
Zero a real vector.
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:185
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Source terms for scalars.
subroutine source_scalar_eval_noforce(f, t)
Eval routine for zero forcing.
subroutine source_scalar_free(f)
Deallocate a source_scalar term f.
subroutine source_scalar_init(f, dm)
Initialize a source_scalar term f.
subroutine source_scalar_set_type(f, f_eval)
Set the eval method for the source_scalar term f.
subroutine source_scalar_set_pw_type(f, f_eval_pw)
Set the pointwise eval method for the source_scalar term f.
subroutine source_scalar_eval_pw(f, t)
Driver for all pointwise source_scalar term evaluatons.
Utilities.
Definition: utils.f90:35
Defines a source term for the scalar transport equation term .