Neko  0.9.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  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, n, glb_n)
204  integer, intent(in) :: n
205  integer, intent(in) :: glb_n
206  real(kind=rp), dimension(n), intent(inout) :: x
207  real(kind=rp) :: rlam
208 
209  rlam = glsum(x, n)/glb_n
210  call cadd(x, -rlam, n)
211 
212  end subroutine ortho
213 
225  subroutine cdtp (dtx, x, dr, ds, dt, coef, es, ee)
226  type(coef_t), intent(in) :: coef
227  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: dtx
228  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(inout) :: x
229  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: dr
230  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: ds
231  real(kind=rp), dimension(coef%Xh%lxyz, coef%msh%nelv), intent(in) :: dt
232  integer, optional :: es, ee
233  integer :: eblk_start, eblk_end
234 
235  if (present(es)) then
236  eblk_start = es
237  else
238  eblk_start = 1
239  end if
240 
241  if (present(ee)) then
242  eblk_end = ee
243  else
244  eblk_end = coef%msh%nelv
245  end if
246 
247  if (neko_bcknd_sx .eq. 1) then
248  call opr_sx_cdtp(dtx, x, dr, ds, dt, coef)
249  else if (neko_bcknd_xsmm .eq. 1) then
250  call opr_xsmm_cdtp(dtx, x, dr, ds, dt, coef)
251  else if (neko_bcknd_device .eq. 1) then
252  call opr_device_cdtp(dtx, x, dr, ds, dt, coef)
253  else
254  call opr_cpu_cdtp(dtx, x, dr, ds, dt, coef, eblk_start, eblk_end)
255  end if
256 
257  end subroutine cdtp
258 
269  subroutine conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
270  type(space_t), intent(inout) :: xh
271  type(coef_t), intent(inout) :: coef
272  real(kind=rp), intent(inout) :: du(xh%lxyz, coef%msh%nelv)
273  real(kind=rp), intent(inout) :: u(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
274  real(kind=rp), intent(inout) :: vx(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
275  real(kind=rp), intent(inout) :: vy(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
276  real(kind=rp), intent(inout) :: vz(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
277  integer, optional :: es, ee
278  integer :: eblk_end, eblk_start
279 
280  associate(nelv => coef%msh%nelv, gdim => coef%msh%gdim)
281  if (present(es)) then
282  eblk_start = es
283  else
284  eblk_start = 1
285  end if
286 
287  if (present(ee)) then
288  eblk_end = ee
289  else
290  eblk_end = coef%msh%nelv
291  end if
292 
293  if (neko_bcknd_sx .eq. 1) then
294  call opr_sx_conv1(du, u, vx, vy, vz, xh, coef, nelv)
295  else if (neko_bcknd_xsmm .eq. 1) then
296  call opr_xsmm_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
297  else if (neko_bcknd_device .eq. 1) then
298  call opr_device_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
299  else
300  call opr_cpu_conv1(du, u, vx, vy, vz, xh, coef, eblk_start, eblk_end)
301  end if
302  end associate
303 
304  end subroutine conv1
305 
321  subroutine convect_scalar(du, u, c, Xh_GLL, Xh_GL, coef_GLL, &
322  coef_GL, GLL_to_GL)
323  type(space_t), intent(in) :: xh_gl
324  type(space_t), intent(in) :: xh_gll
325  type(coef_t), intent(in) :: coef_GLL
326  type(coef_t), intent(in) :: coef_GL
327  type(interpolator_t), intent(inout) :: GLL_to_GL
328  real(kind=rp), intent(inout) :: &
329  du(xh_gll%lx, xh_gll%ly, xh_gll%lz, coef_gl%msh%nelv)
330  real(kind=rp), intent(inout) :: &
331  u(xh_gl%lx, xh_gl%lx, xh_gl%lx, coef_gl%msh%nelv)
332  real(kind=rp), intent(inout) :: c(xh_gl%lxyz, coef_gl%msh%nelv, 3)
333 
334  if (neko_bcknd_sx .eq. 1) then
335  call opr_sx_convect_scalar(du, u, c, xh_gll, xh_gl, &
336  coef_gll, coef_gl, gll_to_gl)
337  else if (neko_bcknd_xsmm .eq. 1) then
338  call opr_xsmm_convect_scalar(du, u, c, xh_gll, xh_gl, &
339  coef_gll, coef_gl, gll_to_gl)
340  else
341  call opr_cpu_convect_scalar(du, u, c, xh_gll, xh_gl, &
342  coef_gll, coef_gl, gll_to_gl)
343  end if
344 
345  end subroutine convect_scalar
346 
347  !! Compute the curl fo a vector field.
348  !! @param w1 Will store the x component of the curl.
349  !! @param w2 Will store the y component of the curl.
350  !! @param w3 Will store the z component of the curl.
351  !! @param u1 The x component of the vector field.
352  !! @param u2 The y component of the vector field.
353  !! @param u3 The z component of the vector field.
354  !! @param work1 A temporary array for computations.
355  !! @param work2 A temporary array for computations.
356  !! @param coef The SEM coefficients.
357  subroutine curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
358  type(field_t), intent(inout) :: w1
359  type(field_t), intent(inout) :: w2
360  type(field_t), intent(inout) :: w3
361  type(field_t), intent(inout) :: u1
362  type(field_t), intent(inout) :: u2
363  type(field_t), intent(inout) :: u3
364  type(field_t), intent(inout) :: work1
365  type(field_t), intent(inout) :: work2
366  type(coef_t), intent(in) :: coef
367 
368  if (neko_bcknd_sx .eq. 1) then
369  call opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
370  else if (neko_bcknd_xsmm .eq. 1) then
371  call opr_xsmm_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
372  else if (neko_bcknd_device .eq. 1) then
373  call opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
374  else
375  call opr_cpu_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
376  end if
377 
378  end subroutine curl
379 
380  !! Compute the CFL number
381  !! @param dt The timestep.
382  !! @param u The x component of velocity.
383  !! @param v The y component of velocity.
384  !! @param w The z component of velocity.
385  !! @param Xh The SEM function space.
386  !! @param coef The SEM coefficients.
387  !! @param nelv The total number of elements.
388  !! @param gdim Number of geometric dimensions.
389  function cfl(dt, u, v, w, Xh, coef, nelv, gdim)
390  type(space_t), intent(in) :: xh
391  type(coef_t), intent(in) :: coef
392  integer, intent(in) :: nelv, gdim
393  real(kind=rp), intent(in) :: dt
394  real(kind=rp), dimension(Xh%lx, Xh%ly, Xh%lz, nelv), intent(in) :: u, v, w
395  real(kind=rp) :: cfl
396  integer :: ierr
397 
398  if (neko_bcknd_sx .eq. 1) then
399  cfl = opr_sx_cfl(dt, u, v, w, xh, coef, nelv)
400  else if (neko_bcknd_device .eq. 1) then
401  cfl = opr_device_cfl(dt, u, v, w, xh, coef, nelv, gdim)
402  else
403  cfl = opr_cpu_cfl(dt, u, v, w, xh, coef, nelv, gdim)
404  end if
405 
406  if (.not. neko_device_mpi) then
407  call mpi_allreduce(mpi_in_place, cfl, 1, &
408  mpi_real_precision, mpi_max, neko_comm, ierr)
409  end if
410 
411  end function cfl
412 
425  subroutine strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
426  type(field_t), intent(in) :: u, v, w
427  type(coef_t), intent(in) :: coef
428  real(kind=rp), intent(inout) :: s11(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
429  real(kind=rp), intent(inout) :: s22(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
430  real(kind=rp), intent(inout) :: s33(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
431  real(kind=rp), intent(inout) :: s12(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
432  real(kind=rp), intent(inout) :: s13(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
433  real(kind=rp), intent(inout) :: s23(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
434 
435  type(c_ptr) :: s11_d, s22_d, s33_d, s12_d, s23_d, s13_d
436 
437  integer :: nelv, lxyz
438 
439  if (neko_bcknd_device .eq. 1) then
440  s11_d = device_get_ptr(s11)
441  s22_d = device_get_ptr(s22)
442  s33_d = device_get_ptr(s33)
443  s12_d = device_get_ptr(s12)
444  s23_d = device_get_ptr(s23)
445  s13_d = device_get_ptr(s13)
446  end if
447 
448  nelv = u%msh%nelv
449  lxyz = u%Xh%lxyz
450 
451  ! we use s11 as a work array here
452  call dudxyz (s12, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
453  call dudxyz (s11, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
454  if (neko_bcknd_device .eq. 1) then
455  call device_add2(s12_d, s11_d, nelv*lxyz)
456  else
457  call add2(s12, s11, nelv*lxyz)
458  end if
459 
460  call dudxyz (s13, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
461  call dudxyz (s11, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
462  if (neko_bcknd_device .eq. 1) then
463  call device_add2(s13_d, s11_d, nelv*lxyz)
464  else
465  call add2(s13, s11, nelv*lxyz)
466  end if
467 
468  call dudxyz (s23, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
469  call dudxyz (s11, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
470  if (neko_bcknd_device .eq. 1) then
471  call device_add2(s23_d, s11_d, nelv*lxyz)
472  else
473  call add2(s23, s11, nelv*lxyz)
474  end if
475 
476  call dudxyz (s11, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
477  call dudxyz (s22, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
478  call dudxyz (s33, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
479 
480  if (neko_bcknd_device .eq. 1) then
481  call device_cmult(s12_d, 0.5_rp, nelv*lxyz)
482  call device_cmult(s13_d, 0.5_rp, nelv*lxyz)
483  call device_cmult(s23_d, 0.5_rp, nelv*lxyz)
484  else
485  call cmult(s12, 0.5_rp, nelv*lxyz)
486  call cmult(s13, 0.5_rp, nelv*lxyz)
487  call cmult(s23, 0.5_rp, nelv*lxyz)
488  end if
489 
490  end subroutine strain_rate
491 
498  subroutine lambda2op(lambda2, u, v, w, coef)
499  type(coef_t), intent(in) :: coef
500  type(field_t), intent(inout) :: lambda2
501  type(field_t), intent(in) :: u, v, w
502 
503  if (neko_bcknd_sx .eq. 1) then
504  call opr_sx_lambda2(lambda2, u, v, w, coef)
505  else if (neko_bcknd_device .eq. 1) then
506  call opr_device_lambda2(lambda2, u, v, w, coef)
507  else
508  call opr_cpu_lambda2(lambda2, u, v, w, coef)
509  end if
510 
511  end subroutine lambda2op
512 
523  subroutine set_convect_rst(cr, cs, ct, cx, cy, cz, Xh, coef)
524  type(space_t), intent(inout) :: xh
525  type(coef_t), intent(inout) :: coef
526  real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
527  intent(inout) :: cr, cs, ct
528  real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
529  intent(in) :: cx, cy, cz
530 
531  if (neko_bcknd_sx .eq. 1) then
532  call opr_sx_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
533  else if (neko_bcknd_xsmm .eq. 1) then
534  call opr_xsmm_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
535  else
536  call opr_cpu_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
537  end if
538 
539  end subroutine set_convect_rst
540 
556  subroutine runge_kutta(phi, c_r1, c_r23, c_r4, Xh_GLL, Xh_GL, coef, &
557  coef_GL, GLL_to_GL, tau, dtau, n, nel, n_GL)
558  type(space_t), intent(inout) :: xh_gll
559  type(space_t), intent(inout) :: xh_gl
560  type(coef_t), intent(inout) :: coef
561  type(coef_t), intent(inout) :: coef_gl
562  type(interpolator_t) :: gll_to_gl
563  real(kind=rp), intent(inout) :: tau, dtau
564  integer, intent(in) :: n, nel, n_gl
565  real(kind=rp), dimension(n), intent(inout) :: phi
566  real(kind=rp), dimension(3 * n_GL), intent(inout) :: c_r1, c_r23, c_r4
567  real(kind=rp) :: c1, c2, c3
568  ! Work Arrays
569  real(kind=rp), dimension(n) :: u1, r1, r2, r3, r4
570  real(kind=rp), dimension(n_GL) :: u1_gl
571  integer :: i, e
572 
573  c1 = 1.
574  c2 = -dtau/2.
575  c3 = -dtau
576 
577  ! Stage 1:
578  call invcol3 (u1, phi, coef%B, n)
579  call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
580  call convect_scalar(r1, u1_gl, c_r1, xh_gll, xh_gl, coef, &
581  coef_gl, gll_to_gl)
582  call col2(r1, coef%B, n)
583 
584  ! Stage 2:
585  call add3s2 (u1, phi, r1, c1, c2, n)
586  call invcol2 (u1, coef%B, n)
587  call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
588  call convect_scalar(r2, u1_gl, c_r23, xh_gll, xh_gl, coef, &
589  coef_gl, gll_to_gl)
590  call col2(r2, coef%B, n)
591 
592  ! Stage 3:
593  call add3s2 (u1, phi, r2, 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(r3, u1_gl, c_r23, xh_gll, xh_gl, coef, &
597  coef_gl, gll_to_gl)
598  call col2(r3, coef%B, n)
599 
600  ! Stage 4:
601  call add3s2 (u1, phi, r3, c1, c3, n)
602  call invcol2 (u1, coef%B, n)
603  call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
604  call convect_scalar(r4, u1_gl, c_r4, xh_gll, xh_gl, coef, &
605  coef_gl, gll_to_gl)
606  call col2(r4, coef%B, n)
607 
608  c1 = -dtau/6.
609  c2 = -dtau/3.
610  do i = 1, n
611  phi(i) = phi(i) + c1 * (r1(i) + r4(i)) + c2 * (r2(i) + r3(i))
612  end do
613 
614  end subroutine runge_kutta
615 
616 end module operators
617 
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 .
Definition: device_math.F90:76
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:311
subroutine, public invcol2(a, b, n)
Vector division .
Definition: math.f90:715
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition: math.f90:323
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition: math.f90:360
subroutine, public add2(a, b, n)
Vector addition .
Definition: math.f90:587
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition: math.f90:487
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition: math.f90:770
subroutine, public col2(a, b, n)
Vector multiplication .
Definition: math.f90:729
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:239
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:195
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:562
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:171
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:528
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:327
subroutine, public grad(ux, uy, uz, u, coef)
Compute the gradient of a scalar field.
Definition: operators.f90:146
subroutine, public curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
Definition: operators.f90:362
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:430
subroutine, public lambda2op(lambda2, u, v, w, coef)
Compute the Lambda2 field for a given velocity field.
Definition: operators.f90:503
subroutine, public ortho(x, n, glb_n)
Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T.
Definition: operators.f90:208
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:230
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:394
subroutine, public conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
Compute the advection term.
Definition: operators.f90:274
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