Neko  0.8.99
A portable framework for high-order spectral element flow simulations
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 !
34 module operators
37  use num_types, only : rp
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
53  use interpolation, only : interpolator_t
54  use math, only : glsum, cmult, add2, add3s2, cadd, copy, col2, invcol2, &
55  invcol3, rzero
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 
66 contains
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  res_d = device_get_ptr(res)
107 
108  call neko_scratch_registry%request_field(work, ind)
109 
110  ! Get dux / dx
111  call dudxyz(res, ux, coef%drdx, coef%dsdx, coef%dtdx, coef)
112 
113  ! Get duy / dy
114  call dudxyz(work%x, uy, coef%drdy, coef%dsdy, coef%dtdy, coef)
115  if (neko_bcknd_device .eq. 1) then
116  call device_add2(res_d, work%x_d, work%size())
117  else
118  call add2(res, work%x, work%size())
119  end if
120 
121  ! Get dux / dz
122  call dudxyz(work%x, uz, coef%drdz, coef%dsdz, coef%dtdz, coef)
123  if (neko_bcknd_device .eq. 1) then
124  call device_add2(res_d, work%x_d, work%size())
125  else
126  call add2(res, work%x, work%size())
127  end if
128 
129  call neko_scratch_registry%relinquish_field(ind)
130 
131  end subroutine div
132 
139  subroutine grad(ux, uy, uz, u, coef)
140  type(coef_t), intent(in) :: coef
141  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: ux
142  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uy
143  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uz
144  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: u
145 
146  call dudxyz(ux, u, coef%drdx, coef%dsdx, coef%dtdx, coef)
147  call dudxyz(uy, u, coef%drdy, coef%dsdy, coef%dtdy, coef)
148  call dudxyz(uz, u, coef%drdz, coef%dsdz, coef%dtdz, coef)
149 
150  end subroutine grad
151 
164  subroutine opgrad(ux, uy, uz, u, coef, es, ee)
165  type(coef_t), intent(in) :: coef
166  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: ux
167  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uy
168  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: uz
169  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: u
170  integer, optional :: es, ee
171  integer :: eblk_start, eblk_end
172 
173  if (present(es)) then
174  eblk_start = es
175  else
176  eblk_start = 1
177  end if
178 
179  if (present(ee)) then
180  eblk_end = ee
181  else
182  eblk_end = coef%msh%nelv
183  end if
184 
185  if (neko_bcknd_sx .eq. 1) then
186  call opr_sx_opgrad(ux, uy, uz, u, coef)
187  else if (neko_bcknd_xsmm .eq. 1) then
188  call opr_xsmm_opgrad(ux, uy, uz, u, coef)
189  else if (neko_bcknd_device .eq. 1) then
190  call opr_device_opgrad(ux, uy, uz, u, coef)
191  else
192  call opr_cpu_opgrad(ux, uy, uz, u, coef, eblk_start, eblk_end)
193  end if
194 
195  end subroutine opgrad
196 
201  subroutine ortho(x, n, glb_n)
202  integer, intent(in) :: n
203  integer, intent(in) :: glb_n
204  real(kind=rp), dimension(n), intent(inout) :: x
205  real(kind=rp) :: rlam
206 
207  rlam = glsum(x, n)/glb_n
208  call cadd(x, -rlam, n)
209 
210  end subroutine ortho
211 
223  subroutine cdtp (dtx, x, dr, ds, dt, coef, es, ee)
224  type(coef_t), intent(in) :: coef
225  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: dtx
226  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: x
227  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: dr
228  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: ds
229  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: dt
230  integer, optional :: es, ee
231  integer :: eblk_start, eblk_end
232 
233  if (present(es)) then
234  eblk_start = es
235  else
236  eblk_start = 1
237  end if
238 
239  if (present(ee)) then
240  eblk_end = ee
241  else
242  eblk_end = coef%msh%nelv
243  end if
244 
245  if (neko_bcknd_sx .eq. 1) then
246  call opr_sx_cdtp(dtx, x, dr, ds, dt, coef)
247  else if (neko_bcknd_xsmm .eq. 1) then
248  call opr_xsmm_cdtp(dtx, x, dr, ds, dt, coef)
249  else if (neko_bcknd_device .eq. 1) then
250  call opr_device_cdtp(dtx, x, dr, ds, dt, coef)
251  else
252  call opr_cpu_cdtp(dtx, x, dr, ds, dt, coef, eblk_start, eblk_end)
253  end if
254 
255  end subroutine cdtp
256 
267  subroutine conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
268  type(space_t), intent(inout) :: xh
269  type(coef_t), intent(inout) :: coef
270  real(kind=rp), intent(inout) :: du(xh%lxyz, coef%msh%nelv)
271  real(kind=rp), intent(inout) :: u(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
272  real(kind=rp), intent(inout) :: vx(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
273  real(kind=rp), intent(inout) :: vy(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
274  real(kind=rp), intent(inout) :: vz(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
275  integer, optional :: es, ee
276  integer :: eblk_end, eblk_start
277 
278  associate(nelv => coef%msh%nelv, gdim => coef%msh%gdim)
279  if (present(es)) then
280  eblk_start = es
281  else
282  eblk_start = 1
283  end if
284 
285  if (present(ee)) then
286  eblk_end = ee
287  else
288  eblk_end = coef%msh%nelv
289  end if
290 
291  if (neko_bcknd_sx .eq. 1) then
292  call opr_sx_conv1(du, u, vx, vy, vz, xh, coef, nelv)
293  else if (neko_bcknd_xsmm .eq. 1) then
294  call opr_xsmm_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
295  else if (neko_bcknd_device .eq. 1) then
296  call opr_device_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
297  else
298  call opr_cpu_conv1(du, u, vx, vy, vz, xh, coef, eblk_start, eblk_end)
299  end if
300  end associate
301 
302  end subroutine conv1
303 
319  subroutine convect_scalar(du, u, c, Xh_GLL, Xh_GL, coef_GLL, &
320  coef_GL, GLL_to_GL)
321  type(space_t), intent(in) :: xh_gl
322  type(space_t), intent(in) :: xh_gll
323  type(coef_t), intent(in) :: coef_GLL
324  type(coef_t), intent(in) :: coef_GL
325  type(interpolator_t), intent(inout) :: GLL_to_GL
326  real(kind=rp), intent(inout) :: &
327  du(xh_gll%lx, xh_gll%ly, xh_gll%lz, coef_gl%msh%nelv)
328  real(kind=rp), intent(inout) :: &
329  u(xh_gl%lx, xh_gl%lx, xh_gl%lx, coef_gl%msh%nelv)
330  real(kind=rp), intent(inout) :: c(xh_gl%lxyz, coef_gl%msh%nelv, 3)
331 
332  if (neko_bcknd_sx .eq. 1) then
333  call opr_sx_convect_scalar(du, u, c, xh_gll, xh_gl, &
334  coef_gll, coef_gl, gll_to_gl)
335  else if (neko_bcknd_xsmm .eq. 1) then
336  call opr_xsmm_convect_scalar(du, u, c, xh_gll, xh_gl, &
337  coef_gll, coef_gl, gll_to_gl)
338  else
339  call opr_cpu_convect_scalar(du, u, c, xh_gll, xh_gl, &
340  coef_gll, coef_gl, gll_to_gl)
341  end if
342 
343  end subroutine convect_scalar
344 
345  !! Compute the curl fo a vector field.
346  !! @param w1 Will store the x component of the curl.
347  !! @param w2 Will store the y component of the curl.
348  !! @param w3 Will store the z component of the curl.
349  !! @param u1 The x component of the vector field.
350  !! @param u2 The y component of the vector field.
351  !! @param u3 The z component of the vector field.
352  !! @param work1 A temporary array for computations.
353  !! @param work2 A temporary array for computations.
354  !! @param coef The SEM coefficients.
355  subroutine curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
356  type(field_t), intent(inout) :: w1
357  type(field_t), intent(inout) :: w2
358  type(field_t), intent(inout) :: w3
359  type(field_t), intent(inout) :: u1
360  type(field_t), intent(inout) :: u2
361  type(field_t), intent(inout) :: u3
362  type(field_t), intent(inout) :: work1
363  type(field_t), intent(inout) :: work2
364  type(coef_t), intent(in) :: coef
365 
366  if (neko_bcknd_sx .eq. 1) then
367  call opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
368  else if (neko_bcknd_xsmm .eq. 1) then
369  call opr_xsmm_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
370  else if (neko_bcknd_device .eq. 1) then
371  call opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
372  else
373  call opr_cpu_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
374  end if
375 
376  end subroutine curl
377 
378  !! Compute the CFL number
379  !! @param dt The timestep.
380  !! @param u The x component of velocity.
381  !! @param v The y component of velocity.
382  !! @param w The z component of velocity.
383  !! @param Xh The SEM function space.
384  !! @param coef The SEM coefficients.
385  !! @param nelv The total number of elements.
386  !! @param gdim Number of geometric dimensions.
387  function cfl(dt, u, v, w, Xh, coef, nelv, gdim)
388  type(space_t), intent(in) :: xh
389  type(coef_t), intent(in) :: coef
390  integer, intent(in) :: nelv, gdim
391  real(kind=rp), intent(in) :: dt
392  real(kind=rp), dimension(Xh%lx, Xh%ly, Xh%lz, nelv), intent(in) :: u, v, w
393  real(kind=rp) :: cfl
394  integer :: ierr
395 
396  if (neko_bcknd_sx .eq. 1) then
397  cfl = opr_sx_cfl(dt, u, v, w, xh, coef, nelv)
398  else if (neko_bcknd_device .eq. 1) then
399  cfl = opr_device_cfl(dt, u, v, w, xh, coef, nelv, gdim)
400  else
401  cfl = opr_cpu_cfl(dt, u, v, w, xh, coef, nelv, gdim)
402  end if
403 
404  if (.not. neko_device_mpi) then
405  call mpi_allreduce(mpi_in_place, cfl, 1, &
406  mpi_real_precision, mpi_max, neko_comm, ierr)
407  end if
408 
409  end function cfl
410 
423  subroutine strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
424  type(field_t), intent(in) :: u, v, w
425  type(coef_t), intent(in) :: coef
426  real(kind=rp), intent(inout) :: s11(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
427  real(kind=rp), intent(inout) :: s22(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
428  real(kind=rp), intent(inout) :: s33(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
429  real(kind=rp), intent(inout) :: s12(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
430  real(kind=rp), intent(inout) :: s13(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
431  real(kind=rp), intent(inout) :: s23(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
432 
433  type(c_ptr) :: s11_d, s22_d, s33_d, s12_d, s23_d, s13_d
434 
435  integer :: nelv, lxyz
436 
437  if (neko_bcknd_device .eq. 1) then
438  s11_d = device_get_ptr(s11)
439  s22_d = device_get_ptr(s22)
440  s33_d = device_get_ptr(s33)
441  s12_d = device_get_ptr(s12)
442  s23_d = device_get_ptr(s23)
443  s13_d = device_get_ptr(s13)
444  end if
445 
446  nelv = u%msh%nelv
447  lxyz = u%Xh%lxyz
448 
449  ! we use s11 as a work array here
450  call dudxyz (s12, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
451  call dudxyz (s11, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
452  if (neko_bcknd_device .eq. 1) then
453  call device_add2(s12_d, s11_d, nelv*lxyz)
454  else
455  call add2(s12, s11, nelv*lxyz)
456  end if
457 
458  call dudxyz (s13, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
459  call dudxyz (s11, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
460  if (neko_bcknd_device .eq. 1) then
461  call device_add2(s13_d, s11_d, nelv*lxyz)
462  else
463  call add2(s13, s11, nelv*lxyz)
464  end if
465 
466  call dudxyz (s23, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
467  call dudxyz (s11, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
468  if (neko_bcknd_device .eq. 1) then
469  call device_add2(s23_d, s11_d, nelv*lxyz)
470  else
471  call add2(s23, s11, nelv*lxyz)
472  end if
473 
474  call dudxyz (s11, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
475  call dudxyz (s22, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
476  call dudxyz (s33, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
477 
478  if (neko_bcknd_device .eq. 1) then
479  call device_cmult(s12_d, 0.5_rp, nelv*lxyz)
480  call device_cmult(s13_d, 0.5_rp, nelv*lxyz)
481  call device_cmult(s23_d, 0.5_rp, nelv*lxyz)
482  else
483  call cmult(s12, 0.5_rp, nelv*lxyz)
484  call cmult(s13, 0.5_rp, nelv*lxyz)
485  call cmult(s23, 0.5_rp, nelv*lxyz)
486  end if
487 
488  end subroutine strain_rate
489 
496  subroutine lambda2op(lambda2, u, v, w, coef)
497  type(coef_t), intent(in) :: coef
498  type(field_t), intent(inout) :: lambda2
499  type(field_t), intent(in) :: u, v, w
500 
501  if (neko_bcknd_sx .eq. 1) then
502  call opr_sx_lambda2(lambda2, u, v, w, coef)
503  else if (neko_bcknd_device .eq. 1) then
504  call opr_device_lambda2(lambda2, u, v, w, coef)
505  else
506  call opr_cpu_lambda2(lambda2, u, v, w, coef)
507  end if
508 
509  end subroutine lambda2op
510 
521  subroutine set_convect_rst(cr, cs, ct, cx, cy, cz, Xh, coef)
522  type(space_t), intent(inout) :: xh
523  type(coef_t), intent(inout) :: coef
524  real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
525  intent(inout) :: cr, cs, ct
526  real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
527  intent(in) :: cx, cy, cz
528 
529  if (neko_bcknd_sx .eq. 1) then
530  call opr_sx_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
531  else if (neko_bcknd_xsmm .eq. 1) then
532  call opr_xsmm_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
533  else
534  call opr_cpu_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
535  end if
536 
537  end subroutine set_convect_rst
538 
554  subroutine runge_kutta(phi, c_r1, c_r23, c_r4, Xh_GLL, Xh_GL, coef, &
555  coef_GL, GLL_to_GL, tau, dtau, n, nel, n_GL)
556  type(space_t), intent(inout) :: xh_gll
557  type(space_t), intent(inout) :: xh_gl
558  type(coef_t), intent(inout) :: coef
559  type(coef_t), intent(inout) :: coef_gl
560  type(interpolator_t) :: gll_to_gl
561  real(kind=rp), intent(inout) :: tau, dtau
562  integer, intent(in) :: n, nel, n_gl
563  real(kind=rp), dimension(n), intent(inout) :: phi
564  real(kind=rp), dimension(3 * n_GL), intent(inout) :: c_r1, c_r23, c_r4
565  real(kind=rp) :: c1, c2, c3
566  ! Work Arrays
567  real(kind=rp), dimension(n) :: u1, r1, r2, r3, r4
568  real(kind=rp), dimension(n_GL) :: u1_gl
569  integer :: i, e
570 
571  c1 = 1.
572  c2 = -dtau/2.
573  c3 = -dtau
574 
575  ! Stage 1:
576  call invcol3 (u1, phi, coef%B, n)
577  call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
578  call convect_scalar(r1, u1_gl, c_r1, xh_gll, xh_gl, coef, &
579  coef_gl, gll_to_gl)
580  call col2(r1, coef%B, n)
581 
582  ! Stage 2:
583  call add3s2 (u1, phi, r1, c1, c2, n)
584  call invcol2 (u1, coef%B, n)
585  call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
586  call convect_scalar(r2, u1_gl, c_r23, xh_gll, xh_gl, coef, &
587  coef_gl, gll_to_gl)
588  call col2(r2, coef%B, n)
589 
590  ! Stage 3:
591  call add3s2 (u1, phi, r2, 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(r3, u1_gl, c_r23, xh_gll, xh_gl, coef, &
595  coef_gl, gll_to_gl)
596  call col2(r3, coef%B, n)
597 
598  ! Stage 4:
599  call add3s2 (u1, phi, r3, c1, c3, n)
600  call invcol2 (u1, coef%B, n)
601  call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
602  call convect_scalar(r4, u1_gl, c_r4, xh_gll, xh_gl, coef, &
603  coef_gl, gll_to_gl)
604  call col2(r4, coef%B, n)
605 
606  c1 = -dtau/6.
607  c2 = -dtau/3.
608  do i = 1, n
609  phi(i) = phi(i) + c1 * (r1(i) + r4(i)) + c2 * (r2(i) + r3(i))
610  end do
611 
612  end subroutine runge_kutta
613 
614 end module operators
615 
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_cmult(a_d, c, n)
Multiplication by constant c .
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:277
subroutine, public invcol2(a, b, n)
Vector division .
Definition: math.f90:675
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition: math.f90:289
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition: math.f90:326
subroutine, public add2(a, b, n)
Vector addition .
Definition: math.f90:547
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition: math.f90:447
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition: math.f90:730
subroutine, public col2(a, b, n)
Vector multiplication .
Definition: math.f90:689
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:228
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:184
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_sx
Definition: neko_config.f90:39
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
logical, parameter neko_device_mpi
Definition: neko_config.f90:46
integer, parameter neko_bcknd_xsmm
Definition: neko_config.f90:40
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.
Definition: operators.f90:560
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.
Definition: operators.f90:169
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.
Definition: operators.f90:526
subroutine, public div(res, ux, uy, uz, coef)
Compute the divergence of a vector field.
Definition: operators.f90:101
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.
Definition: operators.f90:325
subroutine, public grad(ux, uy, uz, u, coef)
Compute the gradient of a scalar field.
Definition: operators.f90:144
subroutine, public curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
Definition: operators.f90:360
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.
Definition: operators.f90:428
subroutine, public lambda2op(lambda2, u, v, w, coef)
Compute the Lambda2 field for a given velocity field.
Definition: operators.f90:501
subroutine, public ortho(x, n, glb_n)
Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T.
Definition: operators.f90:206
subroutine, public cdtp(dtx, x, dr, ds, dt, coef, es, ee)
Apply D^T to a scalar field, where D is the derivative matrix.
Definition: operators.f90:228
subroutine, public dudxyz(du, u, dr, ds, dt, coef)
Compute derivative of a scalar field along a single direction.
Definition: operators.f90:76
real(kind=rp) function, public cfl(dt, u, v, w, Xh, coef, nelv, gdim)
Definition: operators.f90:392
subroutine, public conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
Compute the advection term.
Definition: operators.f90:272
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
subroutine, public opr_cpu_lambda2(lambda2, u, v, w, coef)
Definition: opr_cpu.f90:235
real(kind=rp) function, public opr_cpu_cfl(dt, u, v, w, Xh, coef, nelv, gdim)
Definition: opr_cpu.f90:175
Operators accelerator backends.
Definition: opr_device.F90:34
real(kind=rp) function, public opr_device_cfl(dt, u, v, w, Xh, coef, nelv, gdim)
Definition: opr_device.F90:665
subroutine, public opr_device_cdtp(dtx, x, dr, ds, dt, coef)
Definition: opr_device.F90:432
subroutine, public opr_device_conv1(du, u, vx, vy, vz, Xh, coef, nelv, gdim)
Definition: opr_device.F90:468
subroutine, public opr_device_dudxyz(du, u, dr, ds, dt, coef)
Definition: opr_device.F90:325
subroutine, public opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_Xh)
Definition: opr_device.F90:515
subroutine, public opr_device_opgrad(ux, uy, uz, u, coef)
Definition: opr_device.F90:360
subroutine, public opr_device_lambda2(lambda2, u, v, w, coef)
Definition: opr_device.F90:402
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_convect_scalar(du, u, c, Xh_GLL, Xh_GL, coef_GLL, coef_GL, GLL_to_GL)
Definition: opr_xsmm.F90:371
subroutine, public opr_xsmm_cdtp(dtx, x, dr, ds, dt, coef)
Definition: opr_xsmm.F90:242
subroutine, public opr_xsmm_set_convect_rst(cr, cs, ct, cx, cy, cz, Xh, coef)
Definition: opr_xsmm.F90:467
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_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_Xh)
Definition: opr_xsmm.F90:416
subroutine, public opr_xsmm_conv1(du, u, vx, vy, vz, Xh, coef, nelv, gdim)
Definition: opr_xsmm.F90:297
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