39 module subroutine opr_cpu_conv1(du, u, vx, vy, vz, xh, coef, e_start, e_end)
40 type(space_t),
intent(in) :: Xh
41 type(coef_t),
intent(in) :: coef
42 integer,
intent(in) :: e_start, e_end
43 real(kind=rp),
intent(inout) :: du(xh%lxyz, e_end-e_start+1)
44 real(kind=rp),
intent(inout) :: u(xh%lx, xh%ly, xh%lz, e_end-e_start+1)
45 real(kind=rp),
intent(inout) :: vx(xh%lx, xh%ly, xh%lz, e_end-e_start+1)
46 real(kind=rp),
intent(inout) :: vy(xh%lx, xh%ly, xh%lz, e_end-e_start+1)
47 real(kind=rp),
intent(inout) :: vz(xh%lx, xh%ly, xh%lz, e_end-e_start+1)
50 e_len = e_end-e_start+1
52 if (e_len .eq. 1)
then
53 call opr_cpu_conv1_single(du, u, vx, vy, vz, xh, coef, e_start)
55 call opr_cpu_conv1_many(du, u, vx, vy, vz, xh, coef, e_start, e_len)
57 end subroutine opr_cpu_conv1
59 subroutine opr_cpu_conv1_many(du, u, vx, vy, vz, Xh, coef, e_start, e_len)
60 type(space_t),
intent(in) :: Xh
61 type(coef_t),
intent(in) :: coef
62 integer,
intent(in) :: e_start, e_len
63 real(kind=rp),
intent(inout) :: du(xh%lxyz, e_len)
64 real(kind=rp),
intent(inout) :: u(xh%lx, xh%ly, xh%lz, e_len)
65 real(kind=rp),
intent(inout) :: vx(xh%lx, xh%ly, xh%lz, e_len)
66 real(kind=rp),
intent(inout) :: vy(xh%lx, xh%ly, xh%lz, e_len)
67 real(kind=rp),
intent(inout) :: vz(xh%lx, xh%ly, xh%lz, e_len)
69 associate(drdx => coef%drdx, drdy => coef%drdy, drdz => coef%drdz, &
70 dsdx => coef%dsdx, dsdy => coef%dsdy, dsdz => coef%dsdz, &
71 dtdx => coef%dtdx, dtdy => coef%dtdy, dtdz => coef%dtdz, &
72 jacinv => coef%jacinv)
76 call cpu_conv1_lx14(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
77 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
78 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
79 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
80 jacinv(1,1,1, e_start), e_len)
82 call cpu_conv1_lx13(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
83 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
84 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
85 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
86 jacinv(1,1,1, e_start), e_len)
88 call cpu_conv1_lx12(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
89 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
90 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
91 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
92 jacinv(1,1,1, e_start), e_len)
94 call cpu_conv1_lx11(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
95 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
96 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
97 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
98 jacinv(1,1,1, e_start), e_len)
100 call cpu_conv1_lx10(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
101 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
102 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
103 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
104 jacinv(1,1,1, e_start), e_len)
106 call cpu_conv1_lx9(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
107 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
108 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
109 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
110 jacinv(1,1,1, e_start), e_len)
112 call cpu_conv1_lx8(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
113 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
114 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
115 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
116 jacinv(1,1,1, e_start), e_len)
118 call cpu_conv1_lx7(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
119 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
120 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
121 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
122 jacinv(1,1,1, e_start), e_len)
124 call cpu_conv1_lx6(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
125 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
126 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
127 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
128 jacinv(1,1,1, e_start), e_len)
130 call cpu_conv1_lx5(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
131 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
132 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
133 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
134 jacinv(1,1,1, e_start), e_len)
136 call cpu_conv1_lx4(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
137 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
138 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
139 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
140 jacinv(1,1,1, e_start), e_len)
142 call cpu_conv1_lx3(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
143 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
144 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
145 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
146 jacinv(1,1,1, e_start), e_len)
148 call cpu_conv1_lx2(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
149 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
150 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
151 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
152 jacinv(1,1,1, e_start), e_len)
154 call cpu_conv1_lx(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
155 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
156 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
157 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
158 jacinv(1,1,1, e_start), e_len, xh%lx)
163 end subroutine opr_cpu_conv1_many
165 subroutine cpu_conv1_lx(du, u, vx, vy, vz, dx, dy, dz, &
166 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
168 integer,
intent(in) :: nelv, lx
169 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
170 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
171 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
172 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
173 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
174 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
175 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
176 real(kind=rp),
dimension(lx, lx, lx) :: dudr
177 real(kind=rp),
dimension(lx, lx, lx) :: duds
178 real(kind=rp),
dimension(lx, lx, lx) :: dudt
180 integer :: e, i, j, k, l
187 tmp = tmp + dx(i,k) * u(k,j,1,e)
198 tmp = tmp + dy(j,l) * u(i,l,k,e)
209 tmp = tmp + dz(k,l) * u(i,1,l,e)
215 do i = 1, lx * lx * lx
216 du(i,1,1,e) = jacinv(i,1,1,e) &
218 * ( drdx(i,1,1,e) * dudr(i,1,1) &
219 + dsdx(i,1,1,e) * duds(i,1,1) &
220 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
222 * ( drdy(i,1,1,e) * dudr(i,1,1) &
223 + dsdy(i,1,1,e) * duds(i,1,1) &
224 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
226 * ( drdz(i,1,1,e) * dudr(i,1,1) &
227 + dsdz(i,1,1,e) * duds(i,1,1) &
228 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
232 end subroutine cpu_conv1_lx
234 subroutine cpu_conv1_lx14(du, u, vx, vy, vz, dx, dy, dz, &
235 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
237 integer,
parameter :: lx = 14
238 integer,
intent(in) :: nelv
239 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
240 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
241 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
242 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
243 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
244 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
245 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
246 real(kind=rp),
dimension(lx, lx, lx) :: dudr
247 real(kind=rp),
dimension(lx, lx, lx) :: duds
248 real(kind=rp),
dimension(lx, lx, lx) :: dudt
249 integer :: e, i, j, k
254 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
255 + dx(i,2) * u(2,j,1,e) &
256 + dx(i,3) * u(3,j,1,e) &
257 + dx(i,4) * u(4,j,1,e) &
258 + dx(i,5) * u(5,j,1,e) &
259 + dx(i,6) * u(6,j,1,e) &
260 + dx(i,7) * u(7,j,1,e) &
261 + dx(i,8) * u(8,j,1,e) &
262 + dx(i,9) * u(9,j,1,e) &
263 + dx(i,10) * u(10,j,1,e) &
264 + dx(i,11) * u(11,j,1,e) &
265 + dx(i,12) * u(12,j,1,e) &
266 + dx(i,13) * u(13,j,1,e) &
267 + dx(i,14) * u(14,j,1,e)
274 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
275 + dy(j,2) * u(i,2,k,e) &
276 + dy(j,3) * u(i,3,k,e) &
277 + dy(j,4) * u(i,4,k,e) &
278 + dy(j,5) * u(i,5,k,e) &
279 + dy(j,6) * u(i,6,k,e) &
280 + dy(j,7) * u(i,7,k,e) &
281 + dy(j,8) * u(i,8,k,e) &
282 + dy(j,9) * u(i,9,k,e) &
283 + dy(j,10) * u(i,10,k,e) &
284 + dy(j,11) * u(i,11,k,e) &
285 + dy(j,12) * u(i,12,k,e) &
286 + dy(j,13) * u(i,13,k,e) &
287 + dy(j,14) * u(i,14,k,e)
294 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
295 + dz(k,2) * u(i,1,2,e) &
296 + dz(k,3) * u(i,1,3,e) &
297 + dz(k,4) * u(i,1,4,e) &
298 + dz(k,5) * u(i,1,5,e) &
299 + dz(k,6) * u(i,1,6,e) &
300 + dz(k,7) * u(i,1,7,e) &
301 + dz(k,8) * u(i,1,8,e) &
302 + dz(k,9) * u(i,1,9,e) &
303 + dz(k,10) * u(i,1,10,e) &
304 + dz(k,11) * u(i,1,11,e) &
305 + dz(k,12) * u(i,1,12,e) &
306 + dz(k,13) * u(i,1,13,e) &
307 + dz(k,14) * u(i,1,14,e)
311 do i = 1, lx * lx * lx
312 du(i,1,1,e) = jacinv(i,1,1,e) &
314 * ( drdx(i,1,1,e) * dudr(i,1,1) &
315 + dsdx(i,1,1,e) * duds(i,1,1) &
316 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
318 * ( drdy(i,1,1,e) * dudr(i,1,1) &
319 + dsdy(i,1,1,e) * duds(i,1,1) &
320 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
322 * ( drdz(i,1,1,e) * dudr(i,1,1) &
323 + dsdz(i,1,1,e) * duds(i,1,1) &
324 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
328 end subroutine cpu_conv1_lx14
330 subroutine cpu_conv1_lx13(du, u, vx, vy, vz, dx, dy, dz, &
331 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
333 integer,
parameter :: lx = 13
334 integer,
intent(in) :: nelv
335 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
336 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
337 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
338 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
339 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
340 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
341 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
342 real(kind=rp),
dimension(lx, lx, lx) :: dudr
343 real(kind=rp),
dimension(lx, lx, lx) :: duds
344 real(kind=rp),
dimension(lx, lx, lx) :: dudt
345 integer :: e, i, j, k
350 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
351 + dx(i,2) * u(2,j,1,e) &
352 + dx(i,3) * u(3,j,1,e) &
353 + dx(i,4) * u(4,j,1,e) &
354 + dx(i,5) * u(5,j,1,e) &
355 + dx(i,6) * u(6,j,1,e) &
356 + dx(i,7) * u(7,j,1,e) &
357 + dx(i,8) * u(8,j,1,e) &
358 + dx(i,9) * u(9,j,1,e) &
359 + dx(i,10) * u(10,j,1,e) &
360 + dx(i,11) * u(11,j,1,e) &
361 + dx(i,12) * u(12,j,1,e) &
362 + dx(i,13) * u(13,j,1,e)
369 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
370 + dy(j,2) * u(i,2,k,e) &
371 + dy(j,3) * u(i,3,k,e) &
372 + dy(j,4) * u(i,4,k,e) &
373 + dy(j,5) * u(i,5,k,e) &
374 + dy(j,6) * u(i,6,k,e) &
375 + dy(j,7) * u(i,7,k,e) &
376 + dy(j,8) * u(i,8,k,e) &
377 + dy(j,9) * u(i,9,k,e) &
378 + dy(j,10) * u(i,10,k,e) &
379 + dy(j,11) * u(i,11,k,e) &
380 + dy(j,12) * u(i,12,k,e) &
381 + dy(j,13) * u(i,13,k,e)
388 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
389 + dz(k,2) * u(i,1,2,e) &
390 + dz(k,3) * u(i,1,3,e) &
391 + dz(k,4) * u(i,1,4,e) &
392 + dz(k,5) * u(i,1,5,e) &
393 + dz(k,6) * u(i,1,6,e) &
394 + dz(k,7) * u(i,1,7,e) &
395 + dz(k,8) * u(i,1,8,e) &
396 + dz(k,9) * u(i,1,9,e) &
397 + dz(k,10) * u(i,1,10,e) &
398 + dz(k,11) * u(i,1,11,e) &
399 + dz(k,12) * u(i,1,12,e) &
400 + dz(k,13) * u(i,1,13,e)
404 do i = 1, lx * lx * lx
405 du(i,1,1,e) = jacinv(i,1,1,e) &
407 * ( drdx(i,1,1,e) * dudr(i,1,1) &
408 + dsdx(i,1,1,e) * duds(i,1,1) &
409 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
411 * ( drdy(i,1,1,e) * dudr(i,1,1) &
412 + dsdy(i,1,1,e) * duds(i,1,1) &
413 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
415 * ( drdz(i,1,1,e) * dudr(i,1,1) &
416 + dsdz(i,1,1,e) * duds(i,1,1) &
417 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
421 end subroutine cpu_conv1_lx13
423 subroutine cpu_conv1_lx12(du, u, vx, vy, vz, dx, dy, dz, &
424 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
426 integer,
parameter :: lx = 12
427 integer,
intent(in) :: nelv
428 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
429 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
430 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
431 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
432 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
433 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
434 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
435 real(kind=rp),
dimension(lx, lx, lx) :: dudr
436 real(kind=rp),
dimension(lx, lx, lx) :: duds
437 real(kind=rp),
dimension(lx, lx, lx) :: dudt
438 integer :: e, i, j, k
443 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
444 + dx(i,2) * u(2,j,1,e) &
445 + dx(i,3) * u(3,j,1,e) &
446 + dx(i,4) * u(4,j,1,e) &
447 + dx(i,5) * u(5,j,1,e) &
448 + dx(i,6) * u(6,j,1,e) &
449 + dx(i,7) * u(7,j,1,e) &
450 + dx(i,8) * u(8,j,1,e) &
451 + dx(i,9) * u(9,j,1,e) &
452 + dx(i,10) * u(10,j,1,e) &
453 + dx(i,11) * u(11,j,1,e) &
454 + dx(i,12) * u(12,j,1,e)
461 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
462 + dy(j,2) * u(i,2,k,e) &
463 + dy(j,3) * u(i,3,k,e) &
464 + dy(j,4) * u(i,4,k,e) &
465 + dy(j,5) * u(i,5,k,e) &
466 + dy(j,6) * u(i,6,k,e) &
467 + dy(j,7) * u(i,7,k,e) &
468 + dy(j,8) * u(i,8,k,e) &
469 + dy(j,9) * u(i,9,k,e) &
470 + dy(j,10) * u(i,10,k,e) &
471 + dy(j,11) * u(i,11,k,e) &
472 + dy(j,12) * u(i,12,k,e)
479 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
480 + dz(k,2) * u(i,1,2,e) &
481 + dz(k,3) * u(i,1,3,e) &
482 + dz(k,4) * u(i,1,4,e) &
483 + dz(k,5) * u(i,1,5,e) &
484 + dz(k,6) * u(i,1,6,e) &
485 + dz(k,7) * u(i,1,7,e) &
486 + dz(k,8) * u(i,1,8,e) &
487 + dz(k,9) * u(i,1,9,e) &
488 + dz(k,10) * u(i,1,10,e) &
489 + dz(k,11) * u(i,1,11,e) &
490 + dz(k,12) * u(i,1,12,e)
494 do i = 1, lx * lx * lx
495 du(i,1,1,e) = jacinv(i,1,1,e) &
497 * ( drdx(i,1,1,e) * dudr(i,1,1) &
498 + dsdx(i,1,1,e) * duds(i,1,1) &
499 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
501 * ( drdy(i,1,1,e) * dudr(i,1,1) &
502 + dsdy(i,1,1,e) * duds(i,1,1) &
503 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
505 * ( drdz(i,1,1,e) * dudr(i,1,1) &
506 + dsdz(i,1,1,e) * duds(i,1,1) &
507 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
511 end subroutine cpu_conv1_lx12
513 subroutine cpu_conv1_lx11(du, u, vx, vy, vz, dx, dy, dz, &
514 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
516 integer,
parameter :: lx = 11
517 integer,
intent(in) :: nelv
518 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
519 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
520 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
521 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
522 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
523 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
524 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
525 real(kind=rp),
dimension(lx, lx, lx) :: dudr
526 real(kind=rp),
dimension(lx, lx, lx) :: duds
527 real(kind=rp),
dimension(lx, lx, lx) :: dudt
528 integer :: e, i, j, k
533 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
534 + dx(i,2) * u(2,j,1,e) &
535 + dx(i,3) * u(3,j,1,e) &
536 + dx(i,4) * u(4,j,1,e) &
537 + dx(i,5) * u(5,j,1,e) &
538 + dx(i,6) * u(6,j,1,e) &
539 + dx(i,7) * u(7,j,1,e) &
540 + dx(i,8) * u(8,j,1,e) &
541 + dx(i,9) * u(9,j,1,e) &
542 + dx(i,10) * u(10,j,1,e) &
543 + dx(i,11) * u(11,j,1,e)
550 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
551 + dy(j,2) * u(i,2,k,e) &
552 + dy(j,3) * u(i,3,k,e) &
553 + dy(j,4) * u(i,4,k,e) &
554 + dy(j,5) * u(i,5,k,e) &
555 + dy(j,6) * u(i,6,k,e) &
556 + dy(j,7) * u(i,7,k,e) &
557 + dy(j,8) * u(i,8,k,e) &
558 + dy(j,9) * u(i,9,k,e) &
559 + dy(j,10) * u(i,10,k,e) &
560 + dy(j,11) * u(i,11,k,e)
567 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
568 + dz(k,2) * u(i,1,2,e) &
569 + dz(k,3) * u(i,1,3,e) &
570 + dz(k,4) * u(i,1,4,e) &
571 + dz(k,5) * u(i,1,5,e) &
572 + dz(k,6) * u(i,1,6,e) &
573 + dz(k,7) * u(i,1,7,e) &
574 + dz(k,8) * u(i,1,8,e) &
575 + dz(k,9) * u(i,1,9,e) &
576 + dz(k,10) * u(i,1,10,e) &
577 + dz(k,11) * u(i,1,11,e)
581 do i = 1, lx * lx * lx
582 du(i,1,1,e) = jacinv(i,1,1,e) &
584 * ( drdx(i,1,1,e) * dudr(i,1,1) &
585 + dsdx(i,1,1,e) * duds(i,1,1) &
586 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
588 * ( drdy(i,1,1,e) * dudr(i,1,1) &
589 + dsdy(i,1,1,e) * duds(i,1,1) &
590 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
592 * ( drdz(i,1,1,e) * dudr(i,1,1) &
593 + dsdz(i,1,1,e) * duds(i,1,1) &
594 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
598 end subroutine cpu_conv1_lx11
600 subroutine cpu_conv1_lx10(du, u, vx, vy, vz, dx, dy, dz, &
601 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
603 integer,
parameter :: lx = 10
604 integer,
intent(in) :: nelv
605 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
606 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
607 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
608 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
609 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
610 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
611 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
612 real(kind=rp),
dimension(lx, lx, lx) :: dudr
613 real(kind=rp),
dimension(lx, lx, lx) :: duds
614 real(kind=rp),
dimension(lx, lx, lx) :: dudt
615 integer :: e, i, j, k
620 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
621 + dx(i,2) * u(2,j,1,e) &
622 + dx(i,3) * u(3,j,1,e) &
623 + dx(i,4) * u(4,j,1,e) &
624 + dx(i,5) * u(5,j,1,e) &
625 + dx(i,6) * u(6,j,1,e) &
626 + dx(i,7) * u(7,j,1,e) &
627 + dx(i,8) * u(8,j,1,e) &
628 + dx(i,9) * u(9,j,1,e) &
629 + dx(i,10) * u(10,j,1,e)
636 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
637 + dy(j,2) * u(i,2,k,e) &
638 + dy(j,3) * u(i,3,k,e) &
639 + dy(j,4) * u(i,4,k,e) &
640 + dy(j,5) * u(i,5,k,e) &
641 + dy(j,6) * u(i,6,k,e) &
642 + dy(j,7) * u(i,7,k,e) &
643 + dy(j,8) * u(i,8,k,e) &
644 + dy(j,9) * u(i,9,k,e) &
645 + dy(j,10) * u(i,10,k,e)
652 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
653 + dz(k,2) * u(i,1,2,e) &
654 + dz(k,3) * u(i,1,3,e) &
655 + dz(k,4) * u(i,1,4,e) &
656 + dz(k,5) * u(i,1,5,e) &
657 + dz(k,6) * u(i,1,6,e) &
658 + dz(k,7) * u(i,1,7,e) &
659 + dz(k,8) * u(i,1,8,e) &
660 + dz(k,9) * u(i,1,9,e) &
661 + dz(k,10) * u(i,1,10,e)
665 do i = 1, lx * lx * lx
666 du(i,1,1,e) = jacinv(i,1,1,e) &
668 * ( drdx(i,1,1,e) * dudr(i,1,1) &
669 + dsdx(i,1,1,e) * duds(i,1,1) &
670 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
672 * ( drdy(i,1,1,e) * dudr(i,1,1) &
673 + dsdy(i,1,1,e) * duds(i,1,1) &
674 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
676 * ( drdz(i,1,1,e) * dudr(i,1,1) &
677 + dsdz(i,1,1,e) * duds(i,1,1) &
678 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
682 end subroutine cpu_conv1_lx10
684 subroutine cpu_conv1_lx9(du, u, vx, vy, vz, dx, dy, dz, &
685 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
687 integer,
parameter :: lx = 9
688 integer,
intent(in) :: nelv
689 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
690 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
691 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
692 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
693 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
694 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
695 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
696 real(kind=rp),
dimension(lx, lx, lx) :: dudr
697 real(kind=rp),
dimension(lx, lx, lx) :: duds
698 real(kind=rp),
dimension(lx, lx, lx) :: dudt
699 integer :: e, i, j, k
704 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
705 + dx(i,2) * u(2,j,1,e) &
706 + dx(i,3) * u(3,j,1,e) &
707 + dx(i,4) * u(4,j,1,e) &
708 + dx(i,5) * u(5,j,1,e) &
709 + dx(i,6) * u(6,j,1,e) &
710 + dx(i,7) * u(7,j,1,e) &
711 + dx(i,8) * u(8,j,1,e) &
712 + dx(i,9) * u(9,j,1,e)
719 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
720 + dy(j,2) * u(i,2,k,e) &
721 + dy(j,3) * u(i,3,k,e) &
722 + dy(j,4) * u(i,4,k,e) &
723 + dy(j,5) * u(i,5,k,e) &
724 + dy(j,6) * u(i,6,k,e) &
725 + dy(j,7) * u(i,7,k,e) &
726 + dy(j,8) * u(i,8,k,e) &
727 + dy(j,9) * u(i,9,k,e)
734 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
735 + dz(k,2) * u(i,1,2,e) &
736 + dz(k,3) * u(i,1,3,e) &
737 + dz(k,4) * u(i,1,4,e) &
738 + dz(k,5) * u(i,1,5,e) &
739 + dz(k,6) * u(i,1,6,e) &
740 + dz(k,7) * u(i,1,7,e) &
741 + dz(k,8) * u(i,1,8,e) &
742 + dz(k,9) * u(i,1,9,e)
746 do i = 1, lx * lx * lx
747 du(i,1,1,e) = jacinv(i,1,1,e) &
749 * ( drdx(i,1,1,e) * dudr(i,1,1) &
750 + dsdx(i,1,1,e) * duds(i,1,1) &
751 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
753 * ( drdy(i,1,1,e) * dudr(i,1,1) &
754 + dsdy(i,1,1,e) * duds(i,1,1) &
755 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
757 * ( drdz(i,1,1,e) * dudr(i,1,1) &
758 + dsdz(i,1,1,e) * duds(i,1,1) &
759 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
763 end subroutine cpu_conv1_lx9
765 subroutine cpu_conv1_lx8(du, u, vx, vy, vz, dx, dy, dz, &
766 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
768 integer,
parameter :: lx = 8
769 integer,
intent(in) :: nelv
770 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
771 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
772 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
773 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
774 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
775 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
776 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
777 real(kind=rp),
dimension(lx, lx, lx) :: dudr
778 real(kind=rp),
dimension(lx, lx, lx) :: duds
779 real(kind=rp),
dimension(lx, lx, lx) :: dudt
780 integer :: e, i, j, k
785 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
786 + dx(i,2) * u(2,j,1,e) &
787 + dx(i,3) * u(3,j,1,e) &
788 + dx(i,4) * u(4,j,1,e) &
789 + dx(i,5) * u(5,j,1,e) &
790 + dx(i,6) * u(6,j,1,e) &
791 + dx(i,7) * u(7,j,1,e) &
792 + dx(i,8) * u(8,j,1,e)
799 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
800 + dy(j,2) * u(i,2,k,e) &
801 + dy(j,3) * u(i,3,k,e) &
802 + dy(j,4) * u(i,4,k,e) &
803 + dy(j,5) * u(i,5,k,e) &
804 + dy(j,6) * u(i,6,k,e) &
805 + dy(j,7) * u(i,7,k,e) &
806 + dy(j,8) * u(i,8,k,e)
813 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
814 + dz(k,2) * u(i,1,2,e) &
815 + dz(k,3) * u(i,1,3,e) &
816 + dz(k,4) * u(i,1,4,e) &
817 + dz(k,5) * u(i,1,5,e) &
818 + dz(k,6) * u(i,1,6,e) &
819 + dz(k,7) * u(i,1,7,e) &
820 + dz(k,8) * u(i,1,8,e)
824 do i = 1, lx * lx * lx
825 du(i,1,1,e) = jacinv(i,1,1,e) &
827 * ( drdx(i,1,1,e) * dudr(i,1,1) &
828 + dsdx(i,1,1,e) * duds(i,1,1) &
829 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
831 * ( drdy(i,1,1,e) * dudr(i,1,1) &
832 + dsdy(i,1,1,e) * duds(i,1,1) &
833 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
835 * ( drdz(i,1,1,e) * dudr(i,1,1) &
836 + dsdz(i,1,1,e) * duds(i,1,1) &
837 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
841 end subroutine cpu_conv1_lx8
843 subroutine cpu_conv1_lx7(du, u, vx, vy, vz, dx, dy, dz, &
844 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
846 integer,
parameter :: lx = 7
847 integer,
intent(in) :: nelv
848 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
849 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
850 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
851 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
852 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
853 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
854 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
855 real(kind=rp),
dimension(lx, lx, lx) :: dudr
856 real(kind=rp),
dimension(lx, lx, lx) :: duds
857 real(kind=rp),
dimension(lx, lx, lx) :: dudt
858 integer :: e, i, j, k
863 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
864 + dx(i,2) * u(2,j,1,e) &
865 + dx(i,3) * u(3,j,1,e) &
866 + dx(i,4) * u(4,j,1,e) &
867 + dx(i,5) * u(5,j,1,e) &
868 + dx(i,6) * u(6,j,1,e) &
869 + dx(i,7) * u(7,j,1,e)
876 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
877 + dy(j,2) * u(i,2,k,e) &
878 + dy(j,3) * u(i,3,k,e) &
879 + dy(j,4) * u(i,4,k,e) &
880 + dy(j,5) * u(i,5,k,e) &
881 + dy(j,6) * u(i,6,k,e) &
882 + dy(j,7) * u(i,7,k,e)
889 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
890 + dz(k,2) * u(i,1,2,e) &
891 + dz(k,3) * u(i,1,3,e) &
892 + dz(k,4) * u(i,1,4,e) &
893 + dz(k,5) * u(i,1,5,e) &
894 + dz(k,6) * u(i,1,6,e) &
895 + dz(k,7) * u(i,1,7,e)
899 do i = 1, lx * lx * lx
900 du(i,1,1,e) = jacinv(i,1,1,e) &
902 * ( drdx(i,1,1,e) * dudr(i,1,1) &
903 + dsdx(i,1,1,e) * duds(i,1,1) &
904 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
906 * ( drdy(i,1,1,e) * dudr(i,1,1) &
907 + dsdy(i,1,1,e) * duds(i,1,1) &
908 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
910 * ( drdz(i,1,1,e) * dudr(i,1,1) &
911 + dsdz(i,1,1,e) * duds(i,1,1) &
912 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
916 end subroutine cpu_conv1_lx7
918 subroutine cpu_conv1_lx6(du, u, vx, vy, vz, dx, dy, dz, &
919 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
921 integer,
parameter :: lx = 6
922 integer,
intent(in) :: nelv
923 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
924 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
925 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
926 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
927 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
928 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
929 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
930 real(kind=rp),
dimension(lx, lx, lx) :: dudr
931 real(kind=rp),
dimension(lx, lx, lx) :: duds
932 real(kind=rp),
dimension(lx, lx, lx) :: dudt
933 integer :: e, i, j, k
938 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
939 + dx(i,2) * u(2,j,1,e) &
940 + dx(i,3) * u(3,j,1,e) &
941 + dx(i,4) * u(4,j,1,e) &
942 + dx(i,5) * u(5,j,1,e) &
943 + dx(i,6) * u(6,j,1,e)
950 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
951 + dy(j,2) * u(i,2,k,e) &
952 + dy(j,3) * u(i,3,k,e) &
953 + dy(j,4) * u(i,4,k,e) &
954 + dy(j,5) * u(i,5,k,e) &
955 + dy(j,6) * u(i,6,k,e)
962 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
963 + dz(k,2) * u(i,1,2,e) &
964 + dz(k,3) * u(i,1,3,e) &
965 + dz(k,4) * u(i,1,4,e) &
966 + dz(k,5) * u(i,1,5,e) &
967 + dz(k,6) * u(i,1,6,e)
971 do i = 1, lx * lx * lx
972 du(i,1,1,e) = jacinv(i,1,1,e) &
974 * ( drdx(i,1,1,e) * dudr(i,1,1) &
975 + dsdx(i,1,1,e) * duds(i,1,1) &
976 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
978 * ( drdy(i,1,1,e) * dudr(i,1,1) &
979 + dsdy(i,1,1,e) * duds(i,1,1) &
980 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
982 * ( drdz(i,1,1,e) * dudr(i,1,1) &
983 + dsdz(i,1,1,e) * duds(i,1,1) &
984 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
988 end subroutine cpu_conv1_lx6
990 subroutine cpu_conv1_lx5(du, u, vx, vy, vz, dx, dy, dz, &
991 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
993 integer,
parameter :: lx = 5
994 integer,
intent(in) :: nelv
995 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
996 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
997 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
998 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
999 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
1000 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
1001 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1002 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1003 real(kind=rp),
dimension(lx, lx, lx) :: duds
1004 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1005 integer :: e, i, j, k
1010 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
1011 + dx(i,2) * u(2,j,1,e) &
1012 + dx(i,3) * u(3,j,1,e) &
1013 + dx(i,4) * u(4,j,1,e) &
1014 + dx(i,5) * u(5,j,1,e)
1021 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
1022 + dy(j,2) * u(i,2,k,e) &
1023 + dy(j,3) * u(i,3,k,e) &
1024 + dy(j,4) * u(i,4,k,e) &
1025 + dy(j,5) * u(i,5,k,e)
1032 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
1033 + dz(k,2) * u(i,1,2,e) &
1034 + dz(k,3) * u(i,1,3,e) &
1035 + dz(k,4) * u(i,1,4,e) &
1036 + dz(k,5) * u(i,1,5,e)
1040 do i = 1, lx * lx * lx
1041 du(i,1,1,e) = jacinv(i,1,1,e) &
1043 * ( drdx(i,1,1,e) * dudr(i,1,1) &
1044 + dsdx(i,1,1,e) * duds(i,1,1) &
1045 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
1047 * ( drdy(i,1,1,e) * dudr(i,1,1) &
1048 + dsdy(i,1,1,e) * duds(i,1,1) &
1049 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
1051 * ( drdz(i,1,1,e) * dudr(i,1,1) &
1052 + dsdz(i,1,1,e) * duds(i,1,1) &
1053 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
1057 end subroutine cpu_conv1_lx5
1059 subroutine cpu_conv1_lx4(du, u, vx, vy, vz, dx, dy, dz, &
1060 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1062 integer,
parameter :: lx = 4
1063 integer,
intent(in) :: nelv
1064 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
1065 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
1066 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
1067 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
1068 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
1069 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
1070 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1071 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1072 real(kind=rp),
dimension(lx, lx, lx) :: duds
1073 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1074 integer :: e, i, j, k
1079 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
1080 + dx(i,2) * u(2,j,1,e) &
1081 + dx(i,3) * u(3,j,1,e) &
1082 + dx(i,4) * u(4,j,1,e)
1089 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
1090 + dy(j,2) * u(i,2,k,e) &
1091 + dy(j,3) * u(i,3,k,e) &
1092 + dy(j,4) * u(i,4,k,e)
1099 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
1100 + dz(k,2) * u(i,1,2,e) &
1101 + dz(k,3) * u(i,1,3,e) &
1102 + dz(k,4) * u(i,1,4,e)
1106 do i = 1, lx * lx * lx
1107 du(i,1,1,e) = jacinv(i,1,1,e) &
1109 * ( drdx(i,1,1,e) * dudr(i,1,1) &
1110 + dsdx(i,1,1,e) * duds(i,1,1) &
1111 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
1113 * ( drdy(i,1,1,e) * dudr(i,1,1) &
1114 + dsdy(i,1,1,e) * duds(i,1,1) &
1115 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
1117 * ( drdz(i,1,1,e) * dudr(i,1,1) &
1118 + dsdz(i,1,1,e) * duds(i,1,1) &
1119 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
1123 end subroutine cpu_conv1_lx4
1125 subroutine cpu_conv1_lx3(du, u, vx, vy, vz, dx, dy, dz, &
1126 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1128 integer,
parameter :: lx = 3
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) :: dudr
1138 real(kind=rp),
dimension(lx, lx, lx) :: duds
1139 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1140 integer :: e, i, j, k
1145 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
1146 + dx(i,2) * u(2,j,1,e) &
1147 + dx(i,3) * u(3,j,1,e)
1154 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
1155 + dy(j,2) * u(i,2,k,e) &
1156 + dy(j,3) * u(i,3,k,e)
1163 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
1164 + dz(k,2) * u(i,1,2,e) &
1165 + dz(k,3) * u(i,1,3,e)
1169 do i = 1, lx * lx * lx
1170 du(i,1,1,e) = jacinv(i,1,1,e) &
1172 * ( drdx(i,1,1,e) * dudr(i,1,1) &
1173 + dsdx(i,1,1,e) * duds(i,1,1) &
1174 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
1176 * ( drdy(i,1,1,e) * dudr(i,1,1) &
1177 + dsdy(i,1,1,e) * duds(i,1,1) &
1178 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
1180 * ( drdz(i,1,1,e) * dudr(i,1,1) &
1181 + dsdz(i,1,1,e) * duds(i,1,1) &
1182 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
1186 end subroutine cpu_conv1_lx3
1188 subroutine cpu_conv1_lx2(du, u, vx, vy, vz, dx, dy, dz, &
1189 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1191 integer,
parameter :: lx = 2
1192 integer,
intent(in) :: nelv
1193 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(inout) :: du
1194 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: u, vx, vy, vz
1195 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdx, dsdx, dtdx
1196 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdy, dsdy, dtdy
1197 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: drdz, dsdz, dtdz
1198 real(kind=rp),
dimension(lx, lx, lx, nelv),
intent(in) :: jacinv
1199 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1200 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1201 real(kind=rp),
dimension(lx, lx, lx) :: duds
1202 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1203 integer :: e, i, j, k
1208 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
1209 + dx(i,2) * u(2,j,1,e)
1216 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
1217 + dy(j,2) * u(i,2,k,e)
1224 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
1225 + dz(k,2) * u(i,1,2,e)
1229 do i = 1, lx * lx * lx
1230 du(i,1,1,e) = jacinv(i,1,1,e) &
1232 * ( drdx(i,1,1,e) * dudr(i,1,1) &
1233 + dsdx(i,1,1,e) * duds(i,1,1) &
1234 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
1236 * ( drdy(i,1,1,e) * dudr(i,1,1) &
1237 + dsdy(i,1,1,e) * duds(i,1,1) &
1238 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
1240 * ( drdz(i,1,1,e) * dudr(i,1,1) &
1241 + dsdz(i,1,1,e) * duds(i,1,1) &
1242 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
1246 end subroutine cpu_conv1_lx2
1248 subroutine opr_cpu_conv1_single(du, u, vx, vy, vz, Xh, coef, e)
1249 integer,
parameter :: e_len = 1
1250 type(space_t),
intent(in) :: Xh
1251 type(coef_t),
intent(in) :: coef
1252 integer,
intent(in) :: e
1253 real(kind=rp),
intent(inout) :: du(xh%lxyz, e_len)
1254 real(kind=rp),
intent(inout) :: u(xh%lx, xh%ly, xh%lz, e_len)
1255 real(kind=rp),
intent(inout) :: vx(xh%lx, xh%ly, xh%lz, e_len)
1256 real(kind=rp),
intent(inout) :: vy(xh%lx, xh%ly, xh%lz, e_len)
1257 real(kind=rp),
intent(inout) :: vz(xh%lx, xh%ly, xh%lz, e_len)
1259 associate(drdx => coef%drdx, drdy => coef%drdy, drdz => coef%drdz, &
1260 dsdx => coef%dsdx, dsdy => coef%dsdy, dsdz => coef%dsdz, &
1261 dtdx => coef%dtdx, dtdy => coef%dtdy, dtdz => coef%dtdz, &
1262 jacinv => coef%jacinv)
1265 call cpu_conv1_lx14_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1266 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1267 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1268 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1271 call cpu_conv1_lx13_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1272 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1273 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1274 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1277 call cpu_conv1_lx12_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1278 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1279 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1280 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1283 call cpu_conv1_lx11_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1284 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1285 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1286 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1289 call cpu_conv1_lx10_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1290 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1291 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1292 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1295 call cpu_conv1_lx9_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1296 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1297 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1298 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1301 call cpu_conv1_lx8_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1302 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1303 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1304 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1307 call cpu_conv1_lx7_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1308 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1309 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1310 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1313 call cpu_conv1_lx6_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1314 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1315 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1316 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1319 call cpu_conv1_lx5_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1320 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1321 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1322 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1325 call cpu_conv1_lx4_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1326 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1327 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1328 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1331 call cpu_conv1_lx3_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1332 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1333 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1334 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1337 call cpu_conv1_lx2_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1338 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1339 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1340 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1343 call cpu_conv1_lx_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1344 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1345 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1346 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1347 jacinv(1,1,1,e), xh%lx)
1351 end subroutine opr_cpu_conv1_single
1354 subroutine cpu_conv1_lx_single(du, u, vx, vy, vz, dx, dy, dz, &
1355 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1357 integer,
intent(in) :: lx
1358 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
1359 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
1360 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
1361 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
1362 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
1363 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
1364 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1365 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1366 real(kind=rp),
dimension(lx, lx, lx) :: duds
1367 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1368 real(kind=rp) :: tmp
1369 integer :: i, j, k, l
1375 tmp = tmp + dx(i,k) * u(k,j,1)
1386 tmp = tmp + dy(j,l) * u(i,l,k)
1397 tmp = tmp + dz(k,l) * u(i,1,l)
1403 do i = 1, lx * lx * lx
1404 du(i,1,1) = jacinv(i,1,1) &
1406 * ( drdx(i,1,1) * dudr(i,1,1) &
1407 + dsdx(i,1,1) * duds(i,1,1) &
1408 + dtdx(i,1,1) * dudt(i,1,1) ) &
1410 * ( drdy(i,1,1) * dudr(i,1,1) &
1411 + dsdy(i,1,1) * duds(i,1,1) &
1412 + dtdy(i,1,1) * dudt(i,1,1) ) &
1414 * ( drdz(i,1,1) * dudr(i,1,1) &
1415 + dsdz(i,1,1) * duds(i,1,1) &
1416 + dtdz(i,1,1) * dudt(i,1,1) ) )
1419 end subroutine cpu_conv1_lx_single
1422 subroutine cpu_conv1_lx14_single(du, u, vx, vy, vz, dx, dy, dz, &
1423 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1424 integer,
parameter :: lx = 14
1425 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
1426 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
1427 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
1428 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
1429 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
1430 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
1431 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1432 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1433 real(kind=rp),
dimension(lx, lx, lx) :: duds
1434 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1439 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1440 + dx(i,2) * u(2,j,1) &
1441 + dx(i,3) * u(3,j,1) &
1442 + dx(i,4) * u(4,j,1) &
1443 + dx(i,5) * u(5,j,1) &
1444 + dx(i,6) * u(6,j,1) &
1445 + dx(i,7) * u(7,j,1) &
1446 + dx(i,8) * u(8,j,1) &
1447 + dx(i,9) * u(9,j,1) &
1448 + dx(i,10) * u(10,j,1) &
1449 + dx(i,11) * u(11,j,1) &
1450 + dx(i,12) * u(12,j,1) &
1451 + dx(i,13) * u(13,j,1) &
1452 + dx(i,14) * u(14,j,1)
1459 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1460 + dy(j,2) * u(i,2,k) &
1461 + dy(j,3) * u(i,3,k) &
1462 + dy(j,4) * u(i,4,k) &
1463 + dy(j,5) * u(i,5,k) &
1464 + dy(j,6) * u(i,6,k) &
1465 + dy(j,7) * u(i,7,k) &
1466 + dy(j,8) * u(i,8,k) &
1467 + dy(j,9) * u(i,9,k) &
1468 + dy(j,10) * u(i,10,k) &
1469 + dy(j,11) * u(i,11,k) &
1470 + dy(j,12) * u(i,12,k) &
1471 + dy(j,13) * u(i,13,k) &
1472 + dy(j,14) * u(i,14,k)
1479 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1480 + dz(k,2) * u(i,1,2) &
1481 + dz(k,3) * u(i,1,3) &
1482 + dz(k,4) * u(i,1,4) &
1483 + dz(k,5) * u(i,1,5) &
1484 + dz(k,6) * u(i,1,6) &
1485 + dz(k,7) * u(i,1,7) &
1486 + dz(k,8) * u(i,1,8) &
1487 + dz(k,9) * u(i,1,9) &
1488 + dz(k,10) * u(i,1,10) &
1489 + dz(k,11) * u(i,1,11) &
1490 + dz(k,12) * u(i,1,12) &
1491 + dz(k,13) * u(i,1,13) &
1492 + dz(k,14) * u(i,1,14)
1496 do i = 1, lx * lx * lx
1497 du(i,1,1) = jacinv(i,1,1) &
1499 * ( drdx(i,1,1) * dudr(i,1,1) &
1500 + dsdx(i,1,1) * duds(i,1,1) &
1501 + dtdx(i,1,1) * dudt(i,1,1) ) &
1503 * ( drdy(i,1,1) * dudr(i,1,1) &
1504 + dsdy(i,1,1) * duds(i,1,1) &
1505 + dtdy(i,1,1) * dudt(i,1,1) ) &
1507 * ( drdz(i,1,1) * dudr(i,1,1) &
1508 + dsdz(i,1,1) * duds(i,1,1) &
1509 + dtdz(i,1,1) * dudt(i,1,1) ) )
1512 end subroutine cpu_conv1_lx14_single
1515 subroutine cpu_conv1_lx13_single(du, u, vx, vy, vz, dx, dy, dz, &
1516 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1517 integer,
parameter :: lx = 13
1518 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
1519 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
1520 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
1521 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
1522 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
1523 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
1524 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1525 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1526 real(kind=rp),
dimension(lx, lx, lx) :: duds
1527 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1532 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1533 + dx(i,2) * u(2,j,1) &
1534 + dx(i,3) * u(3,j,1) &
1535 + dx(i,4) * u(4,j,1) &
1536 + dx(i,5) * u(5,j,1) &
1537 + dx(i,6) * u(6,j,1) &
1538 + dx(i,7) * u(7,j,1) &
1539 + dx(i,8) * u(8,j,1) &
1540 + dx(i,9) * u(9,j,1) &
1541 + dx(i,10) * u(10,j,1) &
1542 + dx(i,11) * u(11,j,1) &
1543 + dx(i,12) * u(12,j,1) &
1544 + dx(i,13) * u(13,j,1)
1551 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1552 + dy(j,2) * u(i,2,k) &
1553 + dy(j,3) * u(i,3,k) &
1554 + dy(j,4) * u(i,4,k) &
1555 + dy(j,5) * u(i,5,k) &
1556 + dy(j,6) * u(i,6,k) &
1557 + dy(j,7) * u(i,7,k) &
1558 + dy(j,8) * u(i,8,k) &
1559 + dy(j,9) * u(i,9,k) &
1560 + dy(j,10) * u(i,10,k) &
1561 + dy(j,11) * u(i,11,k) &
1562 + dy(j,12) * u(i,12,k) &
1563 + dy(j,13) * u(i,13,k)
1570 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1571 + dz(k,2) * u(i,1,2) &
1572 + dz(k,3) * u(i,1,3) &
1573 + dz(k,4) * u(i,1,4) &
1574 + dz(k,5) * u(i,1,5) &
1575 + dz(k,6) * u(i,1,6) &
1576 + dz(k,7) * u(i,1,7) &
1577 + dz(k,8) * u(i,1,8) &
1578 + dz(k,9) * u(i,1,9) &
1579 + dz(k,10) * u(i,1,10) &
1580 + dz(k,11) * u(i,1,11) &
1581 + dz(k,12) * u(i,1,12) &
1582 + dz(k,13) * u(i,1,13)
1586 do i = 1, lx * lx * lx
1587 du(i,1,1) = jacinv(i,1,1) &
1589 * ( drdx(i,1,1) * dudr(i,1,1) &
1590 + dsdx(i,1,1) * duds(i,1,1) &
1591 + dtdx(i,1,1) * dudt(i,1,1) ) &
1593 * ( drdy(i,1,1) * dudr(i,1,1) &
1594 + dsdy(i,1,1) * duds(i,1,1) &
1595 + dtdy(i,1,1) * dudt(i,1,1) ) &
1597 * ( drdz(i,1,1) * dudr(i,1,1) &
1598 + dsdz(i,1,1) * duds(i,1,1) &
1599 + dtdz(i,1,1) * dudt(i,1,1) ) )
1602 end subroutine cpu_conv1_lx13_single
1605 subroutine cpu_conv1_lx12_single(du, u, vx, vy, vz, dx, dy, dz, &
1606 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1607 integer,
parameter :: lx = 12
1608 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
1609 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
1610 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
1611 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
1612 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
1613 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
1614 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1615 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1616 real(kind=rp),
dimension(lx, lx, lx) :: duds
1617 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1622 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1623 + dx(i,2) * u(2,j,1) &
1624 + dx(i,3) * u(3,j,1) &
1625 + dx(i,4) * u(4,j,1) &
1626 + dx(i,5) * u(5,j,1) &
1627 + dx(i,6) * u(6,j,1) &
1628 + dx(i,7) * u(7,j,1) &
1629 + dx(i,8) * u(8,j,1) &
1630 + dx(i,9) * u(9,j,1) &
1631 + dx(i,10) * u(10,j,1) &
1632 + dx(i,11) * u(11,j,1) &
1633 + dx(i,12) * u(12,j,1)
1640 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1641 + dy(j,2) * u(i,2,k) &
1642 + dy(j,3) * u(i,3,k) &
1643 + dy(j,4) * u(i,4,k) &
1644 + dy(j,5) * u(i,5,k) &
1645 + dy(j,6) * u(i,6,k) &
1646 + dy(j,7) * u(i,7,k) &
1647 + dy(j,8) * u(i,8,k) &
1648 + dy(j,9) * u(i,9,k) &
1649 + dy(j,10) * u(i,10,k) &
1650 + dy(j,11) * u(i,11,k) &
1651 + dy(j,12) * u(i,12,k)
1658 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1659 + dz(k,2) * u(i,1,2) &
1660 + dz(k,3) * u(i,1,3) &
1661 + dz(k,4) * u(i,1,4) &
1662 + dz(k,5) * u(i,1,5) &
1663 + dz(k,6) * u(i,1,6) &
1664 + dz(k,7) * u(i,1,7) &
1665 + dz(k,8) * u(i,1,8) &
1666 + dz(k,9) * u(i,1,9) &
1667 + dz(k,10) * u(i,1,10) &
1668 + dz(k,11) * u(i,1,11) &
1669 + dz(k,12) * u(i,1,12)
1673 do i = 1, lx * lx * lx
1674 du(i,1,1) = jacinv(i,1,1) &
1676 * ( drdx(i,1,1) * dudr(i,1,1) &
1677 + dsdx(i,1,1) * duds(i,1,1) &
1678 + dtdx(i,1,1) * dudt(i,1,1) ) &
1680 * ( drdy(i,1,1) * dudr(i,1,1) &
1681 + dsdy(i,1,1) * duds(i,1,1) &
1682 + dtdy(i,1,1) * dudt(i,1,1) ) &
1684 * ( drdz(i,1,1) * dudr(i,1,1) &
1685 + dsdz(i,1,1) * duds(i,1,1) &
1686 + dtdz(i,1,1) * dudt(i,1,1) ) )
1689 end subroutine cpu_conv1_lx12_single
1692 subroutine cpu_conv1_lx11_single(du, u, vx, vy, vz, dx, dy, dz, &
1693 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1694 integer,
parameter :: lx = 11
1695 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
1696 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
1697 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
1698 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
1699 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
1700 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
1701 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1702 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1703 real(kind=rp),
dimension(lx, lx, lx) :: duds
1704 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1709 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1710 + dx(i,2) * u(2,j,1) &
1711 + dx(i,3) * u(3,j,1) &
1712 + dx(i,4) * u(4,j,1) &
1713 + dx(i,5) * u(5,j,1) &
1714 + dx(i,6) * u(6,j,1) &
1715 + dx(i,7) * u(7,j,1) &
1716 + dx(i,8) * u(8,j,1) &
1717 + dx(i,9) * u(9,j,1) &
1718 + dx(i,10) * u(10,j,1) &
1719 + dx(i,11) * u(11,j,1)
1726 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1727 + dy(j,2) * u(i,2,k) &
1728 + dy(j,3) * u(i,3,k) &
1729 + dy(j,4) * u(i,4,k) &
1730 + dy(j,5) * u(i,5,k) &
1731 + dy(j,6) * u(i,6,k) &
1732 + dy(j,7) * u(i,7,k) &
1733 + dy(j,8) * u(i,8,k) &
1734 + dy(j,9) * u(i,9,k) &
1735 + dy(j,10) * u(i,10,k) &
1736 + dy(j,11) * u(i,11,k)
1743 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1744 + dz(k,2) * u(i,1,2) &
1745 + dz(k,3) * u(i,1,3) &
1746 + dz(k,4) * u(i,1,4) &
1747 + dz(k,5) * u(i,1,5) &
1748 + dz(k,6) * u(i,1,6) &
1749 + dz(k,7) * u(i,1,7) &
1750 + dz(k,8) * u(i,1,8) &
1751 + dz(k,9) * u(i,1,9) &
1752 + dz(k,10) * u(i,1,10) &
1753 + dz(k,11) * u(i,1,11)
1757 do i = 1, lx * lx * lx
1758 du(i,1,1) = jacinv(i,1,1) &
1760 * ( drdx(i,1,1) * dudr(i,1,1) &
1761 + dsdx(i,1,1) * duds(i,1,1) &
1762 + dtdx(i,1,1) * dudt(i,1,1) ) &
1764 * ( drdy(i,1,1) * dudr(i,1,1) &
1765 + dsdy(i,1,1) * duds(i,1,1) &
1766 + dtdy(i,1,1) * dudt(i,1,1) ) &
1768 * ( drdz(i,1,1) * dudr(i,1,1) &
1769 + dsdz(i,1,1) * duds(i,1,1) &
1770 + dtdz(i,1,1) * dudt(i,1,1) ) )
1773 end subroutine cpu_conv1_lx11_single
1776 subroutine cpu_conv1_lx10_single(du, u, vx, vy, vz, dx, dy, dz, &
1777 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1778 integer,
parameter :: lx = 10
1779 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
1780 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
1781 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
1782 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
1783 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
1784 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
1785 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1786 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1787 real(kind=rp),
dimension(lx, lx, lx) :: duds
1788 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1793 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1794 + dx(i,2) * u(2,j,1) &
1795 + dx(i,3) * u(3,j,1) &
1796 + dx(i,4) * u(4,j,1) &
1797 + dx(i,5) * u(5,j,1) &
1798 + dx(i,6) * u(6,j,1) &
1799 + dx(i,7) * u(7,j,1) &
1800 + dx(i,8) * u(8,j,1) &
1801 + dx(i,9) * u(9,j,1) &
1802 + dx(i,10) * u(10,j,1)
1809 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1810 + dy(j,2) * u(i,2,k) &
1811 + dy(j,3) * u(i,3,k) &
1812 + dy(j,4) * u(i,4,k) &
1813 + dy(j,5) * u(i,5,k) &
1814 + dy(j,6) * u(i,6,k) &
1815 + dy(j,7) * u(i,7,k) &
1816 + dy(j,8) * u(i,8,k) &
1817 + dy(j,9) * u(i,9,k) &
1818 + dy(j,10) * u(i,10,k)
1825 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1826 + dz(k,2) * u(i,1,2) &
1827 + dz(k,3) * u(i,1,3) &
1828 + dz(k,4) * u(i,1,4) &
1829 + dz(k,5) * u(i,1,5) &
1830 + dz(k,6) * u(i,1,6) &
1831 + dz(k,7) * u(i,1,7) &
1832 + dz(k,8) * u(i,1,8) &
1833 + dz(k,9) * u(i,1,9) &
1834 + dz(k,10) * u(i,1,10)
1838 do i = 1, lx * lx * lx
1839 du(i,1,1) = jacinv(i,1,1) &
1841 * ( drdx(i,1,1) * dudr(i,1,1) &
1842 + dsdx(i,1,1) * duds(i,1,1) &
1843 + dtdx(i,1,1) * dudt(i,1,1) ) &
1845 * ( drdy(i,1,1) * dudr(i,1,1) &
1846 + dsdy(i,1,1) * duds(i,1,1) &
1847 + dtdy(i,1,1) * dudt(i,1,1) ) &
1849 * ( drdz(i,1,1) * dudr(i,1,1) &
1850 + dsdz(i,1,1) * duds(i,1,1) &
1851 + dtdz(i,1,1) * dudt(i,1,1) ) )
1854 end subroutine cpu_conv1_lx10_single
1857 subroutine cpu_conv1_lx9_single(du, u, vx, vy, vz, dx, dy, dz, &
1858 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1859 integer,
parameter :: lx = 9
1860 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
1861 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
1862 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
1863 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
1864 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
1865 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
1866 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1867 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1868 real(kind=rp),
dimension(lx, lx, lx) :: duds
1869 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1874 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1875 + dx(i,2) * u(2,j,1) &
1876 + dx(i,3) * u(3,j,1) &
1877 + dx(i,4) * u(4,j,1) &
1878 + dx(i,5) * u(5,j,1) &
1879 + dx(i,6) * u(6,j,1) &
1880 + dx(i,7) * u(7,j,1) &
1881 + dx(i,8) * u(8,j,1) &
1882 + dx(i,9) * u(9,j,1)
1889 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1890 + dy(j,2) * u(i,2,k) &
1891 + dy(j,3) * u(i,3,k) &
1892 + dy(j,4) * u(i,4,k) &
1893 + dy(j,5) * u(i,5,k) &
1894 + dy(j,6) * u(i,6,k) &
1895 + dy(j,7) * u(i,7,k) &
1896 + dy(j,8) * u(i,8,k) &
1897 + dy(j,9) * u(i,9,k)
1904 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1905 + dz(k,2) * u(i,1,2) &
1906 + dz(k,3) * u(i,1,3) &
1907 + dz(k,4) * u(i,1,4) &
1908 + dz(k,5) * u(i,1,5) &
1909 + dz(k,6) * u(i,1,6) &
1910 + dz(k,7) * u(i,1,7) &
1911 + dz(k,8) * u(i,1,8) &
1912 + dz(k,9) * u(i,1,9)
1916 do i = 1, lx * lx * lx
1917 du(i,1,1) = jacinv(i,1,1) &
1919 * ( drdx(i,1,1) * dudr(i,1,1) &
1920 + dsdx(i,1,1) * duds(i,1,1) &
1921 + dtdx(i,1,1) * dudt(i,1,1) ) &
1923 * ( drdy(i,1,1) * dudr(i,1,1) &
1924 + dsdy(i,1,1) * duds(i,1,1) &
1925 + dtdy(i,1,1) * dudt(i,1,1) ) &
1927 * ( drdz(i,1,1) * dudr(i,1,1) &
1928 + dsdz(i,1,1) * duds(i,1,1) &
1929 + dtdz(i,1,1) * dudt(i,1,1) ) )
1932 end subroutine cpu_conv1_lx9_single
1935 subroutine cpu_conv1_lx8_single(du, u, vx, vy, vz, dx, dy, dz, &
1936 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1937 integer,
parameter :: lx = 8
1938 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
1939 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
1940 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
1941 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
1942 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
1943 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
1944 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1945 real(kind=rp),
dimension(lx, lx, lx) :: dudr
1946 real(kind=rp),
dimension(lx, lx, lx) :: duds
1947 real(kind=rp),
dimension(lx, lx, lx) :: dudt
1952 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1953 + dx(i,2) * u(2,j,1) &
1954 + dx(i,3) * u(3,j,1) &
1955 + dx(i,4) * u(4,j,1) &
1956 + dx(i,5) * u(5,j,1) &
1957 + dx(i,6) * u(6,j,1) &
1958 + dx(i,7) * u(7,j,1) &
1959 + dx(i,8) * u(8,j,1)
1966 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1967 + dy(j,2) * u(i,2,k) &
1968 + dy(j,3) * u(i,3,k) &
1969 + dy(j,4) * u(i,4,k) &
1970 + dy(j,5) * u(i,5,k) &
1971 + dy(j,6) * u(i,6,k) &
1972 + dy(j,7) * u(i,7,k) &
1973 + dy(j,8) * u(i,8,k)
1980 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1981 + dz(k,2) * u(i,1,2) &
1982 + dz(k,3) * u(i,1,3) &
1983 + dz(k,4) * u(i,1,4) &
1984 + dz(k,5) * u(i,1,5) &
1985 + dz(k,6) * u(i,1,6) &
1986 + dz(k,7) * u(i,1,7) &
1987 + dz(k,8) * u(i,1,8)
1991 do i = 1, lx * lx * lx
1992 du(i,1,1) = jacinv(i,1,1) &
1994 * ( drdx(i,1,1) * dudr(i,1,1) &
1995 + dsdx(i,1,1) * duds(i,1,1) &
1996 + dtdx(i,1,1) * dudt(i,1,1) ) &
1998 * ( drdy(i,1,1) * dudr(i,1,1) &
1999 + dsdy(i,1,1) * duds(i,1,1) &
2000 + dtdy(i,1,1) * dudt(i,1,1) ) &
2002 * ( drdz(i,1,1) * dudr(i,1,1) &
2003 + dsdz(i,1,1) * duds(i,1,1) &
2004 + dtdz(i,1,1) * dudt(i,1,1) ) )
2007 end subroutine cpu_conv1_lx8_single
2010 subroutine cpu_conv1_lx7_single(du, u, vx, vy, vz, dx, dy, dz, &
2011 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2012 integer,
parameter :: lx = 7
2013 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
2014 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
2015 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
2016 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
2017 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
2018 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
2019 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
2020 real(kind=rp),
dimension(lx, lx, lx) :: dudr
2021 real(kind=rp),
dimension(lx, lx, lx) :: duds
2022 real(kind=rp),
dimension(lx, lx, lx) :: dudt
2027 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2028 + dx(i,2) * u(2,j,1) &
2029 + dx(i,3) * u(3,j,1) &
2030 + dx(i,4) * u(4,j,1) &
2031 + dx(i,5) * u(5,j,1) &
2032 + dx(i,6) * u(6,j,1) &
2033 + dx(i,7) * u(7,j,1)
2040 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2041 + dy(j,2) * u(i,2,k) &
2042 + dy(j,3) * u(i,3,k) &
2043 + dy(j,4) * u(i,4,k) &
2044 + dy(j,5) * u(i,5,k) &
2045 + dy(j,6) * u(i,6,k) &
2046 + dy(j,7) * u(i,7,k)
2053 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2054 + dz(k,2) * u(i,1,2) &
2055 + dz(k,3) * u(i,1,3) &
2056 + dz(k,4) * u(i,1,4) &
2057 + dz(k,5) * u(i,1,5) &
2058 + dz(k,6) * u(i,1,6) &
2059 + dz(k,7) * u(i,1,7)
2063 do i = 1, lx * lx * lx
2064 du(i,1,1) = jacinv(i,1,1) &
2066 * ( drdx(i,1,1) * dudr(i,1,1) &
2067 + dsdx(i,1,1) * duds(i,1,1) &
2068 + dtdx(i,1,1) * dudt(i,1,1) ) &
2070 * ( drdy(i,1,1) * dudr(i,1,1) &
2071 + dsdy(i,1,1) * duds(i,1,1) &
2072 + dtdy(i,1,1) * dudt(i,1,1) ) &
2074 * ( drdz(i,1,1) * dudr(i,1,1) &
2075 + dsdz(i,1,1) * duds(i,1,1) &
2076 + dtdz(i,1,1) * dudt(i,1,1) ) )
2079 end subroutine cpu_conv1_lx7_single
2082 subroutine cpu_conv1_lx6_single(du, u, vx, vy, vz, dx, dy, dz, &
2083 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2084 integer,
parameter :: lx = 6
2085 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
2086 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
2087 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
2088 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
2089 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
2090 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
2091 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
2092 real(kind=rp),
dimension(lx, lx, lx) :: dudr
2093 real(kind=rp),
dimension(lx, lx, lx) :: duds
2094 real(kind=rp),
dimension(lx, lx, lx) :: dudt
2099 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2100 + dx(i,2) * u(2,j,1) &
2101 + dx(i,3) * u(3,j,1) &
2102 + dx(i,4) * u(4,j,1) &
2103 + dx(i,5) * u(5,j,1) &
2104 + dx(i,6) * u(6,j,1)
2111 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2112 + dy(j,2) * u(i,2,k) &
2113 + dy(j,3) * u(i,3,k) &
2114 + dy(j,4) * u(i,4,k) &
2115 + dy(j,5) * u(i,5,k) &
2116 + dy(j,6) * u(i,6,k)
2123 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2124 + dz(k,2) * u(i,1,2) &
2125 + dz(k,3) * u(i,1,3) &
2126 + dz(k,4) * u(i,1,4) &
2127 + dz(k,5) * u(i,1,5) &
2128 + dz(k,6) * u(i,1,6)
2132 do i = 1, lx * lx * lx
2133 du(i,1,1) = jacinv(i,1,1) &
2135 * ( drdx(i,1,1) * dudr(i,1,1) &
2136 + dsdx(i,1,1) * duds(i,1,1) &
2137 + dtdx(i,1,1) * dudt(i,1,1) ) &
2139 * ( drdy(i,1,1) * dudr(i,1,1) &
2140 + dsdy(i,1,1) * duds(i,1,1) &
2141 + dtdy(i,1,1) * dudt(i,1,1) ) &
2143 * ( drdz(i,1,1) * dudr(i,1,1) &
2144 + dsdz(i,1,1) * duds(i,1,1) &
2145 + dtdz(i,1,1) * dudt(i,1,1) ) )
2148 end subroutine cpu_conv1_lx6_single
2151 subroutine cpu_conv1_lx5_single(du, u, vx, vy, vz, dx, dy, dz, &
2152 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2153 integer,
parameter :: lx = 5
2154 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
2155 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
2156 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
2157 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
2158 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
2159 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
2160 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
2161 real(kind=rp),
dimension(lx, lx, lx) :: dudr
2162 real(kind=rp),
dimension(lx, lx, lx) :: duds
2163 real(kind=rp),
dimension(lx, lx, lx) :: dudt
2168 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2169 + dx(i,2) * u(2,j,1) &
2170 + dx(i,3) * u(3,j,1) &
2171 + dx(i,4) * u(4,j,1) &
2172 + dx(i,5) * u(5,j,1)
2179 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2180 + dy(j,2) * u(i,2,k) &
2181 + dy(j,3) * u(i,3,k) &
2182 + dy(j,4) * u(i,4,k) &
2183 + dy(j,5) * u(i,5,k)
2190 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2191 + dz(k,2) * u(i,1,2) &
2192 + dz(k,3) * u(i,1,3) &
2193 + dz(k,4) * u(i,1,4) &
2194 + dz(k,5) * u(i,1,5)
2198 do i = 1, lx * lx * lx
2199 du(i,1,1) = jacinv(i,1,1) &
2201 * ( drdx(i,1,1) * dudr(i,1,1) &
2202 + dsdx(i,1,1) * duds(i,1,1) &
2203 + dtdx(i,1,1) * dudt(i,1,1) ) &
2205 * ( drdy(i,1,1) * dudr(i,1,1) &
2206 + dsdy(i,1,1) * duds(i,1,1) &
2207 + dtdy(i,1,1) * dudt(i,1,1) ) &
2209 * ( drdz(i,1,1) * dudr(i,1,1) &
2210 + dsdz(i,1,1) * duds(i,1,1) &
2211 + dtdz(i,1,1) * dudt(i,1,1) ) )
2214 end subroutine cpu_conv1_lx5_single
2217 subroutine cpu_conv1_lx4_single(du, u, vx, vy, vz, dx, dy, dz, &
2218 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2219 integer,
parameter :: lx = 4
2220 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
2221 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
2222 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
2223 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
2224 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
2225 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
2226 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
2227 real(kind=rp),
dimension(lx, lx, lx) :: dudr
2228 real(kind=rp),
dimension(lx, lx, lx) :: duds
2229 real(kind=rp),
dimension(lx, lx, lx) :: dudt
2234 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2235 + dx(i,2) * u(2,j,1) &
2236 + dx(i,3) * u(3,j,1) &
2237 + dx(i,4) * u(4,j,1)
2244 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2245 + dy(j,2) * u(i,2,k) &
2246 + dy(j,3) * u(i,3,k) &
2247 + dy(j,4) * u(i,4,k)
2254 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2255 + dz(k,2) * u(i,1,2) &
2256 + dz(k,3) * u(i,1,3) &
2257 + dz(k,4) * u(i,1,4)
2261 do i = 1, lx * lx * lx
2262 du(i,1,1) = jacinv(i,1,1) &
2264 * ( drdx(i,1,1) * dudr(i,1,1) &
2265 + dsdx(i,1,1) * duds(i,1,1) &
2266 + dtdx(i,1,1) * dudt(i,1,1) ) &
2268 * ( drdy(i,1,1) * dudr(i,1,1) &
2269 + dsdy(i,1,1) * duds(i,1,1) &
2270 + dtdy(i,1,1) * dudt(i,1,1) ) &
2272 * ( drdz(i,1,1) * dudr(i,1,1) &
2273 + dsdz(i,1,1) * duds(i,1,1) &
2274 + dtdz(i,1,1) * dudt(i,1,1) ) )
2277 end subroutine cpu_conv1_lx4_single
2280 subroutine cpu_conv1_lx3_single(du, u, vx, vy, vz, dx, dy, dz, &
2281 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2282 integer,
parameter :: lx = 3
2283 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
2284 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
2285 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
2286 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
2287 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
2288 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
2289 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
2290 real(kind=rp),
dimension(lx, lx, lx) :: dudr
2291 real(kind=rp),
dimension(lx, lx, lx) :: duds
2292 real(kind=rp),
dimension(lx, lx, lx) :: dudt
2297 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2298 + dx(i,2) * u(2,j,1) &
2299 + dx(i,3) * u(3,j,1)
2306 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2307 + dy(j,2) * u(i,2,k) &
2308 + dy(j,3) * u(i,3,k)
2315 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2316 + dz(k,2) * u(i,1,2) &
2317 + dz(k,3) * u(i,1,3)
2321 do i = 1, lx * lx * lx
2322 du(i,1,1) = jacinv(i,1,1) &
2324 * ( drdx(i,1,1) * dudr(i,1,1) &
2325 + dsdx(i,1,1) * duds(i,1,1) &
2326 + dtdx(i,1,1) * dudt(i,1,1) ) &
2328 * ( drdy(i,1,1) * dudr(i,1,1) &
2329 + dsdy(i,1,1) * duds(i,1,1) &
2330 + dtdy(i,1,1) * dudt(i,1,1) ) &
2332 * ( drdz(i,1,1) * dudr(i,1,1) &
2333 + dsdz(i,1,1) * duds(i,1,1) &
2334 + dtdz(i,1,1) * dudt(i,1,1) ) )
2337 end subroutine cpu_conv1_lx3_single
2340 subroutine cpu_conv1_lx2_single(du, u, vx, vy, vz, dx, dy, dz, &
2341 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2342 integer,
parameter :: lx = 2
2343 real(kind=rp),
dimension(lx, lx, lx),
intent(inout) :: du
2344 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: u, vx, vy, vz
2345 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdx, dsdx, dtdx
2346 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdy, dsdy, dtdy
2347 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: drdz, dsdz, dtdz
2348 real(kind=rp),
dimension(lx, lx, lx),
intent(in) :: jacinv
2349 real(kind=rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
2350 real(kind=rp),
dimension(lx, lx, lx) :: dudr
2351 real(kind=rp),
dimension(lx, lx, lx) :: duds
2352 real(kind=rp),
dimension(lx, lx, lx) :: dudt
2357 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2358 + dx(i,2) * u(2,j,1)
2365 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2366 + dy(j,2) * u(i,2,k)
2373 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2374 + dz(k,2) * u(i,1,2)
2378 do i = 1, lx * lx * lx
2379 du(i,1,1) = jacinv(i,1,1) &
2381 * ( drdx(i,1,1) * dudr(i,1,1) &
2382 + dsdx(i,1,1) * duds(i,1,1) &
2383 + dtdx(i,1,1) * dudt(i,1,1) ) &
2385 * ( drdy(i,1,1) * dudr(i,1,1) &
2386 + dsdy(i,1,1) * duds(i,1,1) &
2387 + dtdy(i,1,1) * dudt(i,1,1) ) &
2389 * ( drdz(i,1,1) * dudr(i,1,1) &
2390 + dsdz(i,1,1) * duds(i,1,1) &
2391 + dtdz(i,1,1) * dudt(i,1,1) ) )
2394 end subroutine cpu_conv1_lx2_single
2396end submodule cpu_conv1