Neko  0.8.99
A portable framework for high-order spectral element flow simulations
adv_no_dealias.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021-2024, 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 advection, only: advection_t
36  use num_types, only: rp
37  use math, only: subcol3, rzero
38  use space, only: space_t
39  use field, only: field_t
40  use coefs, only: coef_t
43  use operators, only: conv1
45  use, intrinsic :: iso_c_binding, only: c_ptr, c_null_ptr, c_associated
46  implicit none
47  private
48 
50  type, public, extends(advection_t) :: adv_no_dealias_t
51  real(kind=rp), allocatable :: temp(:)
52  type(c_ptr) :: temp_d = c_null_ptr
53  contains
55  procedure, pass(this) :: init => init_no_dealias
57  procedure, pass(this) :: free => free_no_dealias
60  procedure, pass(this) :: compute => compute_advection_no_dealias
63  procedure, pass(this) :: compute_scalar => &
65  end type adv_no_dealias_t
66 
67 contains
68 
71  subroutine init_no_dealias(this, coef)
72  class(adv_no_dealias_t), intent(inout) :: this
73  type(coef_t), intent(in) :: coef
74 
75  allocate(this%temp(coef%dof%size()))
76 
77  if (neko_bcknd_device .eq. 1) then
78  call device_map(this%temp, this%temp_d, coef%dof%size())
79  end if
80 
81  end subroutine init_no_dealias
82 
84  subroutine free_no_dealias(this)
85  class(adv_no_dealias_t), intent(inout) :: this
86 
87  if (allocated(this%temp)) then
88  deallocate(this%temp)
89  end if
90  if (c_associated(this%temp_d)) then
91  call device_free(this%temp_d)
92  end if
93  end subroutine free_no_dealias
94 
107  subroutine compute_advection_no_dealias(this, vx, vy, vz, fx, fy, fz, Xh, &
108  coef, n, dt)
109  class(adv_no_dealias_t), intent(inout) :: this
110  type(space_t), intent(inout) :: Xh
111  type(coef_t), intent(inout) :: coef
112  type(field_t), intent(inout) :: vx, vy, vz
113  type(field_t), intent(inout) :: fx, fy, fz
114  integer, intent(in) :: n
115  real(kind=rp), intent(in), optional :: dt
116 
117  if (neko_bcknd_device .eq. 1) then
118 
119  call conv1(this%temp, vx%x, vx%x, vy%x, vz%x, xh, coef)
120  call device_subcol3 (fx%x_d, coef%B_d, this%temp_d, n)
121  call conv1(this%temp, vy%x, vx%x, vy%x, vz%x, xh, coef)
122  call device_subcol3 (fy%x_d, coef%B_d, this%temp_d, n)
123  if (coef%Xh%lz .eq. 1) then
124  call device_rzero (this%temp_d, n)
125  else
126  call conv1(this%temp, vz%x, vx%x, vy%x, vz%x, xh, coef)
127  call device_subcol3(fz%x_d, coef%B_d, this%temp_d, n)
128  end if
129  else
130  call conv1(this%temp, vx%x, vx%x, vy%x, vz%x, xh, coef)
131  call subcol3 (fx%x, coef%B, this%temp, n)
132  call conv1(this%temp, vy%x, vx%x, vy%x, vz%x, xh, coef)
133  call subcol3 (fy%x, coef%B, this%temp, n)
134  if (coef%Xh%lz .eq. 1) then
135  call rzero (this%temp, n)
136  else
137  call conv1(this%temp, vz%x, vx%x, vy%x, vz%x, xh, coef)
138  call subcol3(fz%x, coef%B, this%temp, n)
139  end if
140  end if
141 
142  end subroutine compute_advection_no_dealias
143 
156  subroutine compute_scalar_advection_no_dealias(this, vx, vy, vz, s, fs, Xh, &
157  coef, n, dt)
158  class(adv_no_dealias_t), intent(inout) :: this
159  type(field_t), intent(inout) :: vx, vy, vz
160  type(field_t), intent(inout) :: s
161  type(field_t), intent(inout) :: fs
162  type(space_t), intent(inout) :: Xh
163  type(coef_t), intent(inout) :: coef
164  integer, intent(in) :: n
165  real(kind=rp), intent(in), optional :: dt
166 
167  if (neko_bcknd_device .eq. 1) then
168 
169  call conv1(this%temp, s%x, vx%x, vy%x, vz%x, xh, coef)
170  call device_subcol3 (fs%x_d, coef%B_d, this%temp_d, n)
171  if (coef%Xh%lz .eq. 1) then
172  call device_rzero (this%temp_d, n)
173  end if
174  else
175  ! temp will hold vx*ds/dx + vy*ds/dy + vz*ds/sz
176  call conv1(this%temp, s%x, vx%x, vy%x, vz%x, xh, coef)
177 
178  ! fs = fs - B*temp
179  call subcol3 (fs%x, coef%B, this%temp, n)
180  if (coef%Xh%lz .eq. 1) then
181  call rzero (this%temp, n)
182  end if
183  end if
184 
186 
187 end module adv_no_dealias
Return the device pointer for an associated Fortran array.
Definition: device.F90:81
Map a Fortran array to a device (allocate and associate)
Definition: device.F90:57
Subroutines to add advection terms to the RHS of a transport equation.
subroutine compute_advection_no_dealias(this, vx, vy, vz, fx, fy, fz, Xh, coef, n, dt)
Add the advection term for the fluid, i.e. to the RHS.
subroutine init_no_dealias(this, coef)
Constructor.
subroutine compute_scalar_advection_no_dealias(this, vx, vy, vz, s, fs, Xh, coef, n, dt)
Add the advection term for a scalar, i.e. , to the RHS.
subroutine free_no_dealias(this)
Destructor.
Subroutines to add advection terms to the RHS of a transport equation.
Definition: advection.f90:34
Coefficients.
Definition: coef.f90:34
subroutine, public device_rzero(a_d, n)
Zero a real vector.
subroutine, public device_subcol3(a_d, b_d, c_d, n)
Returns .
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 field.
Definition: field.f90:34
Definition: math.f90:60
subroutine, public subcol3(a, b, c, n)
Returns .
Definition: math.f90:716
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:184
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
Operators.
Definition: operators.f90:34
subroutine, public conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
Compute the advection term.
Definition: operators.f90:272
Defines a function space.
Definition: space.f90:34
Type encapsulating advection routines with no dealiasing applied.
Base abstract type for computing the advection operator.
Definition: advection.f90:46
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:55
The function space for the SEM solution fields.
Definition: space.f90:62