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 : c_ptr, device_get_ptr
59 use comm
60 implicit none
61 private
62
63 public :: dudxyz, opgrad, ortho, cdtp, conv1, curl, cfl, &
65
66contains
67
75 subroutine dudxyz (du, u, dr, ds, dt, coef)
76 type(coef_t), intent(in), target :: coef
77 real(kind=rp), dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, & coef%msh%nelv), intent(inout) :: du
78 real(kind=rp), dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, & coef%msh%nelv), intent(in) :: u, dr, ds, dt
79
80 if (neko_bcknd_sx .eq. 1) then
81 call opr_sx_dudxyz(du, u, dr, ds, dt, coef)
82 else if (neko_bcknd_xsmm .eq. 1) then
83 call opr_xsmm_dudxyz(du, u, dr, ds, dt, coef)
84 else if (neko_bcknd_device .eq. 1) then
85 call opr_device_dudxyz(du, u, dr, ds, dt, coef)
86 else
87 call opr_cpu_dudxyz(du, u, dr, ds, dt, coef)
88 end if
89
90 end subroutine dudxyz
91
98 subroutine div(res, ux, uy, uz, coef)
99 type(coef_t), intent(in), target :: coef
100 real(kind=rp), dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, & coef%msh%nelv), intent(inout) :: res
101 real(kind=rp), dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, & coef%msh%nelv), intent(in) :: ux, uy, uz
102 type(field_t), pointer :: work
103 integer :: ind
104 type(c_ptr) :: res_d
105
106 if (neko_bcknd_device .eq. 1) then
107 res_d = device_get_ptr(res)
108 end if
109
110 call neko_scratch_registry%request_field(work, ind)
111
112 ! Get dux / dx
113 call dudxyz(res, ux, coef%drdx, coef%dsdx, coef%dtdx, coef)
114
115 ! Get duy / dy
116 call dudxyz(work%x, uy, coef%drdy, coef%dsdy, coef%dtdy, coef)
117 if (neko_bcknd_device .eq. 1) then
118 call device_add2(res_d, work%x_d, work%size())
119 else
120 call add2(res, work%x, work%size())
121 end if
122
123 ! Get dux / dz
124 call dudxyz(work%x, uz, coef%drdz, coef%dsdz, coef%dtdz, coef)
125 if (neko_bcknd_device .eq. 1) then
126 call device_add2(res_d, work%x_d, work%size())
127 else
128 call add2(res, work%x, work%size())
129 end if
130
131 call neko_scratch_registry%relinquish_field(ind)
132
133 end subroutine div
134
141 subroutine grad(ux, uy, uz, u, coef)
142 type(coef_t), intent(in) :: coef
143 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: ux
144 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uy
145 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uz
146 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: u
147
148 call dudxyz(ux, u, coef%drdx, coef%dsdx, coef%dtdx, coef)
149 call dudxyz(uy, u, coef%drdy, coef%dsdy, coef%dtdy, coef)
150 call dudxyz(uz, u, coef%drdz, coef%dsdz, coef%dtdz, coef)
151
152 end subroutine grad
153
166 subroutine opgrad(ux, uy, uz, u, coef, es, ee)
167 type(coef_t), intent(in) :: coef
168 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: ux
169 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uy
170 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uz
171 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: u
172 integer, optional :: es, ee
173 integer :: eblk_start, eblk_end
174
175 if (present(es)) then
176 eblk_start = es
177 else
178 eblk_start = 1
179 end if
180
181 if (present(ee)) then
182 eblk_end = ee
183 else
184 eblk_end = coef%msh%nelv
185 end if
186
187 if (neko_bcknd_sx .eq. 1) then
188 call opr_sx_opgrad(ux, uy, uz, u, coef)
189 else if (neko_bcknd_xsmm .eq. 1) then
190 call opr_xsmm_opgrad(ux, uy, uz, u, coef)
191 else if (neko_bcknd_device .eq. 1) then
192 call opr_device_opgrad(ux, uy, uz, u, coef)
193 else
194 call opr_cpu_opgrad(ux, uy, uz, u, coef, eblk_start, eblk_end)
195 end if
196
197 end subroutine opgrad
198
203 subroutine ortho(x, glb_n_points, n)
204 integer, intent(in) :: n
205 integer(kind=i8), intent(in) :: glb_n_points
206 real(kind=rp), dimension(n), intent(inout) :: x
207 real(kind=rp) :: c
208 type(c_ptr) :: x_d
209 if (neko_bcknd_device .eq. 1) then
210 x_d = device_get_ptr(x)
211 c = device_glsum(x_d, n)/glb_n_points
212 call device_cadd(x_d, -c, n)
213 else
214 c = glsum(x, n)/glb_n_points
215 call cadd(x, -c, n)
216 end if
217
218 end subroutine ortho
219
231 subroutine cdtp (dtx, x, dr, ds, dt, coef, es, ee)
232 type(coef_t), intent(in) :: coef
233 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: dtx
234 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: x
235 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: dr
236 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: ds
237 real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: dt
238 integer, optional :: es, ee
239 integer :: eblk_start, eblk_end
240
241 if (present(es)) then
242 eblk_start = es
243 else
244 eblk_start = 1
245 end if
246
247 if (present(ee)) then
248 eblk_end = ee
249 else
250 eblk_end = coef%msh%nelv
251 end if
252
253 if (neko_bcknd_sx .eq. 1) then
254 call opr_sx_cdtp(dtx, x, dr, ds, dt, coef)
255 else if (neko_bcknd_xsmm .eq. 1) then
256 call opr_xsmm_cdtp(dtx, x, dr, ds, dt, coef)
257 else if (neko_bcknd_device .eq. 1) then
258 call opr_device_cdtp(dtx, x, dr, ds, dt, coef)
259 else
260 call opr_cpu_cdtp(dtx, x, dr, ds, dt, coef, eblk_start, eblk_end)
261 end if
262
263 end subroutine cdtp
264
275 subroutine conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
276 type(space_t), intent(in) :: xh
277 type(coef_t), intent(in) :: coef
278 real(kind=rp), intent(inout) :: du(xh%lxyz, coef%msh%nelv)
279 real(kind=rp), intent(inout) :: u(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
280 real(kind=rp), intent(inout) :: vx(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
281 real(kind=rp), intent(inout) :: vy(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
282 real(kind=rp), intent(inout) :: vz(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
283 integer, optional :: es, ee
284 integer :: eblk_end, eblk_start
285
286 associate(nelv => coef%msh%nelv, gdim => coef%msh%gdim)
287 if (present(es)) then
288 eblk_start = es
289 else
290 eblk_start = 1
291 end if
292
293 if (present(ee)) then
294 eblk_end = ee
295 else
296 eblk_end = coef%msh%nelv
297 end if
298
299 if (neko_bcknd_sx .eq. 1) then
300 call opr_sx_conv1(du, u, vx, vy, vz, xh, coef, nelv)
301 else if (neko_bcknd_xsmm .eq. 1) then
302 call opr_xsmm_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
303 else if (neko_bcknd_device .eq. 1) then
304 call opr_device_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
305 else
306 call opr_cpu_conv1(du, u, vx, vy, vz, xh, coef, eblk_start, eblk_end)
307 end if
308 end associate
309
310 end subroutine conv1
311
327 subroutine convect_scalar(du, u, c, Xh_GLL, Xh_GL, coef_GLL, &
328 coef_GL, GLL_to_GL)
329 type(space_t), intent(in) :: xh_gl
330 type(space_t), intent(in) :: xh_gll
331 type(coef_t), intent(in) :: coef_GLL
332 type(coef_t), intent(in) :: coef_GL
333 type(interpolator_t), intent(inout) :: GLL_to_GL
334 real(kind=rp), intent(inout) :: &
335 du(xh_gll%lx, xh_gll%ly, xh_gll%lz, coef_gl%msh%nelv)
336 real(kind=rp), intent(inout) :: &
337 u(xh_gl%lx, xh_gl%lx, xh_gl%lx, coef_gl%msh%nelv)
338 real(kind=rp), intent(inout) :: c(xh_gl%lxyz, coef_gl%msh%nelv, 3)
339
340 if (neko_bcknd_sx .eq. 1) then
341 call opr_sx_convect_scalar(du, u, c, xh_gll, xh_gl, &
342 coef_gll, coef_gl, gll_to_gl)
343 else if (neko_bcknd_xsmm .eq. 1) then
344 call opr_xsmm_convect_scalar(du, u, c, xh_gll, xh_gl, &
345 coef_gll, coef_gl, gll_to_gl)
346 else
347 call opr_cpu_convect_scalar(du, u, c, xh_gll, xh_gl, &
348 coef_gll, coef_gl, gll_to_gl)
349 end if
350
351 end subroutine convect_scalar
352
353 !! Compute the curl fo a vector field.
354 !! @param w1 Will store the x component of the curl.
355 !! @param w2 Will store the y component of the curl.
356 !! @param w3 Will store the z component of the curl.
357 !! @param u1 The x component of the vector field.
358 !! @param u2 The y component of the vector field.
359 !! @param u3 The z component of the vector field.
360 !! @param work1 A temporary array for computations.
361 !! @param work2 A temporary array for computations.
362 !! @param coef The SEM coefficients.
363 subroutine curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
364 type(field_t), intent(inout) :: w1
365 type(field_t), intent(inout) :: w2
366 type(field_t), intent(inout) :: w3
367 type(field_t), intent(in) :: u1
368 type(field_t), intent(in) :: u2
369 type(field_t), intent(in) :: u3
370 type(field_t), intent(inout) :: work1
371 type(field_t), intent(inout) :: work2
372 type(coef_t), intent(in) :: coef
373
374 if (neko_bcknd_sx .eq. 1) then
375 call opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
376 else if (neko_bcknd_xsmm .eq. 1) then
377 call opr_xsmm_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
378 else if (neko_bcknd_device .eq. 1) then
379 call opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
380 else
381 call opr_cpu_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
382 end if
383
384 end subroutine curl
385
386 !! Compute the CFL number
387 !! @param dt The timestep.
388 !! @param u The x component of velocity.
389 !! @param v The y component of velocity.
390 !! @param w The z component of velocity.
391 !! @param Xh The SEM function space.
392 !! @param coef The SEM coefficients.
393 !! @param nelv The total number of elements.
394 !! @param gdim Number of geometric dimensions.
395 function cfl(dt, u, v, w, Xh, coef, nelv, gdim)
396 type(space_t), intent(in) :: xh
397 type(coef_t), intent(in) :: coef
398 integer, intent(in) :: nelv, gdim
399 real(kind=rp), intent(in) :: dt
400 real(kind=rp), dimension(Xh%lx, Xh%ly, Xh%lz, nelv), intent(in) :: u, v, w
401 real(kind=rp) :: cfl
402 integer :: ierr
403
404 if (neko_bcknd_sx .eq. 1) then
405 cfl = opr_sx_cfl(dt, u, v, w, xh, coef, nelv)
406 else if (neko_bcknd_device .eq. 1) then
407 cfl = opr_device_cfl(dt, u, v, w, xh, coef, nelv, gdim)
408 else
409 cfl = opr_cpu_cfl(dt, u, v, w, xh, coef, nelv, gdim)
410 end if
411
412 if (.not. neko_device_mpi) then
413 call mpi_allreduce(mpi_in_place, cfl, 1, &
414 mpi_real_precision, mpi_max, neko_comm, ierr)
415 end if
416
417 end function cfl
418
431 subroutine strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
432 type(field_t), intent(in) :: u, v, w
433 type(coef_t), intent(in) :: coef
434 real(kind=rp), intent(inout) :: s11(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
435 real(kind=rp), intent(inout) :: s22(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
436 real(kind=rp), intent(inout) :: s33(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
437 real(kind=rp), intent(inout) :: s12(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
438 real(kind=rp), intent(inout) :: s13(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
439 real(kind=rp), intent(inout) :: s23(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
440
441 type(c_ptr) :: s11_d, s22_d, s33_d, s12_d, s23_d, s13_d
442
443 integer :: nelv, lxyz
444
445 if (neko_bcknd_device .eq. 1) then
446 s11_d = device_get_ptr(s11)
447 s22_d = device_get_ptr(s22)
448 s33_d = device_get_ptr(s33)
449 s12_d = device_get_ptr(s12)
450 s23_d = device_get_ptr(s23)
451 s13_d = device_get_ptr(s13)
452 end if
453
454 nelv = u%msh%nelv
455 lxyz = u%Xh%lxyz
456
457 ! we use s11 as a work array here
458 call dudxyz (s12, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
459 call dudxyz (s11, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
460 if (neko_bcknd_device .eq. 1) then
461 call device_add2(s12_d, s11_d, nelv*lxyz)
462 else
463 call add2(s12, s11, nelv*lxyz)
464 end if
465
466 call dudxyz (s13, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
467 call dudxyz (s11, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
468 if (neko_bcknd_device .eq. 1) then
469 call device_add2(s13_d, s11_d, nelv*lxyz)
470 else
471 call add2(s13, s11, nelv*lxyz)
472 end if
473
474 call dudxyz (s23, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
475 call dudxyz (s11, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
476 if (neko_bcknd_device .eq. 1) then
477 call device_add2(s23_d, s11_d, nelv*lxyz)
478 else
479 call add2(s23, s11, nelv*lxyz)
480 end if
481
482 call dudxyz (s11, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
483 call dudxyz (s22, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
484 call dudxyz (s33, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
485
486 if (neko_bcknd_device .eq. 1) then
487 call device_cmult(s12_d, 0.5_rp, nelv*lxyz)
488 call device_cmult(s13_d, 0.5_rp, nelv*lxyz)
489 call device_cmult(s23_d, 0.5_rp, nelv*lxyz)
490 else
491 call cmult(s12, 0.5_rp, nelv*lxyz)
492 call cmult(s13, 0.5_rp, nelv*lxyz)
493 call cmult(s23, 0.5_rp, nelv*lxyz)
494 end if
495
496 end subroutine strain_rate
497
504 subroutine lambda2op(lambda2, u, v, w, coef)
505 type(coef_t), intent(in) :: coef
506 type(field_t), intent(inout) :: lambda2
507 type(field_t), intent(in) :: u, v, w
509 if (neko_bcknd_sx .eq. 1) then
510 call opr_sx_lambda2(lambda2, u, v, w, coef)
511 else if (neko_bcknd_device .eq. 1) then
512 call opr_device_lambda2(lambda2, u, v, w, coef)
513 else
514 call opr_cpu_lambda2(lambda2, u, v, w, coef)
515 end if
516
517 end subroutine lambda2op
518
529 subroutine set_convect_rst(cr, cs, ct, cx, cy, cz, Xh, coef)
530 type(space_t), intent(inout) :: xh
531 type(coef_t), intent(inout) :: coef
532 real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
533 intent(inout) :: cr, cs, ct
534 real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
535 intent(in) :: cx, cy, cz
536
537 if (neko_bcknd_sx .eq. 1) then
538 call opr_sx_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
539 else if (neko_bcknd_xsmm .eq. 1) then
540 call opr_xsmm_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
541 else
542 call opr_cpu_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
543 end if
544
545 end subroutine set_convect_rst
546
562 subroutine runge_kutta(phi, c_r1, c_r23, c_r4, Xh_GLL, Xh_GL, coef, &
563 coef_GL, GLL_to_GL, tau, dtau, n, nel, n_GL)
564 type(space_t), intent(in) :: xh_gll
565 type(space_t), intent(inout) :: xh_gl
566 type(coef_t), intent(in) :: coef
567 type(coef_t), intent(inout) :: coef_gl
568 type(interpolator_t) :: gll_to_gl
569 real(kind=rp), intent(inout) :: tau, dtau
570 integer, intent(in) :: n, nel, n_gl
571 real(kind=rp), dimension(n), intent(inout) :: phi
572 real(kind=rp), dimension(3 * n_GL), intent(inout) :: c_r1, c_r23, c_r4
573 real(kind=rp) :: c1, c2, c3
574 ! Work Arrays
575 real(kind=rp), dimension(n) :: u1, r1, r2, r3, r4
576 real(kind=rp), dimension(n_GL) :: u1_gl
577 integer :: i, e
578
579 c1 = 1.
580 c2 = -dtau/2.
581 c3 = -dtau
582
583 ! Stage 1:
584 call invcol3 (u1, phi, coef%B, n)
585 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
586 call convect_scalar(r1, u1_gl, c_r1, xh_gll, xh_gl, coef, &
587 coef_gl, gll_to_gl)
588 call col2(r1, coef%B, n)
589
590 ! Stage 2:
591 call add3s2 (u1, phi, r1, c1, c2, n)
592 call invcol2 (u1, coef%B, n)
593 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
594 call convect_scalar(r2, u1_gl, c_r23, xh_gll, xh_gl, coef, &
595 coef_gl, gll_to_gl)
596 call col2(r2, coef%B, n)
597
598 ! Stage 3:
599 call add3s2 (u1, phi, r2, c1, c2, n)
600 call invcol2 (u1, coef%B, n)
601 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
602 call convect_scalar(r3, u1_gl, c_r23, xh_gll, xh_gl, coef, &
603 coef_gl, gll_to_gl)
604 call col2(r3, coef%B, n)
605
606 ! Stage 4:
607 call add3s2 (u1, phi, r3, c1, c3, n)
608 call invcol2 (u1, coef%B, n)
609 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
610 call convect_scalar(r4, u1_gl, c_r4, xh_gll, xh_gl, coef, &
611 coef_gl, gll_to_gl)
612 call col2(r4, coef%B, n)
613
614 c1 = -dtau/6.
615 c2 = -dtau/3.
616 do i = 1, n
617 phi(i) = phi(i) + c1 * (r1(i) + r4(i)) + c2 * (r2(i) + r3(i))
618 end do
619
620 end subroutine runge_kutta
621
622end module operators
623
Return the device pointer for an associated Fortran array.
Definition device.F90:81
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:76
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