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