45 use,
intrinsic :: iso_c_binding
49 integer,
public,
parameter ::
gl = 0,
gll = 1,
gj = 2
73 real(kind=
rp),
allocatable :: zg(:,:)
75 real(kind=
rp),
allocatable :: dr_inv(:)
76 real(kind=
rp),
allocatable :: ds_inv(:)
77 real(kind=
rp),
allocatable :: dt_inv(:)
79 real(kind=
rp),
allocatable :: wx(:)
80 real(kind=
rp),
allocatable :: wy(:)
81 real(kind=
rp),
allocatable :: wz(:)
83 real(kind=
rp),
allocatable :: w3(:,:,:)
86 real(kind=
rp),
allocatable :: dx(:,:)
88 real(kind=
rp),
allocatable :: dy(:,:)
90 real(kind=
rp),
allocatable :: dz(:,:)
93 real(kind=
rp),
allocatable :: dxt(:,:)
95 real(kind=
rp),
allocatable :: dyt(:,:)
97 real(kind=
rp),
allocatable :: dzt(:,:)
100 real(kind=
rp),
allocatable :: v(:,:)
101 real(kind=
rp),
allocatable :: vt(:,:)
102 real(kind=
rp),
allocatable :: vinv(:,:)
103 real(kind=
rp),
allocatable :: vinvt(:,:)
105 real(kind=
rp),
allocatable :: w(:,:)
110 type(c_ptr) :: dr_inv_d = c_null_ptr
111 type(c_ptr) :: ds_inv_d = c_null_ptr
112 type(c_ptr) :: dt_inv_d = c_null_ptr
113 type(c_ptr) :: dxt_d = c_null_ptr
114 type(c_ptr) :: dyt_d = c_null_ptr
115 type(c_ptr) :: dzt_d = c_null_ptr
116 type(c_ptr) :: dx_d = c_null_ptr
117 type(c_ptr) :: dy_d = c_null_ptr
118 type(c_ptr) :: dz_d = c_null_ptr
119 type(c_ptr) :: wx_d = c_null_ptr
120 type(c_ptr) :: wy_d = c_null_ptr
121 type(c_ptr) :: wz_d = c_null_ptr
122 type(c_ptr) :: zg_d = c_null_ptr
123 type(c_ptr) :: w3_d = c_null_ptr
124 type(c_ptr) :: v_d = c_null_ptr
125 type(c_ptr) :: vt_d = c_null_ptr
126 type(c_ptr) :: vinv_d = c_null_ptr
127 type(c_ptr) :: vinvt_d = c_null_ptr
128 type(c_ptr) :: w_d = c_null_ptr
135 interface operator(.eq.)
137 end interface operator(.eq.)
139 interface operator(.ne.)
141 end interface operator(.ne.)
143 public ::
operator(.eq.),
operator(.ne.)
149 class(
space_t),
intent(inout) :: s
150 integer,
intent(in) :: t
151 integer,
intent(in) :: lx
152 integer,
intent(in) :: ly
153 integer,
optional,
intent(in) :: lz
154 integer :: ix, iy, iz
161 if (
present(lz))
then
164 if (lx .ne. ly .or. lx .ne. lz)
then
165 call neko_error(
"Unsupported polynomial dimension")
170 call neko_error(
"Unsupported polynomial dimension")
177 s%lxyz = s%lx*s%ly*s%lz
179 allocate(s%zg(lx, 3))
185 allocate(s%dr_inv(s%lx))
186 allocate(s%ds_inv(s%ly))
187 allocate(s%dt_inv(s%lz))
189 allocate(s%w3(s%lx, s%ly, s%lz))
191 allocate(s%dx(s%lx, s%lx))
192 allocate(s%dy(s%ly, s%ly))
193 allocate(s%dz(s%lz, s%lz))
195 allocate(s%dxt(s%lx, s%lx))
196 allocate(s%dyt(s%ly, s%ly))
197 allocate(s%dzt(s%lz, s%lz))
199 allocate(s%v(s%lx, s%lx))
200 allocate(s%vt(s%lx, s%lx))
201 allocate(s%vinv(s%lx, s%lx))
202 allocate(s%vinvt(s%lx, s%lx))
203 allocate(s%w(s%lx, s%lx))
207 call zwgll(s%zg(1,1), s%wx, s%lx)
208 call zwgll(s%zg(1,2), s%wy, s%ly)
209 if (s%lz .gt. 1)
then
210 call zwgll(s%zg(1,3), s%wz, s%lz)
215 else if (t .eq.
gl)
then
216 call zwgl(s%zg(1,1), s%wx, s%lx)
217 call zwgl(s%zg(1,2), s%wy, s%ly)
218 if (s%lz .gt. 1)
then
219 call zwgl(s%zg(1,3), s%wz, s%lz)
231 s%w3(ix, iy, iz) = s%wx(ix) * s%wy(iy) * s%wz(iz)
237 call dgll(s%dx, s%dxt, s%zg(1,1), s%lx, s%lx)
238 call dgll(s%dy, s%dyt, s%zg(1,2), s%ly, s%ly)
239 if (s%lz .gt. 1)
then
240 call dgll(s%dz, s%dzt, s%zg(1,3), s%lz, s%lz)
245 else if (t .eq.
gl)
then
246 call setup_intp(s%dx, s%dxt, s%zg(1,1), s%zg(1,1), s%lx, s%lx,1)
247 call setup_intp(s%dy, s%dyt, s%zg(1,2), s%zg(1,2), s%ly, s%ly,1)
248 if (s%lz .gt. 1)
then
249 call setup_intp(s%dz, s%dzt, s%zg(1,3), s%zg(1,3), s%lz, s%lz, 1)
260 if (s%lz .gt. 1)
then
323 class(
space_t),
intent(inout) :: s
325 if (
allocated(s%zg))
then
330 if (
allocated(s%wx))
then
335 if (
allocated(s%wy))
then
340 if (
allocated(s%wz))
then
345 if (
allocated(s%w3))
then
350 if (
allocated(s%dx))
then
355 if (
allocated(s%dy))
then
360 if (
allocated(s%dz))
then
365 if (
allocated(s%dxt))
then
370 if (
allocated(s%dyt))
then
375 if (
allocated(s%dzt))
then
380 if (
allocated(s%dr_inv))
then
385 if (
allocated(s%ds_inv))
then
390 if (
allocated(s%dt_inv))
then
395 if (
allocated(s%v))
then
400 if (
allocated(s%vt))
then
405 if (
allocated(s%vinv))
then
410 if (
allocated(s%vinvt))
then
415 if (
allocated(s%w))
then
425 type(
space_t),
intent(in) :: xh
426 type(
space_t),
intent(in) :: yh
429 if ( (xh%lx .eq. yh%lx) .and. &
430 (xh%ly .eq. yh%ly) .and. &
431 (xh%lz .eq. yh%lz) )
then
442 type(
space_t),
intent(in) :: xh
443 type(
space_t),
intent(in) :: yh
446 if ( (xh%lx .eq. yh%lx) .and. &
447 (xh%ly .eq. yh%ly) .and. &
448 (xh%lz .eq. yh%lz) )
then
457 integer,
intent(in) :: lx
458 real(kind=
rp),
intent(inout) :: dx(lx), x(lx)
462 dx(i) = 0.5*(x(i+1) - x(i-1))
464 dx(lx) = x(lx) - x(lx-1)
466 dx(i) = 1.0_rp / dx(i)
474 type(
space_t),
intent(inout) :: Xh
476 real(kind=
rp) :: l(0:xh%lx-1)
477 real(kind=
rp) :: delta(xh%lx)
478 integer :: i, kj, j, kk
480 logical :: scaled = .false.
482 associate(v=> xh%v, vt => xh%vt, &
483 vinv => xh%vinv, vinvt => xh%vinvt, w => xh%w)
502 delta(i) = 2.0_rp / (2*(i-1)+1)
505 delta(xh%lx) = 2.0_rp / (xh%lx-1)
509 delta(i) = sqrt(1.0_rp / delta(i))
514 v(i,j) = v(i,j) * delta(j)
519 call copy(vt, v, xh%lx * xh%lx)
520 call trsp1(vt, xh%lx)
536 call mxm(vt, xh%lx, w, xh%lx, vinv, xh%lx)
539 call copy(vinvt, vinv, xh%lx * xh%lx)
540 call trsp1(vinvt, xh%lx)
542 call copy(vt, v, xh%lxy)
543 call trsp1(vt, xh%lx)
544 call m%init(xh%lx, xh%lx)
545 call copy(m%x, v, xh%lxy)
546 call m%inverse_on_host()
547 call copy(vinv, m%x, xh%lxy)
548 call copy(vinvt, vinv, xh%lx * xh%lx)
549 call trsp1(vinvt, xh%lx)
Map a Fortran array to a device (allocate and associate)
Copy data between host and device (or device and device)
Synchronize a device or stream.
Unmap a Fortran array from a device (deassociate and free)
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
Fast diagonalization methods from NEKTON.
subroutine, public setup_intp(jh, jht, z_to, z_from, n_to, n_from, derivative)
Compute interpolation weights for points z_to using values at points z_from.
subroutine, public copy(a, b, n)
Copy a vector .
Wrapper for all matrix-matrix product implementations.
subroutine, public mxm(a, n1, b, n2, c, n3)
Compute matrix-matrix product for contiguously packed matrices A,B, and C.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Defines a function space.
pure logical function space_ne(xh, yh)
Check if .
pure logical function space_eq(xh, yh)
Check if .
integer, parameter, public gll
subroutine space_compute_dist(dx, x, lx)
integer, parameter, public gj
subroutine space_free(s)
Deallocate a space s.
subroutine space_init(s, t, lx, ly, lz)
Initialize a function space s with given polynomial dimensions.
integer, parameter, public gl
subroutine space_generate_transformation_matrices(xh)
Generate spectral tranform matrices.
LIBRARY ROUTINES FOR SPECTRAL METHODS.
subroutine dgll(d, dt, z, nz, nzd)
Compute the derivative matrix D and its transpose DT associated with the Nth order Lagrangian interpo...
subroutine zwgll(z, w, np)
Generate NP Gauss-Lobatto Legendre points (Z) and weights (W) associated with Jacobi polynomial P(N)(...
subroutine legendre_poly(l, x, n)
Evaluate Legendre polynomials of degrees 0-N at point x and store in array L.
subroutine zwgl(z, w, np)
Generate NP Gauss Legendre points Z and weights W associated with Jacobi polynomial ....
subroutine, public trsp1(a, n)
In-place transpose of a square tensor.
The function space for the SEM solution fields.