Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
adv_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: vdot3, sub2
38 use space, only: space_t, gl
39 use field, only: field_t
40 use coefs, only: coef_t
44 use operators, only: opgrad
47 use, intrinsic :: iso_c_binding, only: c_ptr, c_null_ptr, c_associated
48 implicit none
49 private
50
52 type, public, extends(advection_t) :: adv_dealias_t
54 type(coef_t) :: coef_gl
56 type(coef_t), pointer :: coef_gll
58 type(interpolator_t) :: gll_to_gl
60 type(space_t) :: xh_gl
62 type(space_t), pointer :: xh_gll
63 real(kind=rp), allocatable :: temp(:), tbf(:)
65 real(kind=rp), allocatable :: tx(:), ty(:), tz(:)
66 real(kind=rp), allocatable :: vr(:), vs(:), vt(:)
68 type(c_ptr) :: temp_d = c_null_ptr
70 type(c_ptr) :: tbf_d = c_null_ptr
72 type(c_ptr) :: tx_d = c_null_ptr
74 type(c_ptr) :: ty_d = c_null_ptr
76 type(c_ptr) :: tz_d = c_null_ptr
78 type(c_ptr) :: vr_d = c_null_ptr
80 type(c_ptr) :: vs_d = c_null_ptr
82 type(c_ptr) :: vt_d = c_null_ptr
83
84 contains
87 procedure, pass(this) :: compute => compute_advection_dealias
90 procedure, pass(this) :: compute_scalar => compute_scalar_advection_dealias
92 procedure, pass(this) :: init => init_dealias
94 procedure, pass(this) :: free => free_dealias
95 end type adv_dealias_t
96
97contains
98
102 subroutine init_dealias(this, lxd, coef)
103 class(adv_dealias_t), target, intent(inout) :: this
104 integer, intent(in) :: lxd
105 type(coef_t), intent(inout), target :: coef
106 integer :: nel, n_GL, n
107
108 call this%Xh_GL%init(gl, lxd, lxd, lxd)
109 this%Xh_GLL => coef%Xh
110 this%coef_GLL => coef
111 call this%GLL_to_GL%init(this%Xh_GL, this%Xh_GLL)
112
113 call this%coef_GL%init(this%Xh_GL, coef%msh)
114
115 nel = coef%msh%nelv
116 n_gl = nel*this%Xh_GL%lxyz
117 n = nel*coef%Xh%lxyz
118 call this%GLL_to_GL%map(this%coef_GL%drdx, coef%drdx, nel, this%Xh_GL)
119 call this%GLL_to_GL%map(this%coef_GL%dsdx, coef%dsdx, nel, this%Xh_GL)
120 call this%GLL_to_GL%map(this%coef_GL%dtdx, coef%dtdx, nel, this%Xh_GL)
121 call this%GLL_to_GL%map(this%coef_GL%drdy, coef%drdy, nel, this%Xh_GL)
122 call this%GLL_to_GL%map(this%coef_GL%dsdy, coef%dsdy, nel, this%Xh_GL)
123 call this%GLL_to_GL%map(this%coef_GL%dtdy, coef%dtdy, nel, this%Xh_GL)
124 call this%GLL_to_GL%map(this%coef_GL%drdz, coef%drdz, nel, this%Xh_GL)
125 call this%GLL_to_GL%map(this%coef_GL%dsdz, coef%dsdz, nel, this%Xh_GL)
126 call this%GLL_to_GL%map(this%coef_GL%dtdz, coef%dtdz, nel, this%Xh_GL)
127 if ((neko_bcknd_hip .eq. 1) .or. (neko_bcknd_cuda .eq. 1) .or. &
128 (neko_bcknd_opencl .eq. 1) .or. (neko_bcknd_sx .eq. 1) .or. &
129 (neko_bcknd_xsmm .eq. 1)) then
130 allocate(this%temp(n_gl))
131 allocate(this%tbf(n_gl))
132 allocate(this%tx(n_gl))
133 allocate(this%ty(n_gl))
134 allocate(this%tz(n_gl))
135 allocate(this%vr(n_gl))
136 allocate(this%vs(n_gl))
137 allocate(this%vt(n_gl))
138 end if
139
140 if (neko_bcknd_device .eq. 1) then
141 call device_map(this%temp, this%temp_d, n_gl)
142 call device_map(this%tbf, this%tbf_d, n_gl)
143 call device_map(this%tx, this%tx_d, n_gl)
144 call device_map(this%ty, this%ty_d, n_gl)
145 call device_map(this%tz, this%tz_d, n_gl)
146 call device_map(this%vr, this%vr_d, n_gl)
147 call device_map(this%vs, this%vs_d, n_gl)
148 call device_map(this%vt, this%vt_d, n_gl)
149 end if
150
151 end subroutine init_dealias
152
154 subroutine free_dealias(this)
155 class(adv_dealias_t), intent(inout) :: this
156
157 if (allocated(this%temp)) then
158 deallocate(this%temp)
159 end if
160
161 if (allocated(this%tbf)) then
162 deallocate(this%tbf)
163 end if
164 if (allocated(this%tx)) then
165 deallocate(this%tx)
166 end if
167 if (allocated(this%ty)) then
168 deallocate(this%ty)
169 end if
170 if (allocated(this%tz)) then
171 deallocate(this%tz)
172 end if
173 if (allocated(this%vr)) then
174 deallocate(this%vr)
175 end if
176 if (allocated(this%vs)) then
177 deallocate(this%vs)
178 end if
179 if (allocated(this%vt)) then
180 deallocate(this%vt)
181 end if
182
183 if (c_associated(this%temp_d)) then
184 call device_free(this%temp_d)
185 end if
186
187 if (c_associated(this%tbf_d)) then
188 call device_free(this%tbf_d)
189 end if
190
191 if (c_associated(this%tx_d)) then
192 call device_free(this%tx_d)
193 end if
194
195 if (c_associated(this%ty_d)) then
196 call device_free(this%ty_d)
197 end if
198
199 if (c_associated(this%tz_d)) then
200 call device_free(this%tz_d)
201 end if
202
203 if (c_associated(this%vr_d)) then
204 call device_free(this%vr_d)
205 end if
206
207 if (c_associated(this%vs_d)) then
208 call device_free(this%vs_d)
209 end if
210
211 if (c_associated(this%vt_d)) then
212 call device_free(this%vt_d)
213 end if
214
215 call this%coef_GL%free()
216 call this%GLL_to_GL%free()
217 call this%Xh_GL%free()
218
219 nullify(this%Xh_GLL)
220 nullify(this%coef_GLL)
221
222 end subroutine free_dealias
223
224
237 subroutine compute_advection_dealias(this, vx, vy, vz, fx, fy, fz, Xh, &
238 coef, n, dt)
239 class(adv_dealias_t), intent(inout) :: this
240 type(space_t), intent(in) :: Xh
241 type(coef_t), intent(in) :: coef
242 type(field_t), intent(inout) :: vx, vy, vz
243 type(field_t), intent(inout) :: fx, fy, fz
244 integer, intent(in) :: n
245 real(kind=rp), intent(in), optional :: dt
246
247 real(kind=rp), dimension(this%Xh_GL%lxyz) :: tx, ty, tz
248 real(kind=rp), dimension(this%Xh_GL%lxyz) :: tfx, tfy, tfz
249 real(kind=rp), dimension(this%Xh_GL%lxyz) :: vr, vs, vt
250 real(kind=rp), dimension(this%Xh_GLL%lxyz) :: tempx, tempy, tempz
251 integer :: e, i, idx, nel, n_GL
252
253 nel = coef%msh%nelv
254 n_gl = nel * this%Xh_GL%lxyz
255
256 !This is extremely primitive and unoptimized on the device //Karp
257 associate(c_gl => this%coef_GL)
258 if (neko_bcknd_device .eq. 1) then
259 call this%GLL_to_GL%map(this%tx, vx%x, nel, this%Xh_GL)
260 call this%GLL_to_GL%map(this%ty, vy%x, nel, this%Xh_GL)
261 call this%GLL_to_GL%map(this%tz, vz%x, nel, this%Xh_GL)
262
263 call opgrad(this%vr, this%vs, this%vt, this%tx, c_gl)
264 call device_vdot3(this%tbf_d, this%vr_d, this%vs_d, this%vt_d, &
265 this%tx_d, this%ty_d, this%tz_d, n_gl)
266 call this%GLL_to_GL%map(this%temp, this%tbf, nel, this%Xh_GLL)
267 call device_sub2(fx%x_d, this%temp_d, n)
268
269
270 call opgrad(this%vr, this%vs, this%vt, this%ty, c_gl)
271 call device_vdot3(this%tbf_d, this%vr_d, this%vs_d, this%vt_d, &
272 this%tx_d, this%ty_d, this%tz_d, n_gl)
273 call this%GLL_to_GL%map(this%temp, this%tbf, nel, this%Xh_GLL)
274 call device_sub2(fy%x_d, this%temp_d, n)
275
276 call opgrad(this%vr, this%vs, this%vt, this%tz, c_gl)
277 call device_vdot3(this%tbf_d, this%vr_d, this%vs_d, this%vt_d, &
278 this%tx_d, this%ty_d, this%tz_d, n_gl)
279 call this%GLL_to_GL%map(this%temp, this%tbf, nel, this%Xh_GLL)
280 call device_sub2(fz%x_d, this%temp_d, n)
281
282 else if ((neko_bcknd_sx .eq. 1) .or. (neko_bcknd_xsmm .eq. 1)) then
283
284 call this%GLL_to_GL%map(this%tx, vx%x, nel, this%Xh_GL)
285 call this%GLL_to_GL%map(this%ty, vy%x, nel, this%Xh_GL)
286 call this%GLL_to_GL%map(this%tz, vz%x, nel, this%Xh_GL)
287
288 call opgrad(this%vr, this%vs, this%vt, this%tx, c_gl)
289 call vdot3(this%tbf, this%vr, this%vs, this%vt, &
290 this%tx, this%ty, this%tz, n_gl)
291 call this%GLL_to_GL%map(this%temp, this%tbf, nel, this%Xh_GLL)
292 call sub2(fx%x, this%temp, n)
293
294
295 call opgrad(this%vr, this%vs, this%vt, this%ty, c_gl)
296 call vdot3(this%tbf, this%vr, this%vs, this%vt, &
297 this%tx, this%ty, this%tz, n_gl)
298 call this%GLL_to_GL%map(this%temp, this%tbf, nel, this%Xh_GLL)
299 call sub2(fy%x, this%temp, n)
300
301 call opgrad(this%vr, this%vs, this%vt, this%tz, c_gl)
302 call vdot3(this%tbf, this%vr, this%vs, this%vt, &
303 this%tx, this%ty, this%tz, n_gl)
304 call this%GLL_to_GL%map(this%temp, this%tbf, nel, this%Xh_GLL)
305 call sub2(fz%x, this%temp, n)
306
307 else
308
309 do e = 1, coef%msh%nelv
310 call this%GLL_to_GL%map(tx, vx%x(1,1,1,e), 1, this%Xh_GL)
311 call this%GLL_to_GL%map(ty, vy%x(1,1,1,e), 1, this%Xh_GL)
312 call this%GLL_to_GL%map(tz, vz%x(1,1,1,e), 1, this%Xh_GL)
313
314 call opgrad(vr, vs, vt, tx, c_gl, e, e)
315 do i = 1, this%Xh_GL%lxyz
316 tfx(i) = tx(i)*vr(i) + ty(i)*vs(i) + tz(i)*vt(i)
317 end do
318
319 call opgrad(vr, vs, vt, ty, c_gl, e, e)
320 do i = 1, this%Xh_GL%lxyz
321 tfy(i) = tx(i)*vr(i) + ty(i)*vs(i) + tz(i)*vt(i)
322 end do
323
324 call opgrad(vr, vs, vt, tz, c_gl, e, e)
325 do i = 1, this%Xh_GL%lxyz
326 tfz(i) = tx(i)*vr(i) + ty(i)*vs(i) + tz(i)*vt(i)
327 end do
328
329 call this%GLL_to_GL%map(tempx, tfx, 1, this%Xh_GLL)
330 call this%GLL_to_GL%map(tempy, tfy, 1, this%Xh_GLL)
331 call this%GLL_to_GL%map(tempz, tfz, 1, this%Xh_GLL)
332
333 idx = (e-1)*this%Xh_GLL%lxyz+1
334 do concurrent(i = 0:this%Xh_GLL%lxyz-1)
335 fx%x(i+idx,1,1,1) = fx%x(i+idx,1,1,1) - tempx(i+1)
336 fy%x(i+idx,1,1,1) = fy%x(i+idx,1,1,1) - tempy(i+1)
337 fz%x(i+idx,1,1,1) = fz%x(i+idx,1,1,1) - tempz(i+1)
338 end do
339 end do
340 end if
341 end associate
342
343 end subroutine compute_advection_dealias
344
357 subroutine compute_scalar_advection_dealias(this, vx, vy, vz, s, fs, Xh, &
358 coef, n, dt)
359 class(adv_dealias_t), intent(inout) :: this
360 type(field_t), intent(inout) :: vx, vy, vz
361 type(field_t), intent(inout) :: s
362 type(field_t), intent(inout) :: fs
363 type(space_t), intent(in) :: Xh
364 type(coef_t), intent(in) :: coef
365 integer, intent(in) :: n
366 real(kind=rp), intent(in), optional :: dt
367
368 real(kind=rp), dimension(this%Xh_GL%lxyz) :: vx_gl, vy_gl, vz_gl, s_gl
369 real(kind=rp), dimension(this%Xh_GL%lxyz) :: dsdx, dsdy, dsdz
370 real(kind=rp), dimension(this%Xh_GL%lxyz) :: f_gl
371 integer :: e, i, idx, nel, n_GL
372 real(kind=rp), dimension(this%Xh_GLL%lxyz) :: temp
373
374 nel = coef%msh%nelv
375 n_gl = nel * this%Xh_GL%lxyz
376
377 associate(c_gl => this%coef_GL)
378 if (neko_bcknd_device .eq. 1) then
379
380 ! Map advecting velocity onto the higher-order space
381 call this%GLL_to_GL%map(this%tx, vx%x, nel, this%Xh_GL)
382 call this%GLL_to_GL%map(this%ty, vy%x, nel, this%Xh_GL)
383 call this%GLL_to_GL%map(this%tz, vz%x, nel, this%Xh_GL)
384
385 ! Map the scalar onto the high-order space
386 call this%GLL_to_GL%map(this%temp, s%x, nel, this%Xh_GL)
387
388 ! Compute the scalar gradient in the high-order space
389 call opgrad(this%vr, this%vs, this%vt, this%temp, c_gl)
390
391 ! Compute the convective term, i.e dot the velocity with the scalar grad
392 call device_vdot3(this%tbf_d, this%vr_d, this%vs_d, this%vt_d, &
393 this%tx_d, this%ty_d, this%tz_d, n_gl)
394
395 ! Map back to the original space (we reuse this%temp)
396 call this%GLL_to_GL%map(this%temp, this%tbf, nel, this%Xh_GLL)
397
398 ! Update the source term
399 call device_sub2(fs%x_d, this%temp_d, n)
400
401 else if ((neko_bcknd_sx .eq. 1) .or. (neko_bcknd_xsmm .eq. 1)) then
402
403 ! Map advecting velocity onto the higher-order space
404 call this%GLL_to_GL%map(this%tx, vx%x, nel, this%Xh_GL)
405 call this%GLL_to_GL%map(this%ty, vy%x, nel, this%Xh_GL)
406 call this%GLL_to_GL%map(this%tz, vz%x, nel, this%Xh_GL)
407
408 ! Map the scalar onto the high-order space
409 call this%GLL_to_GL%map(this%temp, s%x, nel, this%Xh_GL)
410
411 ! Compute the scalar gradient in the high-order space
412 call opgrad(this%vr, this%vs, this%vt, this%temp, c_gl)
413
414 ! Compute the convective term, i.e dot the velocity with the scalar grad
415 call vdot3(this%tbf, this%vr, this%vs, this%vt, &
416 this%tx, this%ty, this%tz, n_gl)
417
418 ! Map back to the original space (we reuse this%temp)
419 call this%GLL_to_GL%map(this%temp, this%tbf, nel, this%Xh_GLL)
420
421 ! Update the source term
422 call sub2(fs%x, this%temp, n)
423
424 else
425 do e = 1, coef%msh%nelv
426 ! Map advecting velocity onto the higher-order space
427 call this%GLL_to_GL%map(vx_gl, vx%x(1,1,1,e), 1, this%Xh_GL)
428 call this%GLL_to_GL%map(vy_gl, vy%x(1,1,1,e), 1, this%Xh_GL)
429 call this%GLL_to_GL%map(vz_gl, vz%x(1,1,1,e), 1, this%Xh_GL)
430
431 ! Map scalar onto the higher-order space
432 call this%GLL_to_GL%map(s_gl, s%x(1,1,1,e), 1, this%Xh_GL)
433
434 ! Gradient of s in the higher-order space
435 call opgrad(dsdx, dsdy, dsdz, s_gl, c_gl, e, e)
436
437 ! vx * ds/dx + vy * ds/dy + vz * ds/dz for each point in the element
438 do i = 1, this%Xh_GL%lxyz
439 f_gl(i) = vx_gl(i)*dsdx(i) + vy_gl(i)*dsdy(i) + vz_gl(i)*dsdz(i)
440 end do
441
442 ! Map back the contructed operator to the original space
443 call this%GLL_to_GL%map(temp, f_gl, 1, this%Xh_GLL)
444
445 idx = (e-1)*this%Xh_GLL%lxyz + 1
446
447 call sub2(fs%x(idx, 1, 1, 1), temp, this%Xh_GLL%lxyz)
448 end do
449 end if
450 end associate
451
453
454end module adv_dealias
Return the device pointer for an associated Fortran array.
Definition device.F90:101
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Subroutines to add advection terms to the RHS of a transport equation.
subroutine init_dealias(this, lxd, coef)
Constructor.
subroutine compute_scalar_advection_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_dealias(this)
Destructor.
subroutine compute_advection_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_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm)
Compute a dot product (3-d version) assuming vector components etc.
subroutine, public device_sub2(a_d, b_d, n, strm)
Vector substraction .
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:219
Defines a field.
Definition field.f90:34
Routines to interpolate between different spaces.
Definition math.f90:60
subroutine, public vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
Definition math.f90:684
subroutine, public sub2(a, b, n)
Vector substraction .
Definition math.f90:768
Build configurations.
integer, parameter neko_bcknd_sx
integer, parameter neko_bcknd_hip
integer, parameter neko_bcknd_device
integer, parameter neko_bcknd_opencl
integer, parameter neko_bcknd_cuda
integer, parameter neko_bcknd_xsmm
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public opgrad(ux, uy, uz, u, coef, es, ee)
Compute the weak gradient of a scalar field, i.e. the gradient multiplied by the mass matrix.
Defines a function space.
Definition space.f90:34
integer, parameter, public gl
Definition space.f90:49
Type encapsulating advection routines with dealiasing.
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:56
Interpolation between two space::space_t.
The function space for the SEM solution fields.
Definition space.f90:63