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