Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
operators.f90
Go to the documentation of this file.
1! Copyright (c) 2020-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!
37 use num_types, only : rp, i8
38 use opr_cpu, only : opr_cpu_cfl, opr_cpu_curl, opr_cpu_opgrad, &
39 opr_cpu_conv1, opr_cpu_convect_scalar, opr_cpu_cdtp, &
40 opr_cpu_dudxyz, opr_cpu_lambda2, opr_cpu_set_convect_rst
41 use opr_sx, only : opr_sx_cfl, opr_sx_curl, opr_sx_opgrad, &
42 opr_sx_conv1, opr_sx_convect_scalar, opr_sx_cdtp, &
43 opr_sx_dudxyz, opr_sx_lambda2, opr_sx_set_convect_rst
50 use space, only : space_t
51 use coefs, only : coef_t
52 use field, only : field_t
54 use math, only : glsum, cmult, add2, add3s2, cadd, copy, col2, invcol2, &
56 use device, only : device_get_ptr
60 use comm
61 use, intrinsic :: iso_c_binding, only : c_ptr
62 implicit none
63 private
64
65 public :: dudxyz, opgrad, ortho, cdtp, conv1, curl, cfl, &
67
68contains
69
77 subroutine dudxyz (du, u, dr, ds, dt, coef)
78 type(coef_t), intent(in), target :: coef
79 real(kind=rp), dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, & coef%msh%nelv), intent(inout) :: du
80 real(kind=rp), dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, & coef%msh%nelv), intent(in) :: u, dr, ds, dt
81
82 if (neko_bcknd_sx .eq. 1) then
83 call opr_sx_dudxyz(du, u, dr, ds, dt, coef)
84 else if (neko_bcknd_xsmm .eq. 1) then
85 call opr_xsmm_dudxyz(du, u, dr, ds, dt, coef)
86 else if (neko_bcknd_device .eq. 1) then
87 call opr_device_dudxyz(du, u, dr, ds, dt, coef)
88 else
89 call opr_cpu_dudxyz(du, u, dr, ds, dt, coef)
90 end if
91
92 end subroutine dudxyz
93
100 subroutine div(res, ux, uy, uz, coef)
101 type(coef_t), intent(in), target :: coef
102 real(kind=rp), dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, & coef%msh%nelv), intent(inout) :: res
103 real(kind=rp), dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, & coef%msh%nelv), intent(in) :: ux, uy, uz
104 type(field_t), pointer :: work
105 integer :: ind
106 type(c_ptr) :: res_d
107
108 if (neko_bcknd_device .eq. 1) then
109 res_d = device_get_ptr(res)
110 end if
111
112 call neko_scratch_registry%request_field(work, ind)
113
114 ! Get dux / dx
115 call dudxyz(res, ux, coef%drdx, coef%dsdx, coef%dtdx, coef)
116
117 ! Get duy / dy
118 call dudxyz(work%x, uy, coef%drdy, coef%dsdy, coef%dtdy, coef)
119 if (neko_bcknd_device .eq. 1) then
120 call device_add2(res_d, work%x_d, work%size())
121 else
122 call add2(res, work%x, work%size())
123 end if
124
125 ! Get dux / dz
126 call dudxyz(work%x, uz, coef%drdz, coef%dsdz, coef%dtdz, coef)
127 if (neko_bcknd_device .eq. 1) then
128 call device_add2(res_d, work%x_d, work%size())
129 else
130 call add2(res, work%x, work%size())
131 end if
132
133 call neko_scratch_registry%relinquish_field(ind)
134
135 end subroutine div
136
143 subroutine grad(ux, uy, uz, u, coef)
144 type(coef_t), intent(in) :: coef
145 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: ux
146 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uy
147 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uz
148 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: u
149
150 call dudxyz(ux, u, coef%drdx, coef%dsdx, coef%dtdx, coef)
151 call dudxyz(uy, u, coef%drdy, coef%dsdy, coef%dtdy, coef)
152 call dudxyz(uz, u, coef%drdz, coef%dsdz, coef%dtdz, coef)
153
154 end subroutine grad
155
168 subroutine opgrad(ux, uy, uz, u, coef, es, ee)
169 type(coef_t), intent(in) :: coef
170 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: ux
171 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uy
172 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uz
173 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: u
174 integer, optional :: es, ee
175 integer :: eblk_start, eblk_end
176
177 if (present(es)) then
178 eblk_start = es
179 else
180 eblk_start = 1
181 end if
182
183 if (present(ee)) then
184 eblk_end = ee
185 else
186 eblk_end = coef%msh%nelv
187 end if
188
189 if (neko_bcknd_sx .eq. 1) then
190 call opr_sx_opgrad(ux, uy, uz, u, coef)
191 else if (neko_bcknd_xsmm .eq. 1) then
192 call opr_xsmm_opgrad(ux, uy, uz, u, coef)
193 else if (neko_bcknd_device .eq. 1) then
194 call opr_device_opgrad(ux, uy, uz, u, coef)
195 else
196 call opr_cpu_opgrad(ux, uy, uz, u, coef, eblk_start, eblk_end)
197 end if
198
199 end subroutine opgrad
200
205 subroutine ortho(x, glb_n_points, n)
206 integer, intent(in) :: n
207 integer(kind=i8), intent(in) :: glb_n_points
208 real(kind=rp), dimension(n), intent(inout) :: x
209 real(kind=rp) :: c
210 type(c_ptr) :: x_d
211 if (neko_bcknd_device .eq. 1) then
212 x_d = device_get_ptr(x)
213 c = device_glsum(x_d, n)/glb_n_points
214 call device_cadd(x_d, -c, n)
215 else
216 c = glsum(x, n)/glb_n_points
217 call cadd(x, -c, n)
218 end if
219
220 end subroutine ortho
221
233 subroutine cdtp (dtx, x, dr, ds, dt, coef, es, ee)
234 type(coef_t), intent(in) :: coef
235 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: dtx
236 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: x
237 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: dr
238 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: ds
239 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: dt
240 integer, optional :: es, ee
241 integer :: eblk_start, eblk_end
242
243 if (present(es)) then
244 eblk_start = es
245 else
246 eblk_start = 1
247 end if
248
249 if (present(ee)) then
250 eblk_end = ee
251 else
252 eblk_end = coef%msh%nelv
253 end if
254
255 if (neko_bcknd_sx .eq. 1) then
256 call opr_sx_cdtp(dtx, x, dr, ds, dt, coef)
257 else if (neko_bcknd_xsmm .eq. 1) then
258 call opr_xsmm_cdtp(dtx, x, dr, ds, dt, coef)
259 else if (neko_bcknd_device .eq. 1) then
260 call opr_device_cdtp(dtx, x, dr, ds, dt, coef)
261 else
262 call opr_cpu_cdtp(dtx, x, dr, ds, dt, coef, eblk_start, eblk_end)
263 end if
264
265 end subroutine cdtp
266
277 subroutine conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
278 type(space_t), intent(in) :: xh
279 type(coef_t), intent(in) :: coef
280 real(kind=rp), intent(inout) :: du(xh%lxyz, coef%msh%nelv)
281 real(kind=rp), intent(inout) :: u(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
282 real(kind=rp), intent(inout) :: vx(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
283 real(kind=rp), intent(inout) :: vy(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
284 real(kind=rp), intent(inout) :: vz(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
285 integer, optional :: es, ee
286 integer :: eblk_end, eblk_start
287
288 associate(nelv => coef%msh%nelv, gdim => coef%msh%gdim)
289 if (present(es)) then
290 eblk_start = es
291 else
292 eblk_start = 1
293 end if
294
295 if (present(ee)) then
296 eblk_end = ee
297 else
298 eblk_end = coef%msh%nelv
299 end if
300
301 if (neko_bcknd_sx .eq. 1) then
302 call opr_sx_conv1(du, u, vx, vy, vz, xh, coef, nelv)
303 else if (neko_bcknd_xsmm .eq. 1) then
304 call opr_xsmm_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
305 else if (neko_bcknd_device .eq. 1) then
306 call opr_device_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
307 else
308 call opr_cpu_conv1(du, u, vx, vy, vz, xh, coef, eblk_start, eblk_end)
309 end if
310 end associate
311
312 end subroutine conv1
313
329 subroutine convect_scalar(du, u, c, Xh_GLL, Xh_GL, coef_GLL, &
330 coef_GL, GLL_to_GL)
331 type(space_t), intent(in) :: xh_gl
332 type(space_t), intent(in) :: xh_gll
333 type(coef_t), intent(in) :: coef_GLL
334 type(coef_t), intent(in) :: coef_GL
335 type(interpolator_t), intent(inout) :: GLL_to_GL
336 real(kind=rp), intent(inout) :: &
337 du(xh_gll%lx, xh_gll%ly, xh_gll%lz, coef_gl%msh%nelv)
338 real(kind=rp), intent(inout) :: &
339 u(xh_gl%lx, xh_gl%lx, xh_gl%lx, coef_gl%msh%nelv)
340 real(kind=rp), intent(inout) :: c(xh_gl%lxyz, coef_gl%msh%nelv, 3)
341
342 if (neko_bcknd_sx .eq. 1) then
343 call opr_sx_convect_scalar(du, u, c, xh_gll, xh_gl, &
344 coef_gll, coef_gl, gll_to_gl)
345 else if (neko_bcknd_xsmm .eq. 1) then
346 call opr_xsmm_convect_scalar(du, u, c, xh_gll, xh_gl, &
347 coef_gll, coef_gl, gll_to_gl)
348 else
349 call opr_cpu_convect_scalar(du, u, c, xh_gll, xh_gl, &
350 coef_gll, coef_gl, gll_to_gl)
351 end if
352
353 end subroutine convect_scalar
354
355 !! Compute the curl fo a vector field.
356 !! @param w1 Will store the x component of the curl.
357 !! @param w2 Will store the y component of the curl.
358 !! @param w3 Will store the z component of the curl.
359 !! @param u1 The x component of the vector field.
360 !! @param u2 The y component of the vector field.
361 !! @param u3 The z component of the vector field.
362 !! @param work1 A temporary array for computations.
363 !! @param work2 A temporary array for computations.
364 !! @param coef The SEM coefficients.
365 subroutine curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
366 type(field_t), intent(inout) :: w1
367 type(field_t), intent(inout) :: w2
368 type(field_t), intent(inout) :: w3
369 type(field_t), intent(in) :: u1
370 type(field_t), intent(in) :: u2
371 type(field_t), intent(in) :: u3
372 type(field_t), intent(inout) :: work1
373 type(field_t), intent(inout) :: work2
374 type(coef_t), intent(in) :: coef
375
376 if (neko_bcknd_sx .eq. 1) then
377 call opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
378 else if (neko_bcknd_xsmm .eq. 1) then
379 call opr_xsmm_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
380 else if (neko_bcknd_device .eq. 1) then
381 call opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
382 else
383 call opr_cpu_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
384 end if
385
386 end subroutine curl
387
388 !! Compute the CFL number
389 !! @param dt The timestep.
390 !! @param u The x component of velocity.
391 !! @param v The y component of velocity.
392 !! @param w The z component of velocity.
393 !! @param Xh The SEM function space.
394 !! @param coef The SEM coefficients.
395 !! @param nelv The total number of elements.
396 !! @param gdim Number of geometric dimensions.
397 function cfl(dt, u, v, w, Xh, coef, nelv, gdim)
398 type(space_t), intent(in) :: xh
399 type(coef_t), intent(in) :: coef
400 integer, intent(in) :: nelv, gdim
401 real(kind=rp), intent(in) :: dt
402 real(kind=rp), dimension(Xh%lx, Xh%ly, Xh%lz, nelv), intent(in) :: u, v, w
403 real(kind=rp) :: cfl
404 integer :: ierr
405
406 if (neko_bcknd_sx .eq. 1) then
407 cfl = opr_sx_cfl(dt, u, v, w, xh, coef, nelv)
408 else if (neko_bcknd_device .eq. 1) then
409 cfl = opr_device_cfl(dt, u, v, w, xh, coef, nelv, gdim)
410 else
411 cfl = opr_cpu_cfl(dt, u, v, w, xh, coef, nelv, gdim)
412 end if
413
414 if (.not. neko_device_mpi) then
415 call mpi_allreduce(mpi_in_place, cfl, 1, &
416 mpi_real_precision, mpi_max, neko_comm, ierr)
417 end if
418
419 end function cfl
420
433 subroutine strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
434 type(field_t), intent(in) :: u, v, w
435 type(coef_t), intent(in) :: coef
436 real(kind=rp), intent(inout) :: s11(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
437 real(kind=rp), intent(inout) :: s22(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
438 real(kind=rp), intent(inout) :: s33(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
439 real(kind=rp), intent(inout) :: s12(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
440 real(kind=rp), intent(inout) :: s13(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
441 real(kind=rp), intent(inout) :: s23(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
442
443 type(c_ptr) :: s11_d, s22_d, s33_d, s12_d, s23_d, s13_d
444
445 integer :: nelv, lxyz
446
447 if (neko_bcknd_device .eq. 1) then
448 s11_d = device_get_ptr(s11)
449 s22_d = device_get_ptr(s22)
450 s33_d = device_get_ptr(s33)
451 s12_d = device_get_ptr(s12)
452 s23_d = device_get_ptr(s23)
453 s13_d = device_get_ptr(s13)
454 end if
455
456 nelv = u%msh%nelv
457 lxyz = u%Xh%lxyz
458
459 ! we use s11 as a work array here
460 call dudxyz (s12, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
461 call dudxyz (s11, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
462 if (neko_bcknd_device .eq. 1) then
463 call device_add2(s12_d, s11_d, nelv*lxyz)
464 else
465 call add2(s12, s11, nelv*lxyz)
466 end if
467
468 call dudxyz (s13, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
469 call dudxyz (s11, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
470 if (neko_bcknd_device .eq. 1) then
471 call device_add2(s13_d, s11_d, nelv*lxyz)
472 else
473 call add2(s13, s11, nelv*lxyz)
474 end if
475
476 call dudxyz (s23, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
477 call dudxyz (s11, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
478 if (neko_bcknd_device .eq. 1) then
479 call device_add2(s23_d, s11_d, nelv*lxyz)
480 else
481 call add2(s23, s11, nelv*lxyz)
482 end if
483
484 call dudxyz (s11, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
485 call dudxyz (s22, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
486 call dudxyz (s33, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
487
488 if (neko_bcknd_device .eq. 1) then
489 call device_cmult(s12_d, 0.5_rp, nelv*lxyz)
490 call device_cmult(s13_d, 0.5_rp, nelv*lxyz)
491 call device_cmult(s23_d, 0.5_rp, nelv*lxyz)
492 else
493 call cmult(s12, 0.5_rp, nelv*lxyz)
494 call cmult(s13, 0.5_rp, nelv*lxyz)
495 call cmult(s23, 0.5_rp, nelv*lxyz)
496 end if
497
498 end subroutine strain_rate
499
506 subroutine lambda2op(lambda2, u, v, w, coef)
507 type(coef_t), intent(in) :: coef
508 type(field_t), intent(inout) :: lambda2
509 type(field_t), intent(in) :: u, v, w
511 if (neko_bcknd_sx .eq. 1) then
512 call opr_sx_lambda2(lambda2, u, v, w, coef)
513 else if (neko_bcknd_device .eq. 1) then
514 call opr_device_lambda2(lambda2, u, v, w, coef)
515 else
516 call opr_cpu_lambda2(lambda2, u, v, w, coef)
517 end if
518
519 end subroutine lambda2op
520
531 subroutine set_convect_rst(cr, cs, ct, cx, cy, cz, Xh, coef)
532 type(space_t), intent(inout) :: xh
533 type(coef_t), intent(inout) :: coef
534 real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
535 intent(inout) :: cr, cs, ct
536 real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
537 intent(in) :: cx, cy, cz
538
539 if (neko_bcknd_sx .eq. 1) then
540 call opr_sx_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
541 else if (neko_bcknd_xsmm .eq. 1) then
542 call opr_xsmm_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
543 else
544 call opr_cpu_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
545 end if
546
547 end subroutine set_convect_rst
548
564 subroutine runge_kutta(phi, c_r1, c_r23, c_r4, Xh_GLL, Xh_GL, coef, &
565 coef_GL, GLL_to_GL, tau, dtau, n, nel, n_GL)
566 type(space_t), intent(in) :: xh_gll
567 type(space_t), intent(inout) :: xh_gl
568 type(coef_t), intent(in) :: coef
569 type(coef_t), intent(inout) :: coef_gl
570 type(interpolator_t) :: gll_to_gl
571 real(kind=rp), intent(inout) :: tau, dtau
572 integer, intent(in) :: n, nel, n_gl
573 real(kind=rp), dimension(n), intent(inout) :: phi
574 real(kind=rp), dimension(3 * n_GL), intent(inout) :: c_r1, c_r23, c_r4
575 real(kind=rp) :: c1, c2, c3
576 ! Work Arrays
577 real(kind=rp), dimension(n) :: u1, r1, r2, r3, r4
578 real(kind=rp), dimension(n_GL) :: u1_gl
579 integer :: i, e
580
581 c1 = 1.
582 c2 = -dtau/2.
583 c3 = -dtau
584
585 ! Stage 1:
586 call invcol3 (u1, phi, coef%B, n)
587 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
588 call convect_scalar(r1, u1_gl, c_r1, xh_gll, xh_gl, coef, &
589 coef_gl, gll_to_gl)
590 call col2(r1, coef%B, n)
591
592 ! Stage 2:
593 call add3s2 (u1, phi, r1, c1, c2, n)
594 call invcol2 (u1, coef%B, n)
595 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
596 call convect_scalar(r2, u1_gl, c_r23, xh_gll, xh_gl, coef, &
597 coef_gl, gll_to_gl)
598 call col2(r2, coef%B, n)
599
600 ! Stage 3:
601 call add3s2 (u1, phi, r2, c1, c2, n)
602 call invcol2 (u1, coef%B, n)
603 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
604 call convect_scalar(r3, u1_gl, c_r23, xh_gll, xh_gl, coef, &
605 coef_gl, gll_to_gl)
606 call col2(r3, coef%B, n)
607
608 ! Stage 4:
609 call add3s2 (u1, phi, r3, c1, c3, n)
610 call invcol2 (u1, coef%B, n)
611 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
612 call convect_scalar(r4, u1_gl, c_r4, xh_gll, xh_gl, coef, &
613 coef_gl, gll_to_gl)
614 call col2(r4, coef%B, n)
615
616 c1 = -dtau/6.
617 c2 = -dtau/3.
618 do i = 1, n
619 phi(i) = phi(i) + c1 * (r1(i) + r4(i)) + c2 * (r2(i) + r3(i))
620 end do
621
622 end subroutine runge_kutta
623
624end module operators
625
Return the device pointer for an associated Fortran array.
Definition device.F90:92
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
subroutine, public device_add2(a_d, b_d, n)
Vector addition .
subroutine, public device_cadd(a_d, c, n)
Add a scalar to vector .
subroutine, public device_cmult(a_d, c, n)
Multiplication by constant c .
real(kind=rp) function, public device_glsum(a_d, n)
Sum a vector of length n.
subroutine, public device_copy(a_d, b_d, n)
Copy a vector .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
Defines a field.
Definition field.f90:34
Routines to interpolate between different spaces.
A simulation component that computes lambda2 The values are stored in the field registry under the na...
Definition lambda2.f90:37
Definition math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition math.f90:310
subroutine, public invcol2(a, b, n)
Vector division .
Definition math.f90:714
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:322
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:359
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:586
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition math.f90:486
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition math.f90:769
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:728
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:238
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:194
Build configurations.
integer, parameter neko_bcknd_sx
integer, parameter neko_bcknd_device
logical, parameter neko_device_mpi
integer, parameter neko_bcknd_xsmm
integer, parameter, public i8
Definition num_types.f90:7
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public runge_kutta(phi, c_r1, c_r23, c_r4, xh_gll, xh_gl, coef, coef_gl, gll_to_gl, tau, dtau, n, nel, n_gl)
Compute one step of Runge Kutta time interpolation for OIFS scheme.
subroutine, public set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
Transforms the convecting velocity field to the rst form of the GL space.
subroutine, public ortho(x, glb_n_points, n)
Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T.
subroutine convect_scalar(du, u, c, xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
Apply the convecting velocity c to the to the scalar field u, used in the OIFS scheme.
real(kind=rp) function, public cfl(dt, u, v, w, xh, coef, nelv, gdim)
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.
subroutine, public div(res, ux, uy, uz, coef)
Compute the divergence of a vector field.
subroutine, public conv1(du, u, vx, vy, vz, xh, coef, es, ee)
Compute the advection term.
subroutine, public grad(ux, uy, uz, u, coef)
Compute the gradient of a scalar field.
subroutine, public curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
subroutine, public strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
Compute the strain rate tensor, i.e 0.5 * du_i/dx_j + du_j/dx_i.
subroutine, public lambda2op(lambda2, u, v, w, coef)
Compute the Lambda2 field for a given velocity field.
subroutine, public cdtp(dtx, x, dr, ds, dt, coef, es, ee)
Apply D^T to a scalar field, where D is the derivative matrix.
subroutine, public dudxyz(du, u, dr, ds, dt, coef)
Compute derivative of a scalar field along a single direction.
Definition operators.f90:78
Operators CPU backend.
Definition opr_cpu.f90:34
subroutine, public opr_cpu_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_xh)
Definition opr_cpu.f90:124
real(kind=rp) function, public opr_cpu_cfl(dt, u, v, w, xh, coef, nelv, gdim)
Definition opr_cpu.f90:175
subroutine, public opr_cpu_lambda2(lambda2, u, v, w, coef)
Definition opr_cpu.f90:235
Operators accelerator backends.
subroutine, public opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_xh)
subroutine, public opr_device_cdtp(dtx, x, dr, ds, dt, coef)
real(kind=rp) function, public opr_device_cfl(dt, u, v, w, xh, coef, nelv, gdim)
subroutine, public opr_device_dudxyz(du, u, dr, ds, dt, coef)
subroutine, public opr_device_opgrad(ux, uy, uz, u, coef)
subroutine, public opr_device_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
subroutine, public opr_device_lambda2(lambda2, u, v, w, coef)
Operators SX-Aurora backend.
Definition opr_sx.f90:2
subroutine, public opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_xh)
Definition opr_sx.f90:100
Operators libxsmm backend.
Definition opr_xsmm.F90:61
subroutine, public opr_xsmm_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
Definition opr_xsmm.F90:297
subroutine, public opr_xsmm_cdtp(dtx, x, dr, ds, dt, coef)
Definition opr_xsmm.F90:242
subroutine, public opr_xsmm_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_xh)
Definition opr_xsmm.F90:416
subroutine, public opr_xsmm_dudxyz(du, u, dr, ds, dt, coef)
Definition opr_xsmm.F90:92
subroutine, public opr_xsmm_opgrad(ux, uy, uz, u, coef)
Definition opr_xsmm.F90:146
subroutine, public opr_xsmm_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
Definition opr_xsmm.F90:467
subroutine, public opr_xsmm_convect_scalar(du, u, c, xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
Definition opr_xsmm.F90:371
Defines a registry for storing and requesting temporary fields This can be used when you have a funct...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Defines a function space.
Definition space.f90:34
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Interpolation between two space::space_t.
The function space for the SEM solution fields.
Definition space.f90:62