Neko  0.9.0
A portable framework for high-order spectral element flow simulations
sx_conv1.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021-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 submodule(opr_sx) sx_conv1
35  implicit none
36 
37 contains
38 
39  module subroutine opr_sx_conv1(du, u, vx, vy, vz, xh, coef, nelv)
40  type(space_t), intent(inout) :: Xh
41  type(coef_t), intent(inout) :: coef
42  integer, intent(in) :: nelv
43  real(kind=rp), intent(inout) :: du(xh%lxyz, nelv)
44  real(kind=rp), intent(inout) :: u(xh%lx, xh%ly, xh%lz, nelv)
45  real(kind=rp), intent(inout) :: vx(xh%lx, xh%ly, xh%lz, nelv)
46  real(kind=rp), intent(inout) :: vy(xh%lx, xh%ly, xh%lz, nelv)
47  real(kind=rp), intent(inout) :: vz(xh%lx, xh%ly, xh%lz, nelv)
48 
49  select case (xh%lx)
50  case (14)
51  call sx_conv1_lx14(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
52  coef%drdx, coef%dsdx, coef%dtdx, &
53  coef%drdy, coef%dsdy, coef%dtdy, &
54  coef%drdz, coef%dsdz, coef%dtdz, &
55  coef%jacinv, nelv)
56  case (13)
57  call sx_conv1_lx13(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
58  coef%drdx, coef%dsdx, coef%dtdx, &
59  coef%drdy, coef%dsdy, coef%dtdy, &
60  coef%drdz, coef%dsdz, coef%dtdz, &
61  coef%jacinv, nelv)
62  case (12)
63  call sx_conv1_lx12(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
64  coef%drdx, coef%dsdx, coef%dtdx, &
65  coef%drdy, coef%dsdy, coef%dtdy, &
66  coef%drdz, coef%dsdz, coef%dtdz, &
67  coef%jacinv, nelv)
68  case (11)
69  call sx_conv1_lx11(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
70  coef%drdx, coef%dsdx, coef%dtdx, &
71  coef%drdy, coef%dsdy, coef%dtdy, &
72  coef%drdz, coef%dsdz, coef%dtdz, &
73  coef%jacinv, nelv)
74  case (10)
75  call sx_conv1_lx10(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
76  coef%drdx, coef%dsdx, coef%dtdx, &
77  coef%drdy, coef%dsdy, coef%dtdy, &
78  coef%drdz, coef%dsdz, coef%dtdz, &
79  coef%jacinv, nelv)
80  case (9)
81  call sx_conv1_lx9(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
82  coef%drdx, coef%dsdx, coef%dtdx, &
83  coef%drdy, coef%dsdy, coef%dtdy, &
84  coef%drdz, coef%dsdz, coef%dtdz, &
85  coef%jacinv, nelv)
86  case (8)
87  call sx_conv1_lx8(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
88  coef%drdx, coef%dsdx, coef%dtdx, &
89  coef%drdy, coef%dsdy, coef%dtdy, &
90  coef%drdz, coef%dsdz, coef%dtdz, &
91  coef%jacinv, nelv)
92  case (7)
93  call sx_conv1_lx7(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
94  coef%drdx, coef%dsdx, coef%dtdx, &
95  coef%drdy, coef%dsdy, coef%dtdy, &
96  coef%drdz, coef%dsdz, coef%dtdz, &
97  coef%jacinv, nelv)
98  case (6)
99  call sx_conv1_lx6(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
100  coef%drdx, coef%dsdx, coef%dtdx, &
101  coef%drdy, coef%dsdy, coef%dtdy, &
102  coef%drdz, coef%dsdz, coef%dtdz, &
103  coef%jacinv, nelv)
104  case (5)
105  call sx_conv1_lx5(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
106  coef%drdx, coef%dsdx, coef%dtdx, &
107  coef%drdy, coef%dsdy, coef%dtdy, &
108  coef%drdz, coef%dsdz, coef%dtdz, &
109  coef%jacinv, nelv)
110  case (4)
111  call sx_conv1_lx4(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
112  coef%drdx, coef%dsdx, coef%dtdx, &
113  coef%drdy, coef%dsdy, coef%dtdy, &
114  coef%drdz, coef%dsdz, coef%dtdz, &
115  coef%jacinv, nelv)
116  case (3)
117  call sx_conv1_lx3(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
118  coef%drdx, coef%dsdx, coef%dtdx, &
119  coef%drdy, coef%dsdy, coef%dtdy, &
120  coef%drdz, coef%dsdz, coef%dtdz, &
121  coef%jacinv, nelv)
122  case (2)
123  call sx_conv1_lx2(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
124  coef%drdx, coef%dsdx, coef%dtdx, &
125  coef%drdy, coef%dsdy, coef%dtdy, &
126  coef%drdz, coef%dsdz, coef%dtdz, &
127  coef%jacinv, nelv)
128  case default
129  call sx_conv1_lx(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
130  coef%drdx, coef%dsdx, coef%dtdx, &
131  coef%drdy, coef%dsdy, coef%dtdy, &
132  coef%drdz, coef%dsdz, coef%dtdz, &
133  coef%jacinv, nelv, xh%lx)
134  end select
135 
136  end subroutine opr_sx_conv1
137 
138  subroutine sx_conv1_lx(du, u, vx, vy, vz, dx, dy, dz, &
139  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
140  jacinv, nelv, lx)
141  integer, intent(in) :: nelv, lx
142  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
143  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
144  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
145  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
146  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
147  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
148  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
149  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
150  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
151  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
152  real(kind=rp) :: wr, ws, wt
153  integer :: e, i, j, k, jj, kk
154 
155  do i = 1, lx
156  do jj = 1, lx * lx * nelv
157  wr = 0d0
158  do kk = 1, lx
159  wr = wr + dx(i, kk) * u(kk, jj,1,1)
160  end do
161  dudr(i, jj,1,1) = wr
162  end do
163  end do
164 
165  do k = 1, lx
166  do i = 1, lx
167  do j = 1, lx
168  do e = 1, nelv
169  ws = 0d0
170  !NEC$ unroll_completely
171  do kk = 1, lx
172  ws = ws + dy(j, kk) * u(i, kk,k,e)
173  end do
174  duds(i,j,k,e) = ws
175  end do
176  end do
177  end do
178  end do
179 
180  do j = 1, lx
181  do i = 1, lx
182  do k = 1, lx
183  do e = 1, nelv
184  wt = 0d0
185  !NEC$ unroll_completely
186  do kk = 1, lx
187  wt = wt + dz(k, kk) * u(i,j, kk,e)
188  end do
189  dudt(i,j,k,e) = wt
190  end do
191  end do
192  end do
193  end do
194 
195  do i = 1, nelv * lx * lx * lx
196  du(i,1,1,1) = jacinv(i,1,1,1) * &
197  (vx(i,1,1,1) * ( &
198  drdx(i,1,1,1)*dudr(i,1,1,1) &
199  + dsdx(i,1,1,1)*duds(i,1,1,1) &
200  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
201  + vy(i,1,1,1)*( &
202  drdy(i,1,1,1)*dudr(i,1,1,1) &
203  + dsdy(i,1,1,1)*duds(i,1,1,1) &
204  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
205  + vz(i,1,1,1)*( &
206  drdz(i,1,1,1)*dudr(i,1,1,1) &
207  + dsdz(i,1,1,1)*duds(i,1,1,1) &
208  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
209  end do
210 
211  end subroutine sx_conv1_lx
212 
213  subroutine sx_conv1_lx14(du, u, vx, vy, vz, dx, dy, dz, &
214  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
215  jacinv, nelv)
216  integer, parameter :: lx = 14
217  integer, intent(in) :: nelv
218  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
219  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
220  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
221  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
222  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
223  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
224  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
225  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
226  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
227  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
228  real(kind=rp) :: wr, ws, wt
229  integer :: e, i, j, k, jj, kk
230 
231  do i = 1, lx
232  do jj = 1, lx * lx * nelv
233  wr = 0d0
234  do kk = 1, lx
235  wr = wr + dx(i, kk) * u(kk, jj,1,1)
236  end do
237  dudr(i, jj,1,1) = wr
238  end do
239  end do
240 
241  do k = 1, lx
242  do i = 1, lx
243  do j = 1, lx
244  do e = 1, nelv
245  ws = 0d0
246  !NEC$ unroll_completely
247  do kk = 1, lx
248  ws = ws + dy(j, kk) * u(i, kk,k,e)
249  end do
250  duds(i,j,k,e) = ws
251  end do
252  end do
253  end do
254  end do
255 
256  do j = 1, lx
257  do i = 1, lx
258  do k = 1, lx
259  do e = 1, nelv
260  wt = 0d0
261  !NEC$ unroll_completely
262  do kk = 1, lx
263  wt = wt + dz(k, kk) * u(i,j, kk,e)
264  end do
265  dudt(i,j,k,e) = wt
266  end do
267  end do
268  end do
269  end do
270 
271  do i = 1, nelv * lx * lx * lx
272  du(i,1,1,1) = jacinv(i,1,1,1) * &
273  (vx(i,1,1,1) * ( &
274  drdx(i,1,1,1)*dudr(i,1,1,1) &
275  + dsdx(i,1,1,1)*duds(i,1,1,1) &
276  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
277  + vy(i,1,1,1)*( &
278  drdy(i,1,1,1)*dudr(i,1,1,1) &
279  + dsdy(i,1,1,1)*duds(i,1,1,1) &
280  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
281  + vz(i,1,1,1)*( &
282  drdz(i,1,1,1)*dudr(i,1,1,1) &
283  + dsdz(i,1,1,1)*duds(i,1,1,1) &
284  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
285  end do
286 
287  end subroutine sx_conv1_lx14
288 
289  subroutine sx_conv1_lx13(du, u, vx, vy, vz, dx, dy, dz, &
290  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
291  jacinv, nelv)
292  integer, parameter :: lx = 13
293  integer, intent(in) :: nelv
294  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
295  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
296  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
297  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
298  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
299  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
300  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
301  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
302  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
303  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
304  real(kind=rp) :: wr, ws, wt
305  integer :: e, i, j, k, jj, kk
306 
307  do i = 1, lx
308  do jj = 1, lx * lx * nelv
309  wr = 0d0
310  do kk = 1, lx
311  wr = wr + dx(i, kk) * u(kk, jj,1,1)
312  end do
313  dudr(i, jj,1,1) = wr
314  end do
315  end do
316 
317  do k = 1, lx
318  do i = 1, lx
319  do j = 1, lx
320  do e = 1, nelv
321  ws = 0d0
322  !NEC$ unroll_completely
323  do kk = 1, lx
324  ws = ws + dy(j, kk) * u(i, kk,k,e)
325  end do
326  duds(i,j,k,e) = ws
327  end do
328  end do
329  end do
330  end do
331 
332  do j = 1, lx
333  do i = 1, lx
334  do k = 1, lx
335  do e = 1, nelv
336  wt = 0d0
337  !NEC$ unroll_completely
338  do kk = 1, lx
339  wt = wt + dz(k, kk) * u(i,j, kk,e)
340  end do
341  dudt(i,j,k,e) = wt
342  end do
343  end do
344  end do
345  end do
346 
347  do i = 1, nelv * lx * lx * lx
348  du(i,1,1,1) = jacinv(i,1,1,1) * &
349  (vx(i,1,1,1) * ( &
350  drdx(i,1,1,1)*dudr(i,1,1,1) &
351  + dsdx(i,1,1,1)*duds(i,1,1,1) &
352  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
353  + vy(i,1,1,1)*( &
354  drdy(i,1,1,1)*dudr(i,1,1,1) &
355  + dsdy(i,1,1,1)*duds(i,1,1,1) &
356  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
357  + vz(i,1,1,1)*( &
358  drdz(i,1,1,1)*dudr(i,1,1,1) &
359  + dsdz(i,1,1,1)*duds(i,1,1,1) &
360  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
361  end do
362 
363  end subroutine sx_conv1_lx13
364 
365  subroutine sx_conv1_lx12(du, u, vx, vy, vz, dx, dy, dz, &
366  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
367  jacinv, nelv)
368  integer, parameter :: lx = 12
369  integer, intent(in) :: nelv
370  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
371  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
372  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
373  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
374  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
375  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
376  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
377  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
378  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
379  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
380  real(kind=rp) :: wr, ws, wt
381  integer :: e, i, j, k, jj, kk
382 
383  do i = 1, lx
384  do jj = 1, lx * lx * nelv
385  wr = 0d0
386  do kk = 1, lx
387  wr = wr + dx(i, kk) * u(kk, jj,1,1)
388  end do
389  dudr(i, jj,1,1) = wr
390  end do
391  end do
392 
393  do k = 1, lx
394  do i = 1, lx
395  do j = 1, lx
396  do e = 1, nelv
397  ws = 0d0
398  !NEC$ unroll_completely
399  do kk = 1, lx
400  ws = ws + dy(j, kk) * u(i, kk,k,e)
401  end do
402  duds(i,j,k,e) = ws
403  end do
404  end do
405  end do
406  end do
407 
408  do j = 1, lx
409  do i = 1, lx
410  do k = 1, lx
411  do e = 1, nelv
412  wt = 0d0
413  !NEC$ unroll_completely
414  do kk = 1, lx
415  wt = wt + dz(k, kk) * u(i,j, kk,e)
416  end do
417  dudt(i,j,k,e) = wt
418  end do
419  end do
420  end do
421  end do
422 
423  do i = 1, nelv * lx * lx * lx
424  du(i,1,1,1) = jacinv(i,1,1,1) * &
425  (vx(i,1,1,1) * ( &
426  drdx(i,1,1,1)*dudr(i,1,1,1) &
427  + dsdx(i,1,1,1)*duds(i,1,1,1) &
428  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
429  + vy(i,1,1,1)*( &
430  drdy(i,1,1,1)*dudr(i,1,1,1) &
431  + dsdy(i,1,1,1)*duds(i,1,1,1) &
432  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
433  + vz(i,1,1,1)*( &
434  drdz(i,1,1,1)*dudr(i,1,1,1) &
435  + dsdz(i,1,1,1)*duds(i,1,1,1) &
436  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
437  end do
438 
439  end subroutine sx_conv1_lx12
440 
441  subroutine sx_conv1_lx11(du, u, vx, vy, vz, dx, dy, dz, &
442  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
443  jacinv, nelv)
444  integer, parameter :: lx = 11
445  integer, intent(in) :: nelv
446  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
447  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
448  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
449  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
450  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
451  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
452  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
453  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
454  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
455  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
456  real(kind=rp) :: wr, ws, wt
457  integer :: e, i, j, k, jj, kk
458 
459  do i = 1, lx
460  do jj = 1, lx * lx * nelv
461  wr = 0d0
462  do kk = 1, lx
463  wr = wr + dx(i, kk) * u(kk, jj,1,1)
464  end do
465  dudr(i, jj,1,1) = wr
466  end do
467  end do
468 
469  do k = 1, lx
470  do i = 1, lx
471  do j = 1, lx
472  do e = 1, nelv
473  ws = 0d0
474  !NEC$ unroll_completely
475  do kk = 1, lx
476  ws = ws + dy(j, kk) * u(i, kk,k,e)
477  end do
478  duds(i,j,k,e) = ws
479  end do
480  end do
481  end do
482  end do
483 
484  do j = 1, lx
485  do i = 1, lx
486  do k = 1, lx
487  do e = 1, nelv
488  wt = 0d0
489  !NEC$ unroll_completely
490  do kk = 1, lx
491  wt = wt + dz(k, kk) * u(i,j, kk,e)
492  end do
493  dudt(i,j,k,e) = wt
494  end do
495  end do
496  end do
497  end do
498 
499  do i = 1, nelv * lx * lx * lx
500  du(i,1,1,1) = jacinv(i,1,1,1) * &
501  (vx(i,1,1,1) * ( &
502  drdx(i,1,1,1)*dudr(i,1,1,1) &
503  + dsdx(i,1,1,1)*duds(i,1,1,1) &
504  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
505  + vy(i,1,1,1)*( &
506  drdy(i,1,1,1)*dudr(i,1,1,1) &
507  + dsdy(i,1,1,1)*duds(i,1,1,1) &
508  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
509  + vz(i,1,1,1)*( &
510  drdz(i,1,1,1)*dudr(i,1,1,1) &
511  + dsdz(i,1,1,1)*duds(i,1,1,1) &
512  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
513  end do
514 
515  end subroutine sx_conv1_lx11
516 
517  subroutine sx_conv1_lx10(du, u, vx, vy, vz, dx, dy, dz, &
518  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
519  jacinv, nelv)
520  integer, parameter :: lx = 10
521  integer, intent(in) :: nelv
522  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
523  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
524  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
525  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
526  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
527  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
528  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
529  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
530  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
531  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
532  real(kind=rp) :: wr, ws, wt
533  integer :: e, i, j, k, jj, kk
534 
535  do i = 1, lx
536  do jj = 1, lx * lx * nelv
537  wr = 0d0
538  do kk = 1, lx
539  wr = wr + dx(i, kk) * u(kk, jj,1,1)
540  end do
541  dudr(i, jj,1,1) = wr
542  end do
543  end do
544 
545  do k = 1, lx
546  do i = 1, lx
547  do j = 1, lx
548  do e = 1, nelv
549  ws = 0d0
550  !NEC$ unroll_completely
551  do kk = 1, lx
552  ws = ws + dy(j, kk) * u(i, kk,k,e)
553  end do
554  duds(i,j,k,e) = ws
555  end do
556  end do
557  end do
558  end do
559 
560  do j = 1, lx
561  do i = 1, lx
562  do k = 1, lx
563  do e = 1, nelv
564  wt = 0d0
565  !NEC$ unroll_completely
566  do kk = 1, lx
567  wt = wt + dz(k, kk) * u(i,j, kk,e)
568  end do
569  dudt(i,j,k,e) = wt
570  end do
571  end do
572  end do
573  end do
574 
575  do i = 1, nelv * lx * lx * lx
576  du(i,1,1,1) = jacinv(i,1,1,1) * &
577  (vx(i,1,1,1) * ( &
578  drdx(i,1,1,1)*dudr(i,1,1,1) &
579  + dsdx(i,1,1,1)*duds(i,1,1,1) &
580  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
581  + vy(i,1,1,1)*( &
582  drdy(i,1,1,1)*dudr(i,1,1,1) &
583  + dsdy(i,1,1,1)*duds(i,1,1,1) &
584  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
585  + vz(i,1,1,1)*( &
586  drdz(i,1,1,1)*dudr(i,1,1,1) &
587  + dsdz(i,1,1,1)*duds(i,1,1,1) &
588  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
589  end do
590 
591  end subroutine sx_conv1_lx10
592 
593  subroutine sx_conv1_lx9(du, u, vx, vy, vz, dx, dy, dz, &
594  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
595  jacinv, nelv)
596  integer, parameter :: lx = 9
597  integer, intent(in) :: nelv
598  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
599  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
600  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
601  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
602  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
603  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
604  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
605  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
606  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
607  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
608  real(kind=rp) :: wr, ws, wt
609  integer :: e, i, j, k, jj, kk
610 
611  do i = 1, lx
612  do jj = 1, lx * lx * nelv
613  wr = 0d0
614  do kk = 1, lx
615  wr = wr + dx(i, kk) * u(kk, jj,1,1)
616  end do
617  dudr(i, jj,1,1) = wr
618  end do
619  end do
620 
621  do k = 1, lx
622  do i = 1, lx
623  do j = 1, lx
624  do e = 1, nelv
625  ws = 0d0
626  !NEC$ unroll_completely
627  do kk = 1, lx
628  ws = ws + dy(j, kk) * u(i, kk,k,e)
629  end do
630  duds(i,j,k,e) = ws
631  end do
632  end do
633  end do
634  end do
635 
636  do j = 1, lx
637  do i = 1, lx
638  do k = 1, lx
639  do e = 1, nelv
640  wt = 0d0
641  !NEC$ unroll_completely
642  do kk = 1, lx
643  wt = wt + dz(k, kk) * u(i,j, kk,e)
644  end do
645  dudt(i,j,k,e) = wt
646  end do
647  end do
648  end do
649  end do
650 
651  do i = 1, nelv * lx * lx * lx
652  du(i,1,1,1) = jacinv(i,1,1,1) * &
653  (vx(i,1,1,1) * ( &
654  drdx(i,1,1,1)*dudr(i,1,1,1) &
655  + dsdx(i,1,1,1)*duds(i,1,1,1) &
656  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
657  + vy(i,1,1,1)*( &
658  drdy(i,1,1,1)*dudr(i,1,1,1) &
659  + dsdy(i,1,1,1)*duds(i,1,1,1) &
660  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
661  + vz(i,1,1,1)*( &
662  drdz(i,1,1,1)*dudr(i,1,1,1) &
663  + dsdz(i,1,1,1)*duds(i,1,1,1) &
664  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
665  end do
666 
667  end subroutine sx_conv1_lx9
668 
669  subroutine sx_conv1_lx8(du, u, vx, vy, vz, dx, dy, dz, &
670  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
671  jacinv, nelv)
672  integer, parameter :: lx = 8
673  integer, intent(in) :: nelv
674  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
675  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
676  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
677  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
678  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
679  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
680  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
681  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
682  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
683  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
684  real(kind=rp) :: wr, ws, wt
685  integer :: e, i, j, k, jj, kk
686 
687  do i = 1, lx
688  do jj = 1, lx * lx * nelv
689  wr = 0d0
690  do kk = 1, lx
691  wr = wr + dx(i, kk) * u(kk, jj,1,1)
692  end do
693  dudr(i, jj,1,1) = wr
694  end do
695  end do
696 
697  do k = 1, lx
698  do i = 1, lx
699  do j = 1, lx
700  do e = 1, nelv
701  ws = 0d0
702  !NEC$ unroll_completely
703  do kk = 1, lx
704  ws = ws + dy(j, kk) * u(i, kk,k,e)
705  end do
706  duds(i,j,k,e) = ws
707  end do
708  end do
709  end do
710  end do
711 
712  do j = 1, lx
713  do i = 1, lx
714  do k = 1, lx
715  do e = 1, nelv
716  wt = 0d0
717  !NEC$ unroll_completely
718  do kk = 1, lx
719  wt = wt + dz(k, kk) * u(i,j, kk,e)
720  end do
721  dudt(i,j,k,e) = wt
722  end do
723  end do
724  end do
725  end do
726 
727  do i = 1, nelv * lx * lx * lx
728  du(i,1,1,1) = jacinv(i,1,1,1) * &
729  (vx(i,1,1,1) * ( &
730  drdx(i,1,1,1)*dudr(i,1,1,1) &
731  + dsdx(i,1,1,1)*duds(i,1,1,1) &
732  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
733  + vy(i,1,1,1)*( &
734  drdy(i,1,1,1)*dudr(i,1,1,1) &
735  + dsdy(i,1,1,1)*duds(i,1,1,1) &
736  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
737  + vz(i,1,1,1)*( &
738  drdz(i,1,1,1)*dudr(i,1,1,1) &
739  + dsdz(i,1,1,1)*duds(i,1,1,1) &
740  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
741  end do
742 
743  end subroutine sx_conv1_lx8
744 
745  subroutine sx_conv1_lx7(du, u, vx, vy, vz, dx, dy, dz, &
746  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
747  jacinv, nelv)
748  integer, parameter :: lx = 7
749  integer, intent(in) :: nelv
750  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
751  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
752  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
753  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
754  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
755  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
756  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
757  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
758  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
759  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
760  real(kind=rp) :: wr, ws, wt
761  integer :: e, i, j, k, jj, kk
762 
763  do i = 1, lx
764  do jj = 1, lx * lx * nelv
765  wr = 0d0
766  do kk = 1, lx
767  wr = wr + dx(i, kk) * u(kk, jj,1,1)
768  end do
769  dudr(i, jj,1,1) = wr
770  end do
771  end do
772 
773  do k = 1, lx
774  do i = 1, lx
775  do j = 1, lx
776  do e = 1, nelv
777  ws = 0d0
778  !NEC$ unroll_completely
779  do kk = 1, lx
780  ws = ws + dy(j, kk) * u(i, kk,k,e)
781  end do
782  duds(i,j,k,e) = ws
783  end do
784  end do
785  end do
786  end do
787 
788  do j = 1, lx
789  do i = 1, lx
790  do k = 1, lx
791  do e = 1, nelv
792  wt = 0d0
793  !NEC$ unroll_completely
794  do kk = 1, lx
795  wt = wt + dz(k, kk) * u(i,j, kk,e)
796  end do
797  dudt(i,j,k,e) = wt
798  end do
799  end do
800  end do
801  end do
802 
803  do i = 1, nelv * lx * lx * lx
804  du(i,1,1,1) = jacinv(i,1,1,1) * &
805  (vx(i,1,1,1) * ( &
806  drdx(i,1,1,1)*dudr(i,1,1,1) &
807  + dsdx(i,1,1,1)*duds(i,1,1,1) &
808  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
809  + vy(i,1,1,1)*( &
810  drdy(i,1,1,1)*dudr(i,1,1,1) &
811  + dsdy(i,1,1,1)*duds(i,1,1,1) &
812  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
813  + vz(i,1,1,1)*( &
814  drdz(i,1,1,1)*dudr(i,1,1,1) &
815  + dsdz(i,1,1,1)*duds(i,1,1,1) &
816  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
817  end do
818 
819  end subroutine sx_conv1_lx7
820 
821  subroutine sx_conv1_lx6(du, u, vx, vy, vz, dx, dy, dz, &
822  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
823  jacinv, nelv)
824  integer, parameter :: lx = 6
825  integer, intent(in) :: nelv
826  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
827  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
828  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
829  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
830  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
831  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
832  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
833  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
834  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
835  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
836  real(kind=rp) :: wr, ws, wt
837  integer :: e, i, j, k, jj, kk
838 
839  do i = 1, lx
840  do jj = 1, lx * lx * nelv
841  wr = 0d0
842  do kk = 1, lx
843  wr = wr + dx(i, kk) * u(kk, jj,1,1)
844  end do
845  dudr(i, jj,1,1) = wr
846  end do
847  end do
848 
849  do k = 1, lx
850  do i = 1, lx
851  do j = 1, lx
852  do e = 1, nelv
853  ws = 0d0
854  !NEC$ unroll_completely
855  do kk = 1, lx
856  ws = ws + dy(j, kk) * u(i, kk,k,e)
857  end do
858  duds(i,j,k,e) = ws
859  end do
860  end do
861  end do
862  end do
863 
864  do j = 1, lx
865  do i = 1, lx
866  do k = 1, lx
867  do e = 1, nelv
868  wt = 0d0
869  !NEC$ unroll_completely
870  do kk = 1, lx
871  wt = wt + dz(k, kk) * u(i,j, kk,e)
872  end do
873  dudt(i,j,k,e) = wt
874  end do
875  end do
876  end do
877  end do
878 
879  do i = 1, nelv * lx * lx * lx
880  du(i,1,1,1) = jacinv(i,1,1,1) * &
881  (vx(i,1,1,1) * ( &
882  drdx(i,1,1,1)*dudr(i,1,1,1) &
883  + dsdx(i,1,1,1)*duds(i,1,1,1) &
884  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
885  + vy(i,1,1,1)*( &
886  drdy(i,1,1,1)*dudr(i,1,1,1) &
887  + dsdy(i,1,1,1)*duds(i,1,1,1) &
888  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
889  + vz(i,1,1,1)*( &
890  drdz(i,1,1,1)*dudr(i,1,1,1) &
891  + dsdz(i,1,1,1)*duds(i,1,1,1) &
892  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
893  end do
894 
895  end subroutine sx_conv1_lx6
896 
897  subroutine sx_conv1_lx5(du, u, vx, vy, vz, dx, dy, dz, &
898  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
899  jacinv, nelv)
900  integer, parameter :: lx = 5
901  integer, intent(in) :: nelv
902  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
903  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
904  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
905  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
906  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
907  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
908  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
909  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
910  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
911  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
912  real(kind=rp) :: wr, ws, wt
913  integer :: e, i, j, k, jj, kk
914 
915  do i = 1, lx
916  do jj = 1, lx * lx * nelv
917  wr = 0d0
918  do kk = 1, lx
919  wr = wr + dx(i, kk) * u(kk, jj,1,1)
920  end do
921  dudr(i, jj,1,1) = wr
922  end do
923  end do
924 
925  do k = 1, lx
926  do i = 1, lx
927  do j = 1, lx
928  do e = 1, nelv
929  ws = 0d0
930  !NEC$ unroll_completely
931  do kk = 1, lx
932  ws = ws + dy(j, kk) * u(i, kk,k,e)
933  end do
934  duds(i,j,k,e) = ws
935  end do
936  end do
937  end do
938  end do
939 
940  do j = 1, lx
941  do i = 1, lx
942  do k = 1, lx
943  do e = 1, nelv
944  wt = 0d0
945  !NEC$ unroll_completely
946  do kk = 1, lx
947  wt = wt + dz(k, kk) * u(i,j, kk,e)
948  end do
949  dudt(i,j,k,e) = wt
950  end do
951  end do
952  end do
953  end do
954 
955  do i = 1, nelv * lx * lx * lx
956  du(i,1,1,1) = jacinv(i,1,1,1) * &
957  (vx(i,1,1,1) * ( &
958  drdx(i,1,1,1)*dudr(i,1,1,1) &
959  + dsdx(i,1,1,1)*duds(i,1,1,1) &
960  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
961  + vy(i,1,1,1)*( &
962  drdy(i,1,1,1)*dudr(i,1,1,1) &
963  + dsdy(i,1,1,1)*duds(i,1,1,1) &
964  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
965  + vz(i,1,1,1)*( &
966  drdz(i,1,1,1)*dudr(i,1,1,1) &
967  + dsdz(i,1,1,1)*duds(i,1,1,1) &
968  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
969  end do
970 
971  end subroutine sx_conv1_lx5
972 
973  subroutine sx_conv1_lx4(du, u, vx, vy, vz, dx, dy, dz, &
974  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
975  jacinv, nelv)
976  integer, parameter :: lx = 4
977  integer, intent(in) :: nelv
978  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
979  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
980  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
981  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
982  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
983  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
984  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
985  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
986  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
987  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
988  real(kind=rp) :: wr, ws, wt
989  integer :: e, i, j, k, jj, kk
990 
991  do i = 1, lx
992  do jj = 1, lx * lx * nelv
993  wr = 0d0
994  do kk = 1, lx
995  wr = wr + dx(i, kk) * u(kk, jj,1,1)
996  end do
997  dudr(i, jj,1,1) = wr
998  end do
999  end do
1000 
1001  do k = 1, lx
1002  do i = 1, lx
1003  do j = 1, lx
1004  do e = 1, nelv
1005  ws = 0d0
1006  !NEC$ unroll_completely
1007  do kk = 1, lx
1008  ws = ws + dy(j, kk) * u(i, kk,k,e)
1009  end do
1010  duds(i,j,k,e) = ws
1011  end do
1012  end do
1013  end do
1014  end do
1015 
1016  do j = 1, lx
1017  do i = 1, lx
1018  do k = 1, lx
1019  do e = 1, nelv
1020  wt = 0d0
1021  !NEC$ unroll_completely
1022  do kk = 1, lx
1023  wt = wt + dz(k, kk) * u(i,j, kk,e)
1024  end do
1025  dudt(i,j,k,e) = wt
1026  end do
1027  end do
1028  end do
1029  end do
1030 
1031  do i = 1, nelv * lx * lx * lx
1032  du(i,1,1,1) = jacinv(i,1,1,1) * &
1033  (vx(i,1,1,1) * ( &
1034  drdx(i,1,1,1)*dudr(i,1,1,1) &
1035  + dsdx(i,1,1,1)*duds(i,1,1,1) &
1036  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
1037  + vy(i,1,1,1)*( &
1038  drdy(i,1,1,1)*dudr(i,1,1,1) &
1039  + dsdy(i,1,1,1)*duds(i,1,1,1) &
1040  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
1041  + vz(i,1,1,1)*( &
1042  drdz(i,1,1,1)*dudr(i,1,1,1) &
1043  + dsdz(i,1,1,1)*duds(i,1,1,1) &
1044  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
1045  end do
1046 
1047  end subroutine sx_conv1_lx4
1048 
1049  subroutine sx_conv1_lx3(du, u, vx, vy, vz, dx, dy, dz, &
1050  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1051  jacinv, nelv)
1052  integer, parameter :: lx = 3
1053  integer, intent(in) :: nelv
1054  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
1055  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
1056  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
1057  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
1058  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
1059  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
1060  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1061  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
1062  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
1063  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
1064  real(kind=rp) :: wr, ws, wt
1065  integer :: e, i, j, k, jj, kk
1066 
1067  do i = 1, lx
1068  do jj = 1, lx * lx * nelv
1069  wr = 0d0
1070  do kk = 1, lx
1071  wr = wr + dx(i, kk) * u(kk, jj,1,1)
1072  end do
1073  dudr(i, jj,1,1) = wr
1074  end do
1075  end do
1076 
1077  do k = 1, lx
1078  do i = 1, lx
1079  do j = 1, lx
1080  do e = 1, nelv
1081  ws = 0d0
1082  !NEC$ unroll_completely
1083  do kk = 1, lx
1084  ws = ws + dy(j, kk) * u(i, kk,k,e)
1085  end do
1086  duds(i,j,k,e) = ws
1087  end do
1088  end do
1089  end do
1090  end do
1091 
1092  do j = 1, lx
1093  do i = 1, lx
1094  do k = 1, lx
1095  do e = 1, nelv
1096  wt = 0d0
1097  !NEC$ unroll_completely
1098  do kk = 1, lx
1099  wt = wt + dz(k, kk) * u(i,j, kk,e)
1100  end do
1101  dudt(i,j,k,e) = wt
1102  end do
1103  end do
1104  end do
1105  end do
1106 
1107  do i = 1, nelv * lx * lx * lx
1108  du(i,1,1,1) = jacinv(i,1,1,1) * &
1109  (vx(i,1,1,1) * ( &
1110  drdx(i,1,1,1)*dudr(i,1,1,1) &
1111  + dsdx(i,1,1,1)*duds(i,1,1,1) &
1112  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
1113  + vy(i,1,1,1)*( &
1114  drdy(i,1,1,1)*dudr(i,1,1,1) &
1115  + dsdy(i,1,1,1)*duds(i,1,1,1) &
1116  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
1117  + vz(i,1,1,1)*( &
1118  drdz(i,1,1,1)*dudr(i,1,1,1) &
1119  + dsdz(i,1,1,1)*duds(i,1,1,1) &
1120  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
1121  end do
1122 
1123  end subroutine sx_conv1_lx3
1124 
1125  subroutine sx_conv1_lx2(du, u, vx, vy, vz, dx, dy, dz, &
1126  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1127  jacinv, nelv)
1128  integer, parameter :: lx = 2
1129  integer, intent(in) :: nelv
1130  real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
1131  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
1132  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
1133  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
1134  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
1135  real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
1136  real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1137  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudr
1138  real(kind=rp), dimension(lx, lx, lx, nelv) :: duds
1139  real(kind=rp), dimension(lx, lx, lx, nelv) :: dudt
1140  real(kind=rp) :: wr, ws, wt
1141  integer :: e, i, j, k, jj, kk
1142 
1143  do i = 1, lx
1144  do jj = 1, lx * lx * nelv
1145  wr = 0d0
1146  do kk = 1, lx
1147  wr = wr + dx(i, kk) * u(kk, jj,1,1)
1148  end do
1149  dudr(i, jj,1,1) = wr
1150  end do
1151  end do
1152 
1153  do k = 1, lx
1154  do i = 1, lx
1155  do j = 1, lx
1156  do e = 1, nelv
1157  ws = 0d0
1158  !NEC$ unroll_completely
1159  do kk = 1, lx
1160  ws = ws + dy(j, kk) * u(i, kk,k,e)
1161  end do
1162  duds(i,j,k,e) = ws
1163  end do
1164  end do
1165  end do
1166  end do
1167 
1168  do j = 1, lx
1169  do i = 1, lx
1170  do k = 1, lx
1171  do e = 1, nelv
1172  wt = 0d0
1173  !NEC$ unroll_completely
1174  do kk = 1, lx
1175  wt = wt + dz(k, kk) * u(i,j, kk,e)
1176  end do
1177  dudt(i,j,k,e) = wt
1178  end do
1179  end do
1180  end do
1181  end do
1182 
1183  do i = 1, nelv * lx * lx * lx
1184  du(i,1,1,1) = jacinv(i,1,1,1) * &
1185  (vx(i,1,1,1) * ( &
1186  drdx(i,1,1,1)*dudr(i,1,1,1) &
1187  + dsdx(i,1,1,1)*duds(i,1,1,1) &
1188  + dtdx(i,1,1,1)*dudt(i,1,1,1)) &
1189  + vy(i,1,1,1)*( &
1190  drdy(i,1,1,1)*dudr(i,1,1,1) &
1191  + dsdy(i,1,1,1)*duds(i,1,1,1) &
1192  + dtdy(i,1,1,1)*dudt(i,1,1,1)) &
1193  + vz(i,1,1,1)*( &
1194  drdz(i,1,1,1)*dudr(i,1,1,1) &
1195  + dsdz(i,1,1,1)*duds(i,1,1,1) &
1196  + dtdz(i,1,1,1)*dudt(i,1,1,1)))
1197  end do
1198 
1199  end subroutine sx_conv1_lx2
1200 
1201 end submodule sx_conv1
Operators SX-Aurora backend.
Definition: opr_sx.f90:2