Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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, addcol3
38 use field_math, only : field_col3
39 use space, only : space_t
40 use field, only : field_t
41 use coefs, only : coef_t
44 use operators, only : conv1, div
46 implicit none
47 private
48
50 type, public, extends(advection_t) :: adv_no_dealias_t
51 contains
53 procedure, pass(this) :: init => init_no_dealias
55 procedure, pass(this) :: free => free_no_dealias
58 procedure, pass(this) :: compute => compute_advection_no_dealias
61 procedure, pass(this) :: compute_scalar => &
65 procedure, pass(this) :: compute_ale => compute_ale_advection_no_dealias
67 procedure, pass(this) :: recompute_metrics => recompute_metrics_no_dealias
68 end type adv_no_dealias_t
69
70contains
71
74 subroutine init_no_dealias(this, coef)
75 class(adv_no_dealias_t), intent(inout) :: this
76 type(coef_t), intent(in) :: coef
77
78 end subroutine init_no_dealias
79
81 subroutine free_no_dealias(this)
82 class(adv_no_dealias_t), intent(inout) :: this
83
84 end subroutine free_no_dealias
85
98 subroutine compute_advection_no_dealias(this, vx, vy, vz, fx, fy, fz, Xh, &
99 coef, n, dt)
100 class(adv_no_dealias_t), intent(inout) :: this
101 type(space_t), intent(in) :: Xh
102 type(coef_t), intent(in) :: coef
103 type(field_t), intent(inout) :: vx, vy, vz
104 type(field_t), intent(inout) :: fx, fy, fz
105 integer, intent(in) :: n
106 real(kind=rp), intent(in), optional :: dt
107
108 type(field_t), pointer :: temp
109 integer :: id_temp
110
111 call neko_scratch_registry%request_field(temp, id_temp, .false.)
112
113 if (neko_bcknd_device .eq. 1) then
114
115 call conv1(temp%x, vx%x, vx%x, vy%x, vz%x, xh, coef)
116 call device_subcol3 (fx%x_d, coef%B_d, temp%x_d, n)
117 call conv1(temp%x, vy%x, vx%x, vy%x, vz%x, xh, coef)
118 call device_subcol3 (fy%x_d, coef%B_d, temp%x_d, n)
119 if (coef%Xh%lz .eq. 1) then
120 call device_rzero (temp%x_d, n)
121 else
122 call conv1(temp%x, vz%x, vx%x, vy%x, vz%x, xh, coef)
123 call device_subcol3(fz%x_d, coef%B_d, temp%x_d, n)
124 end if
125 else
126 call conv1(temp%x, vx%x, vx%x, vy%x, vz%x, xh, coef)
127 call subcol3 (fx%x, coef%B, temp%x, n)
128 call conv1(temp%x, vy%x, vx%x, vy%x, vz%x, xh, coef)
129 call subcol3 (fy%x, coef%B, temp%x, n)
130 if (coef%Xh%lz .eq. 1) then
131 call rzero (temp%x, n)
132 else
133 call conv1(temp%x, vz%x, vx%x, vy%x, vz%x, xh, coef)
134 call subcol3(fz%x, coef%B, temp%x, n)
135 end if
136 end if
137
138 call neko_scratch_registry%relinquish_field(id_temp)
139 end subroutine compute_advection_no_dealias
140
153 subroutine compute_scalar_advection_no_dealias(this, vx, vy, vz, s, fs, Xh, &
154 coef, n, dt)
155 class(adv_no_dealias_t), intent(inout) :: this
156 type(field_t), intent(inout) :: vx, vy, vz
157 type(field_t), intent(inout) :: s
158 type(field_t), intent(inout) :: fs
159 type(space_t), intent(in) :: Xh
160 type(coef_t), intent(in) :: coef
161 integer, intent(in) :: n
162 real(kind=rp), intent(in), optional :: dt
163
164 type(field_t), pointer :: temp
165 integer :: id_temp
166 call neko_scratch_registry%request_field(temp, id_temp, .false.)
167
168 if (neko_bcknd_device .eq. 1) then
169
170 call conv1(temp%x, s%x, vx%x, vy%x, vz%x, xh, coef)
171 call device_subcol3 (fs%x_d, coef%B_d, temp%x_d, n)
172 if (coef%Xh%lz .eq. 1) then
173 call device_rzero (temp%x_d, n)
174 end if
175 else
176 ! temp will hold vx*ds/dx + vy*ds/dy + vz*ds/sz
177 call conv1(temp%x, s%x, vx%x, vy%x, vz%x, xh, coef)
178
179 ! fs = fs - B*temp
180 call subcol3 (fs%x, coef%B, temp%x, n)
181 if (coef%Xh%lz .eq. 1) then
182 call rzero (temp%x, n)
183 end if
184 end if
185
186 call neko_scratch_registry%relinquish_field(id_temp)
188
189 subroutine recompute_metrics_no_dealias(this, coef, moving_boundary)
190 class(adv_no_dealias_t), intent(inout) :: this
191 type(coef_t), intent(in) :: coef
192 logical, intent(in) :: moving_boundary
193 ! no-op
194 end subroutine recompute_metrics_no_dealias
195
196 subroutine compute_ale_advection_no_dealias(this, vx, vy, vz, &
197 wm_x, wm_y, wm_z, fx, fy, fz, Xh, coef, n, dt)
198 class(adv_no_dealias_t), intent(inout) :: this
199 type(field_t), intent(inout) :: vx, vy, vz
200 type(field_t), intent(inout) :: wm_x, wm_y, wm_z
201 type(field_t), intent(inout) :: fx, fy, fz
202 type(space_t), intent(in) :: Xh
203 type(coef_t), intent(in) :: coef
204 integer, intent(in) :: n
205 real(kind=rp), intent(in), optional :: dt
206 type(field_t), pointer :: temp, work_x, work_y, work_z
207 integer :: id(4)
208
209 call neko_scratch_registry%request_field(temp, id(1), .false.)
210 call neko_scratch_registry%request_field(work_x, id(2), .false.)
211 call neko_scratch_registry%request_field(work_y, id(3), .false.)
212 call neko_scratch_registry%request_field(work_z, id(4), .false.)
213
214 ! u.wm_*
215 call field_col3(work_x, vx, wm_x)
216 call field_col3(work_y, vx, wm_y)
217 call field_col3(work_z, vx, wm_z)
218 call div(temp%x_d, work_x%x_d, work_y%x_d, work_z%x_d, coef)
219
220 if (neko_bcknd_device .eq. 1) then
221 call device_addcol3(fx%x_d, coef%B_d, temp%x_d, n)
222 else
223 call addcol3(fx%x, coef%B, temp%x, n)
224 end if
225
226 ! v.wm_*
227 call field_col3(work_x, vy, wm_x)
228 call field_col3(work_y, vy, wm_y)
229 call field_col3(work_z, vy, wm_z)
230 call div(temp%x_d, work_x%x_d, work_y%x_d, work_z%x_d, coef)
231
232 if (neko_bcknd_device .eq. 1) then
233 call device_addcol3(fy%x_d, coef%B_d, temp%x_d, n)
234 else
235 call addcol3(fy%x, coef%B, temp%x, n)
236 end if
237
238 ! w.wm_*
239 call field_col3(work_x, vz, wm_x)
240 call field_col3(work_y, vz, wm_y)
241 call field_col3(work_z, vz, wm_z)
242 call div(temp%x_d, work_x%x_d, work_y%x_d, work_z%x_d, coef)
243
244 if (neko_bcknd_device .eq. 1) then
245 call device_addcol3(fz%x_d, coef%B_d, temp%x_d, n)
246 else
247 call addcol3(fz%x, coef%B, temp%x, n)
248 end if
249
250 call neko_scratch_registry%relinquish_field(id)
251
253end module adv_no_dealias
Compute the divergence of a vector field.
Definition operators.f90:85
Subroutines to add advection terms to the RHS of a transport equation.
subroutine init_no_dealias(this, coef)
Constructor.
subroutine free_no_dealias(this)
Destructor.
subroutine compute_ale_advection_no_dealias(this, vx, vy, vz, wm_x, wm_y, wm_z, fx, fy, fz, xh, coef, n, dt)
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 recompute_metrics_no_dealias(this, coef, moving_boundary)
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.
Subroutines to add advection terms to the RHS of a transport equation.
Definition advection.f90:34
Coefficients.
Definition coef.f90:34
subroutine, public device_addcol3(a_d, b_d, c_d, n, strm)
Returns .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_subcol3(a_d, b_d, c_d, n, strm)
Returns .
subroutine, public field_col3(a, b, c, n)
Vector multiplication with 3 vectors .
Defines a field.
Definition field.f90:34
Definition math.f90:60
subroutine, public subcol3(a, b, c, n)
Returns .
Definition math.f90:1075
subroutine, public addcol3(a, b, c, n)
Returns .
Definition math.f90:1162
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:233
Build configurations.
integer, parameter neko_bcknd_device
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.
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
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:63
The function space for the SEM solution fields.
Definition space.f90:63