Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
cpu_conv1.f90
Go to the documentation of this file.
1! Copyright (c) 2021-2024, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
34submodule(opr_cpu) cpu_conv1
35 implicit none
36
37contains
38
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)
48 integer :: e_len
49
50 e_len = e_end-e_start+1
51
52 if (e_len .eq. 1) then
53 call opr_cpu_conv1_single(du, u, vx, vy, vz, xh, coef, e_start)
54 else
55 call opr_cpu_conv1_many(du, u, vx, vy, vz, xh, coef, e_start, e_len)
56 end if
57 end subroutine opr_cpu_conv1
58
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)
68
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)
73 select case (xh%lx)
74 case (14)
75 call cpu_conv1_lx14(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
76 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
77 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
78 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
79 jacinv(1,1,1, e_start), e_len)
80 case (13)
81 call cpu_conv1_lx13(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
82 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
83 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
84 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
85 jacinv(1,1,1, e_start), e_len)
86 case (12)
87 call cpu_conv1_lx12(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
88 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
89 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
90 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
91 jacinv(1,1,1, e_start), e_len)
92 case (11)
93 call cpu_conv1_lx11(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
94 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
95 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
96 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
97 jacinv(1,1,1, e_start), e_len)
98 case (10)
99 call cpu_conv1_lx10(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
100 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
101 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
102 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
103 jacinv(1,1,1, e_start), e_len)
104 case (9)
105 call cpu_conv1_lx9(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
106 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
107 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
108 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
109 jacinv(1,1,1, e_start), e_len)
110 case (8)
111 call cpu_conv1_lx8(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
112 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
113 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
114 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
115 jacinv(1,1,1, e_start), e_len)
116 case (7)
117 call cpu_conv1_lx7(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
118 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
119 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
120 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
121 jacinv(1,1,1, e_start), e_len)
122 case (6)
123 call cpu_conv1_lx6(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
124 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
125 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
126 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
127 jacinv(1,1,1, e_start), e_len)
128 case (5)
129 call cpu_conv1_lx5(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
130 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
131 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
132 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
133 jacinv(1,1,1, e_start), e_len)
134 case (4)
135 call cpu_conv1_lx4(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
136 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
137 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
138 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
139 jacinv(1,1,1, e_start), e_len)
140 case (3)
141 call cpu_conv1_lx3(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
142 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
143 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
144 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
145 jacinv(1,1,1, e_start), e_len)
146 case (2)
147 call cpu_conv1_lx2(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
148 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
149 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
150 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
151 jacinv(1,1,1, e_start), e_len)
152 case default
153 call cpu_conv1_lx(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
154 drdx(1,1,1, e_start), dsdx(1,1,1, e_start), dtdx(1,1,1, e_start),&
155 drdy(1,1,1, e_start), dsdy(1,1,1, e_start), dtdy(1,1,1, e_start),&
156 drdz(1,1,1, e_start), dsdz(1,1,1, e_start), dtdz(1,1,1, e_start),&
157 jacinv(1,1,1, e_start), e_len, xh%lx)
158 end select
159 end associate
160
161 end subroutine opr_cpu_conv1_many
162
163 subroutine cpu_conv1_lx(du, u, vx, vy, vz, dx, dy, dz, &
164 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
165 jacinv, nelv, lx)
166 integer, intent(in) :: nelv, lx
167 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
168 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
169 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
170 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
171 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
172 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
173 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
174 real(kind=rp), dimension(lx, lx, lx) :: dudr
175 real(kind=rp), dimension(lx, lx, lx) :: duds
176 real(kind=rp), dimension(lx, lx, lx) :: dudt
177 real(kind=rp) :: tmp
178 integer :: e, i, j, k, l
179
180 do e = 1, nelv
181 do j = 1, lx * lx
182 do i = 1, lx
183 tmp = 0.0_rp
184 do k = 1, lx
185 tmp = tmp + dx(i,k) * u(k,j,1,e)
186 end do
187 dudr(i,j,1) = tmp
188 end do
189 end do
190
191 do k = 1, lx
192 do j = 1, lx
193 do i = 1, lx
194 tmp = 0.0_rp
195 do l = 1, lx
196 tmp = tmp + dy(j,l) * u(i,l,k,e)
197 end do
198 duds(i,j,k) = tmp
199 end do
200 end do
201 end do
202
203 do k = 1, lx
204 do i = 1, lx*lx
205 tmp = 0.0_rp
206 do l = 1, lx
207 tmp = tmp + dz(k,l) * u(i,1,l,e)
208 end do
209 dudt(i,1,k) = tmp
210 end do
211 end do
212
213 do i = 1, lx * lx * lx
214 du(i,1,1,e) = jacinv(i,1,1,e) &
215 * ( vx(i,1,1,e) &
216 * ( drdx(i,1,1,e) * dudr(i,1,1) &
217 + dsdx(i,1,1,e) * duds(i,1,1) &
218 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
219 + vy(i,1,1,e) &
220 * ( drdy(i,1,1,e) * dudr(i,1,1) &
221 + dsdy(i,1,1,e) * duds(i,1,1) &
222 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
223 + vz(i,1,1,e) &
224 * ( drdz(i,1,1,e) * dudr(i,1,1) &
225 + dsdz(i,1,1,e) * duds(i,1,1) &
226 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
227 end do
228 end do
229
230 end subroutine cpu_conv1_lx
231
232 subroutine cpu_conv1_lx14(du, u, vx, vy, vz, dx, dy, dz, &
233 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
234 jacinv, nelv)
235 integer, parameter :: lx = 14
236 integer, intent(in) :: nelv
237 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
238 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
239 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
240 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
241 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
242 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
243 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
244 real(kind=rp), dimension(lx, lx, lx) :: dudr
245 real(kind=rp), dimension(lx, lx, lx) :: duds
246 real(kind=rp), dimension(lx, lx, lx) :: dudt
247 integer :: e, i, j, k
248
249 do e = 1, nelv
250 do j = 1, lx * lx
251 do i = 1, lx
252 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
253 + dx(i,2) * u(2,j,1,e) &
254 + dx(i,3) * u(3,j,1,e) &
255 + dx(i,4) * u(4,j,1,e) &
256 + dx(i,5) * u(5,j,1,e) &
257 + dx(i,6) * u(6,j,1,e) &
258 + dx(i,7) * u(7,j,1,e) &
259 + dx(i,8) * u(8,j,1,e) &
260 + dx(i,9) * u(9,j,1,e) &
261 + dx(i,10) * u(10,j,1,e) &
262 + dx(i,11) * u(11,j,1,e) &
263 + dx(i,12) * u(12,j,1,e) &
264 + dx(i,13) * u(13,j,1,e) &
265 + dx(i,14) * u(14,j,1,e)
266 end do
267 end do
268
269 do k = 1, lx
270 do j = 1, lx
271 do i = 1, lx
272 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
273 + dy(j,2) * u(i,2,k,e) &
274 + dy(j,3) * u(i,3,k,e) &
275 + dy(j,4) * u(i,4,k,e) &
276 + dy(j,5) * u(i,5,k,e) &
277 + dy(j,6) * u(i,6,k,e) &
278 + dy(j,7) * u(i,7,k,e) &
279 + dy(j,8) * u(i,8,k,e) &
280 + dy(j,9) * u(i,9,k,e) &
281 + dy(j,10) * u(i,10,k,e) &
282 + dy(j,11) * u(i,11,k,e) &
283 + dy(j,12) * u(i,12,k,e) &
284 + dy(j,13) * u(i,13,k,e) &
285 + dy(j,14) * u(i,14,k,e)
286 end do
287 end do
288 end do
289
290 do k = 1, lx
291 do i = 1, lx*lx
292 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
293 + dz(k,2) * u(i,1,2,e) &
294 + dz(k,3) * u(i,1,3,e) &
295 + dz(k,4) * u(i,1,4,e) &
296 + dz(k,5) * u(i,1,5,e) &
297 + dz(k,6) * u(i,1,6,e) &
298 + dz(k,7) * u(i,1,7,e) &
299 + dz(k,8) * u(i,1,8,e) &
300 + dz(k,9) * u(i,1,9,e) &
301 + dz(k,10) * u(i,1,10,e) &
302 + dz(k,11) * u(i,1,11,e) &
303 + dz(k,12) * u(i,1,12,e) &
304 + dz(k,13) * u(i,1,13,e) &
305 + dz(k,14) * u(i,1,14,e)
306 end do
307 end do
308
309 do i = 1, lx * lx * lx
310 du(i,1,1,e) = jacinv(i,1,1,e) &
311 * ( vx(i,1,1,e) &
312 * ( drdx(i,1,1,e) * dudr(i,1,1) &
313 + dsdx(i,1,1,e) * duds(i,1,1) &
314 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
315 + vy(i,1,1,e) &
316 * ( drdy(i,1,1,e) * dudr(i,1,1) &
317 + dsdy(i,1,1,e) * duds(i,1,1) &
318 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
319 + vz(i,1,1,e) &
320 * ( drdz(i,1,1,e) * dudr(i,1,1) &
321 + dsdz(i,1,1,e) * duds(i,1,1) &
322 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
323 end do
324 end do
325
326 end subroutine cpu_conv1_lx14
327
328 subroutine cpu_conv1_lx13(du, u, vx, vy, vz, dx, dy, dz, &
329 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
330 jacinv, nelv)
331 integer, parameter :: lx = 13
332 integer, intent(in) :: nelv
333 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
334 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
335 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
336 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
337 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
338 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
339 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
340 real(kind=rp), dimension(lx, lx, lx) :: dudr
341 real(kind=rp), dimension(lx, lx, lx) :: duds
342 real(kind=rp), dimension(lx, lx, lx) :: dudt
343 integer :: e, i, j, k
344
345 do e = 1, nelv
346 do j = 1, lx * lx
347 do i = 1, lx
348 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
349 + dx(i,2) * u(2,j,1,e) &
350 + dx(i,3) * u(3,j,1,e) &
351 + dx(i,4) * u(4,j,1,e) &
352 + dx(i,5) * u(5,j,1,e) &
353 + dx(i,6) * u(6,j,1,e) &
354 + dx(i,7) * u(7,j,1,e) &
355 + dx(i,8) * u(8,j,1,e) &
356 + dx(i,9) * u(9,j,1,e) &
357 + dx(i,10) * u(10,j,1,e) &
358 + dx(i,11) * u(11,j,1,e) &
359 + dx(i,12) * u(12,j,1,e) &
360 + dx(i,13) * u(13,j,1,e)
361 end do
362 end do
363
364 do k = 1, lx
365 do j = 1, lx
366 do i = 1, lx
367 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
368 + dy(j,2) * u(i,2,k,e) &
369 + dy(j,3) * u(i,3,k,e) &
370 + dy(j,4) * u(i,4,k,e) &
371 + dy(j,5) * u(i,5,k,e) &
372 + dy(j,6) * u(i,6,k,e) &
373 + dy(j,7) * u(i,7,k,e) &
374 + dy(j,8) * u(i,8,k,e) &
375 + dy(j,9) * u(i,9,k,e) &
376 + dy(j,10) * u(i,10,k,e) &
377 + dy(j,11) * u(i,11,k,e) &
378 + dy(j,12) * u(i,12,k,e) &
379 + dy(j,13) * u(i,13,k,e)
380 end do
381 end do
382 end do
383
384 do k = 1, lx
385 do i = 1, lx*lx
386 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
387 + dz(k,2) * u(i,1,2,e) &
388 + dz(k,3) * u(i,1,3,e) &
389 + dz(k,4) * u(i,1,4,e) &
390 + dz(k,5) * u(i,1,5,e) &
391 + dz(k,6) * u(i,1,6,e) &
392 + dz(k,7) * u(i,1,7,e) &
393 + dz(k,8) * u(i,1,8,e) &
394 + dz(k,9) * u(i,1,9,e) &
395 + dz(k,10) * u(i,1,10,e) &
396 + dz(k,11) * u(i,1,11,e) &
397 + dz(k,12) * u(i,1,12,e) &
398 + dz(k,13) * u(i,1,13,e)
399 end do
400 end do
401
402 do i = 1, lx * lx * lx
403 du(i,1,1,e) = jacinv(i,1,1,e) &
404 * ( vx(i,1,1,e) &
405 * ( drdx(i,1,1,e) * dudr(i,1,1) &
406 + dsdx(i,1,1,e) * duds(i,1,1) &
407 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
408 + vy(i,1,1,e) &
409 * ( drdy(i,1,1,e) * dudr(i,1,1) &
410 + dsdy(i,1,1,e) * duds(i,1,1) &
411 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
412 + vz(i,1,1,e) &
413 * ( drdz(i,1,1,e) * dudr(i,1,1) &
414 + dsdz(i,1,1,e) * duds(i,1,1) &
415 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
416 end do
417 end do
418
419 end subroutine cpu_conv1_lx13
420
421 subroutine cpu_conv1_lx12(du, u, vx, vy, vz, dx, dy, dz, &
422 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
423 jacinv, nelv)
424 integer, parameter :: lx = 12
425 integer, intent(in) :: nelv
426 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
427 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
428 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
429 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
430 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
431 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
432 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
433 real(kind=rp), dimension(lx, lx, lx) :: dudr
434 real(kind=rp), dimension(lx, lx, lx) :: duds
435 real(kind=rp), dimension(lx, lx, lx) :: dudt
436 integer :: e, i, j, k
437
438 do e = 1, nelv
439 do j = 1, lx * lx
440 do i = 1, lx
441 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
442 + dx(i,2) * u(2,j,1,e) &
443 + dx(i,3) * u(3,j,1,e) &
444 + dx(i,4) * u(4,j,1,e) &
445 + dx(i,5) * u(5,j,1,e) &
446 + dx(i,6) * u(6,j,1,e) &
447 + dx(i,7) * u(7,j,1,e) &
448 + dx(i,8) * u(8,j,1,e) &
449 + dx(i,9) * u(9,j,1,e) &
450 + dx(i,10) * u(10,j,1,e) &
451 + dx(i,11) * u(11,j,1,e) &
452 + dx(i,12) * u(12,j,1,e)
453 end do
454 end do
455
456 do k = 1, lx
457 do j = 1, lx
458 do i = 1, lx
459 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
460 + dy(j,2) * u(i,2,k,e) &
461 + dy(j,3) * u(i,3,k,e) &
462 + dy(j,4) * u(i,4,k,e) &
463 + dy(j,5) * u(i,5,k,e) &
464 + dy(j,6) * u(i,6,k,e) &
465 + dy(j,7) * u(i,7,k,e) &
466 + dy(j,8) * u(i,8,k,e) &
467 + dy(j,9) * u(i,9,k,e) &
468 + dy(j,10) * u(i,10,k,e) &
469 + dy(j,11) * u(i,11,k,e) &
470 + dy(j,12) * u(i,12,k,e)
471 end do
472 end do
473 end do
474
475 do k = 1, lx
476 do i = 1, lx*lx
477 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
478 + dz(k,2) * u(i,1,2,e) &
479 + dz(k,3) * u(i,1,3,e) &
480 + dz(k,4) * u(i,1,4,e) &
481 + dz(k,5) * u(i,1,5,e) &
482 + dz(k,6) * u(i,1,6,e) &
483 + dz(k,7) * u(i,1,7,e) &
484 + dz(k,8) * u(i,1,8,e) &
485 + dz(k,9) * u(i,1,9,e) &
486 + dz(k,10) * u(i,1,10,e) &
487 + dz(k,11) * u(i,1,11,e) &
488 + dz(k,12) * u(i,1,12,e)
489 end do
490 end do
491
492 do i = 1, lx * lx * lx
493 du(i,1,1,e) = jacinv(i,1,1,e) &
494 * ( vx(i,1,1,e) &
495 * ( drdx(i,1,1,e) * dudr(i,1,1) &
496 + dsdx(i,1,1,e) * duds(i,1,1) &
497 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
498 + vy(i,1,1,e) &
499 * ( drdy(i,1,1,e) * dudr(i,1,1) &
500 + dsdy(i,1,1,e) * duds(i,1,1) &
501 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
502 + vz(i,1,1,e) &
503 * ( drdz(i,1,1,e) * dudr(i,1,1) &
504 + dsdz(i,1,1,e) * duds(i,1,1) &
505 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
506 end do
507 end do
508
509 end subroutine cpu_conv1_lx12
510
511 subroutine cpu_conv1_lx11(du, u, vx, vy, vz, dx, dy, dz, &
512 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
513 jacinv, nelv)
514 integer, parameter :: lx = 11
515 integer, intent(in) :: nelv
516 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
517 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
518 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
519 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
520 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
521 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
522 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
523 real(kind=rp), dimension(lx, lx, lx) :: dudr
524 real(kind=rp), dimension(lx, lx, lx) :: duds
525 real(kind=rp), dimension(lx, lx, lx) :: dudt
526 integer :: e, i, j, k
527
528 do e = 1, nelv
529 do j = 1, lx * lx
530 do i = 1, lx
531 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
532 + dx(i,2) * u(2,j,1,e) &
533 + dx(i,3) * u(3,j,1,e) &
534 + dx(i,4) * u(4,j,1,e) &
535 + dx(i,5) * u(5,j,1,e) &
536 + dx(i,6) * u(6,j,1,e) &
537 + dx(i,7) * u(7,j,1,e) &
538 + dx(i,8) * u(8,j,1,e) &
539 + dx(i,9) * u(9,j,1,e) &
540 + dx(i,10) * u(10,j,1,e) &
541 + dx(i,11) * u(11,j,1,e)
542 end do
543 end do
544
545 do k = 1, lx
546 do j = 1, lx
547 do i = 1, lx
548 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
549 + dy(j,2) * u(i,2,k,e) &
550 + dy(j,3) * u(i,3,k,e) &
551 + dy(j,4) * u(i,4,k,e) &
552 + dy(j,5) * u(i,5,k,e) &
553 + dy(j,6) * u(i,6,k,e) &
554 + dy(j,7) * u(i,7,k,e) &
555 + dy(j,8) * u(i,8,k,e) &
556 + dy(j,9) * u(i,9,k,e) &
557 + dy(j,10) * u(i,10,k,e) &
558 + dy(j,11) * u(i,11,k,e)
559 end do
560 end do
561 end do
562
563 do k = 1, lx
564 do i = 1, lx*lx
565 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
566 + dz(k,2) * u(i,1,2,e) &
567 + dz(k,3) * u(i,1,3,e) &
568 + dz(k,4) * u(i,1,4,e) &
569 + dz(k,5) * u(i,1,5,e) &
570 + dz(k,6) * u(i,1,6,e) &
571 + dz(k,7) * u(i,1,7,e) &
572 + dz(k,8) * u(i,1,8,e) &
573 + dz(k,9) * u(i,1,9,e) &
574 + dz(k,10) * u(i,1,10,e) &
575 + dz(k,11) * u(i,1,11,e)
576 end do
577 end do
578
579 do i = 1, lx * lx * lx
580 du(i,1,1,e) = jacinv(i,1,1,e) &
581 * ( vx(i,1,1,e) &
582 * ( drdx(i,1,1,e) * dudr(i,1,1) &
583 + dsdx(i,1,1,e) * duds(i,1,1) &
584 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
585 + vy(i,1,1,e) &
586 * ( drdy(i,1,1,e) * dudr(i,1,1) &
587 + dsdy(i,1,1,e) * duds(i,1,1) &
588 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
589 + vz(i,1,1,e) &
590 * ( drdz(i,1,1,e) * dudr(i,1,1) &
591 + dsdz(i,1,1,e) * duds(i,1,1) &
592 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
593 end do
594 end do
595
596 end subroutine cpu_conv1_lx11
597
598 subroutine cpu_conv1_lx10(du, u, vx, vy, vz, dx, dy, dz, &
599 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
600 jacinv, nelv)
601 integer, parameter :: lx = 10
602 integer, intent(in) :: nelv
603 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
604 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
605 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
606 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
607 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
608 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
609 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
610 real(kind=rp), dimension(lx, lx, lx) :: dudr
611 real(kind=rp), dimension(lx, lx, lx) :: duds
612 real(kind=rp), dimension(lx, lx, lx) :: dudt
613 integer :: e, i, j, k
614
615 do e = 1, nelv
616 do j = 1, lx * lx
617 do i = 1, lx
618 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
619 + dx(i,2) * u(2,j,1,e) &
620 + dx(i,3) * u(3,j,1,e) &
621 + dx(i,4) * u(4,j,1,e) &
622 + dx(i,5) * u(5,j,1,e) &
623 + dx(i,6) * u(6,j,1,e) &
624 + dx(i,7) * u(7,j,1,e) &
625 + dx(i,8) * u(8,j,1,e) &
626 + dx(i,9) * u(9,j,1,e) &
627 + dx(i,10) * u(10,j,1,e)
628 end do
629 end do
630
631 do k = 1, lx
632 do j = 1, lx
633 do i = 1, lx
634 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
635 + dy(j,2) * u(i,2,k,e) &
636 + dy(j,3) * u(i,3,k,e) &
637 + dy(j,4) * u(i,4,k,e) &
638 + dy(j,5) * u(i,5,k,e) &
639 + dy(j,6) * u(i,6,k,e) &
640 + dy(j,7) * u(i,7,k,e) &
641 + dy(j,8) * u(i,8,k,e) &
642 + dy(j,9) * u(i,9,k,e) &
643 + dy(j,10) * u(i,10,k,e)
644 end do
645 end do
646 end do
647
648 do k = 1, lx
649 do i = 1, lx*lx
650 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
651 + dz(k,2) * u(i,1,2,e) &
652 + dz(k,3) * u(i,1,3,e) &
653 + dz(k,4) * u(i,1,4,e) &
654 + dz(k,5) * u(i,1,5,e) &
655 + dz(k,6) * u(i,1,6,e) &
656 + dz(k,7) * u(i,1,7,e) &
657 + dz(k,8) * u(i,1,8,e) &
658 + dz(k,9) * u(i,1,9,e) &
659 + dz(k,10) * u(i,1,10,e)
660 end do
661 end do
662
663 do i = 1, lx * lx * lx
664 du(i,1,1,e) = jacinv(i,1,1,e) &
665 * ( vx(i,1,1,e) &
666 * ( drdx(i,1,1,e) * dudr(i,1,1) &
667 + dsdx(i,1,1,e) * duds(i,1,1) &
668 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
669 + vy(i,1,1,e) &
670 * ( drdy(i,1,1,e) * dudr(i,1,1) &
671 + dsdy(i,1,1,e) * duds(i,1,1) &
672 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
673 + vz(i,1,1,e) &
674 * ( drdz(i,1,1,e) * dudr(i,1,1) &
675 + dsdz(i,1,1,e) * duds(i,1,1) &
676 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
677 end do
678 end do
679
680 end subroutine cpu_conv1_lx10
681
682 subroutine cpu_conv1_lx9(du, u, vx, vy, vz, dx, dy, dz, &
683 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
684 jacinv, nelv)
685 integer, parameter :: lx = 9
686 integer, intent(in) :: nelv
687 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
688 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
689 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
690 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
691 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
692 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
693 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
694 real(kind=rp), dimension(lx, lx, lx) :: dudr
695 real(kind=rp), dimension(lx, lx, lx) :: duds
696 real(kind=rp), dimension(lx, lx, lx) :: dudt
697 integer :: e, i, j, k
698
699 do e = 1, nelv
700 do j = 1, lx * lx
701 do i = 1, lx
702 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
703 + dx(i,2) * u(2,j,1,e) &
704 + dx(i,3) * u(3,j,1,e) &
705 + dx(i,4) * u(4,j,1,e) &
706 + dx(i,5) * u(5,j,1,e) &
707 + dx(i,6) * u(6,j,1,e) &
708 + dx(i,7) * u(7,j,1,e) &
709 + dx(i,8) * u(8,j,1,e) &
710 + dx(i,9) * u(9,j,1,e)
711 end do
712 end do
713
714 do k = 1, lx
715 do j = 1, lx
716 do i = 1, lx
717 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
718 + dy(j,2) * u(i,2,k,e) &
719 + dy(j,3) * u(i,3,k,e) &
720 + dy(j,4) * u(i,4,k,e) &
721 + dy(j,5) * u(i,5,k,e) &
722 + dy(j,6) * u(i,6,k,e) &
723 + dy(j,7) * u(i,7,k,e) &
724 + dy(j,8) * u(i,8,k,e) &
725 + dy(j,9) * u(i,9,k,e)
726 end do
727 end do
728 end do
729
730 do k = 1, lx
731 do i = 1, lx*lx
732 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
733 + dz(k,2) * u(i,1,2,e) &
734 + dz(k,3) * u(i,1,3,e) &
735 + dz(k,4) * u(i,1,4,e) &
736 + dz(k,5) * u(i,1,5,e) &
737 + dz(k,6) * u(i,1,6,e) &
738 + dz(k,7) * u(i,1,7,e) &
739 + dz(k,8) * u(i,1,8,e) &
740 + dz(k,9) * u(i,1,9,e)
741 end do
742 end do
743
744 do i = 1, lx * lx * lx
745 du(i,1,1,e) = jacinv(i,1,1,e) &
746 * ( vx(i,1,1,e) &
747 * ( drdx(i,1,1,e) * dudr(i,1,1) &
748 + dsdx(i,1,1,e) * duds(i,1,1) &
749 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
750 + vy(i,1,1,e) &
751 * ( drdy(i,1,1,e) * dudr(i,1,1) &
752 + dsdy(i,1,1,e) * duds(i,1,1) &
753 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
754 + vz(i,1,1,e) &
755 * ( drdz(i,1,1,e) * dudr(i,1,1) &
756 + dsdz(i,1,1,e) * duds(i,1,1) &
757 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
758 end do
759 end do
760
761 end subroutine cpu_conv1_lx9
762
763 subroutine cpu_conv1_lx8(du, u, vx, vy, vz, dx, dy, dz, &
764 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
765 jacinv, nelv)
766 integer, parameter :: lx = 8
767 integer, intent(in) :: nelv
768 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
769 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
770 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
771 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
772 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
773 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
774 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
775 real(kind=rp), dimension(lx, lx, lx) :: dudr
776 real(kind=rp), dimension(lx, lx, lx) :: duds
777 real(kind=rp), dimension(lx, lx, lx) :: dudt
778 integer :: e, i, j, k
779
780 do e = 1, nelv
781 do j = 1, lx * lx
782 do i = 1, lx
783 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
784 + dx(i,2) * u(2,j,1,e) &
785 + dx(i,3) * u(3,j,1,e) &
786 + dx(i,4) * u(4,j,1,e) &
787 + dx(i,5) * u(5,j,1,e) &
788 + dx(i,6) * u(6,j,1,e) &
789 + dx(i,7) * u(7,j,1,e) &
790 + dx(i,8) * u(8,j,1,e)
791 end do
792 end do
793
794 do k = 1, lx
795 do j = 1, lx
796 do i = 1, lx
797 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
798 + dy(j,2) * u(i,2,k,e) &
799 + dy(j,3) * u(i,3,k,e) &
800 + dy(j,4) * u(i,4,k,e) &
801 + dy(j,5) * u(i,5,k,e) &
802 + dy(j,6) * u(i,6,k,e) &
803 + dy(j,7) * u(i,7,k,e) &
804 + dy(j,8) * u(i,8,k,e)
805 end do
806 end do
807 end do
808
809 do k = 1, lx
810 do i = 1, lx*lx
811 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
812 + dz(k,2) * u(i,1,2,e) &
813 + dz(k,3) * u(i,1,3,e) &
814 + dz(k,4) * u(i,1,4,e) &
815 + dz(k,5) * u(i,1,5,e) &
816 + dz(k,6) * u(i,1,6,e) &
817 + dz(k,7) * u(i,1,7,e) &
818 + dz(k,8) * u(i,1,8,e)
819 end do
820 end do
821
822 do i = 1, lx * lx * lx
823 du(i,1,1,e) = jacinv(i,1,1,e) &
824 * ( vx(i,1,1,e) &
825 * ( drdx(i,1,1,e) * dudr(i,1,1) &
826 + dsdx(i,1,1,e) * duds(i,1,1) &
827 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
828 + vy(i,1,1,e) &
829 * ( drdy(i,1,1,e) * dudr(i,1,1) &
830 + dsdy(i,1,1,e) * duds(i,1,1) &
831 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
832 + vz(i,1,1,e) &
833 * ( drdz(i,1,1,e) * dudr(i,1,1) &
834 + dsdz(i,1,1,e) * duds(i,1,1) &
835 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
836 end do
837 end do
838
839 end subroutine cpu_conv1_lx8
840
841 subroutine cpu_conv1_lx7(du, u, vx, vy, vz, dx, dy, dz, &
842 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
843 jacinv, nelv)
844 integer, parameter :: lx = 7
845 integer, intent(in) :: nelv
846 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
847 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
848 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
849 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
850 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
851 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
852 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
853 real(kind=rp), dimension(lx, lx, lx) :: dudr
854 real(kind=rp), dimension(lx, lx, lx) :: duds
855 real(kind=rp), dimension(lx, lx, lx) :: dudt
856 integer :: e, i, j, k
857
858 do e = 1, nelv
859 do j = 1, lx * lx
860 do i = 1, lx
861 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
862 + dx(i,2) * u(2,j,1,e) &
863 + dx(i,3) * u(3,j,1,e) &
864 + dx(i,4) * u(4,j,1,e) &
865 + dx(i,5) * u(5,j,1,e) &
866 + dx(i,6) * u(6,j,1,e) &
867 + dx(i,7) * u(7,j,1,e)
868 end do
869 end do
870
871 do k = 1, lx
872 do j = 1, lx
873 do i = 1, lx
874 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
875 + dy(j,2) * u(i,2,k,e) &
876 + dy(j,3) * u(i,3,k,e) &
877 + dy(j,4) * u(i,4,k,e) &
878 + dy(j,5) * u(i,5,k,e) &
879 + dy(j,6) * u(i,6,k,e) &
880 + dy(j,7) * u(i,7,k,e)
881 end do
882 end do
883 end do
884
885 do k = 1, lx
886 do i = 1, lx*lx
887 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
888 + dz(k,2) * u(i,1,2,e) &
889 + dz(k,3) * u(i,1,3,e) &
890 + dz(k,4) * u(i,1,4,e) &
891 + dz(k,5) * u(i,1,5,e) &
892 + dz(k,6) * u(i,1,6,e) &
893 + dz(k,7) * u(i,1,7,e)
894 end do
895 end do
896
897 do i = 1, lx * lx * lx
898 du(i,1,1,e) = jacinv(i,1,1,e) &
899 * ( vx(i,1,1,e) &
900 * ( drdx(i,1,1,e) * dudr(i,1,1) &
901 + dsdx(i,1,1,e) * duds(i,1,1) &
902 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
903 + vy(i,1,1,e) &
904 * ( drdy(i,1,1,e) * dudr(i,1,1) &
905 + dsdy(i,1,1,e) * duds(i,1,1) &
906 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
907 + vz(i,1,1,e) &
908 * ( drdz(i,1,1,e) * dudr(i,1,1) &
909 + dsdz(i,1,1,e) * duds(i,1,1) &
910 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
911 end do
912 end do
913
914 end subroutine cpu_conv1_lx7
915
916 subroutine cpu_conv1_lx6(du, u, vx, vy, vz, dx, dy, dz, &
917 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
918 jacinv, nelv)
919 integer, parameter :: lx = 6
920 integer, intent(in) :: nelv
921 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
922 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
923 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
924 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
925 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
926 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
927 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
928 real(kind=rp), dimension(lx, lx, lx) :: dudr
929 real(kind=rp), dimension(lx, lx, lx) :: duds
930 real(kind=rp), dimension(lx, lx, lx) :: dudt
931 integer :: e, i, j, k
932
933 do e = 1, nelv
934 do j = 1, lx * lx
935 do i = 1, lx
936 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
937 + dx(i,2) * u(2,j,1,e) &
938 + dx(i,3) * u(3,j,1,e) &
939 + dx(i,4) * u(4,j,1,e) &
940 + dx(i,5) * u(5,j,1,e) &
941 + dx(i,6) * u(6,j,1,e)
942 end do
943 end do
944
945 do k = 1, lx
946 do j = 1, lx
947 do i = 1, lx
948 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
949 + dy(j,2) * u(i,2,k,e) &
950 + dy(j,3) * u(i,3,k,e) &
951 + dy(j,4) * u(i,4,k,e) &
952 + dy(j,5) * u(i,5,k,e) &
953 + dy(j,6) * u(i,6,k,e)
954 end do
955 end do
956 end do
957
958 do k = 1, lx
959 do i = 1, lx*lx
960 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
961 + dz(k,2) * u(i,1,2,e) &
962 + dz(k,3) * u(i,1,3,e) &
963 + dz(k,4) * u(i,1,4,e) &
964 + dz(k,5) * u(i,1,5,e) &
965 + dz(k,6) * u(i,1,6,e)
966 end do
967 end do
968
969 do i = 1, lx * lx * lx
970 du(i,1,1,e) = jacinv(i,1,1,e) &
971 * ( vx(i,1,1,e) &
972 * ( drdx(i,1,1,e) * dudr(i,1,1) &
973 + dsdx(i,1,1,e) * duds(i,1,1) &
974 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
975 + vy(i,1,1,e) &
976 * ( drdy(i,1,1,e) * dudr(i,1,1) &
977 + dsdy(i,1,1,e) * duds(i,1,1) &
978 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
979 + vz(i,1,1,e) &
980 * ( drdz(i,1,1,e) * dudr(i,1,1) &
981 + dsdz(i,1,1,e) * duds(i,1,1) &
982 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
983 end do
984 end do
985
986 end subroutine cpu_conv1_lx6
987
988 subroutine cpu_conv1_lx5(du, u, vx, vy, vz, dx, dy, dz, &
989 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
990 jacinv, nelv)
991 integer, parameter :: lx = 5
992 integer, intent(in) :: nelv
993 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
994 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
995 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
996 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
997 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
998 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
999 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1000 real(kind=rp), dimension(lx, lx, lx) :: dudr
1001 real(kind=rp), dimension(lx, lx, lx) :: duds
1002 real(kind=rp), dimension(lx, lx, lx) :: dudt
1003 integer :: e, i, j, k
1004
1005 do e = 1, nelv
1006 do j = 1, lx * lx
1007 do i = 1, lx
1008 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
1009 + dx(i,2) * u(2,j,1,e) &
1010 + dx(i,3) * u(3,j,1,e) &
1011 + dx(i,4) * u(4,j,1,e) &
1012 + dx(i,5) * u(5,j,1,e)
1013 end do
1014 end do
1015
1016 do k = 1, lx
1017 do j = 1, lx
1018 do i = 1, lx
1019 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
1020 + dy(j,2) * u(i,2,k,e) &
1021 + dy(j,3) * u(i,3,k,e) &
1022 + dy(j,4) * u(i,4,k,e) &
1023 + dy(j,5) * u(i,5,k,e)
1024 end do
1025 end do
1026 end do
1027
1028 do k = 1, lx
1029 do i = 1, lx*lx
1030 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
1031 + dz(k,2) * u(i,1,2,e) &
1032 + dz(k,3) * u(i,1,3,e) &
1033 + dz(k,4) * u(i,1,4,e) &
1034 + dz(k,5) * u(i,1,5,e)
1035 end do
1036 end do
1037
1038 do i = 1, lx * lx * lx
1039 du(i,1,1,e) = jacinv(i,1,1,e) &
1040 * ( vx(i,1,1,e) &
1041 * ( drdx(i,1,1,e) * dudr(i,1,1) &
1042 + dsdx(i,1,1,e) * duds(i,1,1) &
1043 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
1044 + vy(i,1,1,e) &
1045 * ( drdy(i,1,1,e) * dudr(i,1,1) &
1046 + dsdy(i,1,1,e) * duds(i,1,1) &
1047 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
1048 + vz(i,1,1,e) &
1049 * ( drdz(i,1,1,e) * dudr(i,1,1) &
1050 + dsdz(i,1,1,e) * duds(i,1,1) &
1051 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
1052 end do
1053 end do
1054
1055 end subroutine cpu_conv1_lx5
1056
1057 subroutine cpu_conv1_lx4(du, u, vx, vy, vz, dx, dy, dz, &
1058 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1059 jacinv, nelv)
1060 integer, parameter :: lx = 4
1061 integer, intent(in) :: nelv
1062 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
1063 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
1064 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
1065 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
1066 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
1067 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
1068 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1069 real(kind=rp), dimension(lx, lx, lx) :: dudr
1070 real(kind=rp), dimension(lx, lx, lx) :: duds
1071 real(kind=rp), dimension(lx, lx, lx) :: dudt
1072 integer :: e, i, j, k
1073
1074 do e = 1, nelv
1075 do j = 1, lx * lx
1076 do i = 1, lx
1077 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
1078 + dx(i,2) * u(2,j,1,e) &
1079 + dx(i,3) * u(3,j,1,e) &
1080 + dx(i,4) * u(4,j,1,e)
1081 end do
1082 end do
1083
1084 do k = 1, lx
1085 do j = 1, lx
1086 do i = 1, lx
1087 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
1088 + dy(j,2) * u(i,2,k,e) &
1089 + dy(j,3) * u(i,3,k,e) &
1090 + dy(j,4) * u(i,4,k,e)
1091 end do
1092 end do
1093 end do
1094
1095 do k = 1, lx
1096 do i = 1, lx*lx
1097 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
1098 + dz(k,2) * u(i,1,2,e) &
1099 + dz(k,3) * u(i,1,3,e) &
1100 + dz(k,4) * u(i,1,4,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_lx4
1122
1123 subroutine cpu_conv1_lx3(du, u, vx, vy, vz, dx, dy, dz, &
1124 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1125 jacinv, nelv)
1126 integer, parameter :: lx = 3
1127 integer, intent(in) :: nelv
1128 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
1129 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
1130 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
1131 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
1132 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
1133 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
1134 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1135 real(kind=rp), dimension(lx, lx, lx) :: dudr
1136 real(kind=rp), dimension(lx, lx, lx) :: duds
1137 real(kind=rp), dimension(lx, lx, lx) :: dudt
1138 integer :: e, i, j, k
1139
1140 do e = 1, nelv
1141 do j = 1, lx * lx
1142 do i = 1, lx
1143 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
1144 + dx(i,2) * u(2,j,1,e) &
1145 + dx(i,3) * u(3,j,1,e)
1146 end do
1147 end do
1148
1149 do k = 1, lx
1150 do j = 1, lx
1151 do i = 1, lx
1152 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
1153 + dy(j,2) * u(i,2,k,e) &
1154 + dy(j,3) * u(i,3,k,e)
1155 end do
1156 end do
1157 end do
1158
1159 do k = 1, lx
1160 do i = 1, lx*lx
1161 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
1162 + dz(k,2) * u(i,1,2,e) &
1163 + dz(k,3) * u(i,1,3,e)
1164 end do
1165 end do
1166
1167 do i = 1, lx * lx * lx
1168 du(i,1,1,e) = jacinv(i,1,1,e) &
1169 * ( vx(i,1,1,e) &
1170 * ( drdx(i,1,1,e) * dudr(i,1,1) &
1171 + dsdx(i,1,1,e) * duds(i,1,1) &
1172 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
1173 + vy(i,1,1,e) &
1174 * ( drdy(i,1,1,e) * dudr(i,1,1) &
1175 + dsdy(i,1,1,e) * duds(i,1,1) &
1176 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
1177 + vz(i,1,1,e) &
1178 * ( drdz(i,1,1,e) * dudr(i,1,1) &
1179 + dsdz(i,1,1,e) * duds(i,1,1) &
1180 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
1181 end do
1182 end do
1183
1184 end subroutine cpu_conv1_lx3
1185
1186 subroutine cpu_conv1_lx2(du, u, vx, vy, vz, dx, dy, dz, &
1187 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1188 jacinv, nelv)
1189 integer, parameter :: lx = 2
1190 integer, intent(in) :: nelv
1191 real(kind=rp), dimension(lx, lx, lx, nelv), intent(inout) :: du
1192 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, vx, vy, vz
1193 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
1194 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
1195 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
1196 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
1197 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1198 real(kind=rp), dimension(lx, lx, lx) :: dudr
1199 real(kind=rp), dimension(lx, lx, lx) :: duds
1200 real(kind=rp), dimension(lx, lx, lx) :: dudt
1201 integer :: e, i, j, k
1202
1203 do e = 1, nelv
1204 do j = 1, lx * lx
1205 do i = 1, lx
1206 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
1207 + dx(i,2) * u(2,j,1,e)
1208 end do
1209 end do
1210
1211 do k = 1, lx
1212 do j = 1, lx
1213 do i = 1, lx
1214 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
1215 + dy(j,2) * u(i,2,k,e)
1216 end do
1217 end do
1218 end do
1219
1220 do k = 1, lx
1221 do i = 1, lx*lx
1222 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
1223 + dz(k,2) * u(i,1,2,e)
1224 end do
1225 end do
1226
1227 do i = 1, lx * lx * lx
1228 du(i,1,1,e) = jacinv(i,1,1,e) &
1229 * ( vx(i,1,1,e) &
1230 * ( drdx(i,1,1,e) * dudr(i,1,1) &
1231 + dsdx(i,1,1,e) * duds(i,1,1) &
1232 + dtdx(i,1,1,e) * dudt(i,1,1) ) &
1233 + vy(i,1,1,e) &
1234 * ( drdy(i,1,1,e) * dudr(i,1,1) &
1235 + dsdy(i,1,1,e) * duds(i,1,1) &
1236 + dtdy(i,1,1,e) * dudt(i,1,1) ) &
1237 + vz(i,1,1,e) &
1238 * ( drdz(i,1,1,e) * dudr(i,1,1) &
1239 + dsdz(i,1,1,e) * duds(i,1,1) &
1240 + dtdz(i,1,1,e) * dudt(i,1,1) ) )
1241 end do
1242 end do
1243
1244 end subroutine cpu_conv1_lx2
1245
1246 subroutine opr_cpu_conv1_single(du, u, vx, vy, vz, Xh, coef, e)
1247 integer, parameter :: e_len = 1
1248 type(space_t), intent(in) :: Xh
1249 type(coef_t), intent(in) :: coef
1250 integer, intent(in) :: e
1251 real(kind=rp), intent(inout) :: du(xh%lxyz, e_len)
1252 real(kind=rp), intent(inout) :: u(xh%lx, xh%ly, xh%lz, e_len)
1253 real(kind=rp), intent(inout) :: vx(xh%lx, xh%ly, xh%lz, e_len)
1254 real(kind=rp), intent(inout) :: vy(xh%lx, xh%ly, xh%lz, e_len)
1255 real(kind=rp), intent(inout) :: vz(xh%lx, xh%ly, xh%lz, e_len)
1256
1257 associate(drdx => coef%drdx, drdy => coef%drdy, drdz => coef%drdz, &
1258 dsdx => coef%dsdx, dsdy => coef%dsdy, dsdz => coef%dsdz, &
1259 dtdx => coef%dtdx, dtdy => coef%dtdy, dtdz => coef%dtdz, &
1260 jacinv => coef%jacinv)
1261 select case (xh%lx)
1262 case (14)
1263 call cpu_conv1_lx14_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1264 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1265 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1266 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1267 jacinv(1,1,1,e))
1268 case (13)
1269 call cpu_conv1_lx13_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1270 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1271 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1272 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1273 jacinv(1,1,1,e))
1274 case (12)
1275 call cpu_conv1_lx12_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1276 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1277 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1278 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1279 jacinv(1,1,1,e))
1280 case (11)
1281 call cpu_conv1_lx11_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1282 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1283 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1284 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1285 jacinv(1,1,1,e))
1286 case (10)
1287 call cpu_conv1_lx10_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1288 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1289 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1290 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1291 jacinv(1,1,1,e))
1292 case (9)
1293 call cpu_conv1_lx9_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1294 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1295 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1296 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1297 jacinv(1,1,1,e))
1298 case (8)
1299 call cpu_conv1_lx8_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1300 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1301 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1302 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1303 jacinv(1,1,1,e))
1304 case (7)
1305 call cpu_conv1_lx7_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1306 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1307 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1308 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1309 jacinv(1,1,1,e))
1310 case (6)
1311 call cpu_conv1_lx6_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1312 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1313 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1314 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1315 jacinv(1,1,1,e))
1316 case (5)
1317 call cpu_conv1_lx5_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1318 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1319 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1320 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1321 jacinv(1,1,1,e))
1322 case (4)
1323 call cpu_conv1_lx4_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1324 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1325 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1326 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1327 jacinv(1,1,1,e))
1328 case (3)
1329 call cpu_conv1_lx3_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1330 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1331 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1332 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1333 jacinv(1,1,1,e))
1334 case (2)
1335 call cpu_conv1_lx2_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1336 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1337 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1338 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1339 jacinv(1,1,1,e))
1340 case default
1341 call cpu_conv1_lx_single(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
1342 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e),&
1343 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e),&
1344 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e),&
1345 jacinv(1,1,1,e), xh%lx)
1346 end select
1347 end associate
1348
1349 end subroutine opr_cpu_conv1_single
1350
1351 subroutine cpu_conv1_lx_single(du, u, vx, vy, vz, dx, dy, dz, &
1352 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1353 jacinv, lx)
1354 integer, intent(in) :: lx
1355 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
1356 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
1357 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
1358 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
1359 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
1360 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
1361 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1362 real(kind=rp), dimension(lx, lx, lx) :: dudr
1363 real(kind=rp), dimension(lx, lx, lx) :: duds
1364 real(kind=rp), dimension(lx, lx, lx) :: dudt
1365 real(kind=rp) :: tmp
1366 integer :: i, j, k, l
1367
1368 do j = 1, lx * lx
1369 do i = 1, lx
1370 tmp = 0.0_rp
1371 do k = 1, lx
1372 tmp = tmp + dx(i,k) * u(k,j,1)
1373 end do
1374 dudr(i,j,1) = tmp
1375 end do
1376 end do
1377
1378 do k = 1, lx
1379 do j = 1, lx
1380 do i = 1, lx
1381 tmp = 0.0_rp
1382 do l = 1, lx
1383 tmp = tmp + dy(j,l) * u(i,l,k)
1384 end do
1385 duds(i,j,k) = tmp
1386 end do
1387 end do
1388 end do
1389
1390 do k = 1, lx
1391 do i = 1, lx*lx
1392 tmp = 0.0_rp
1393 do l = 1, lx
1394 tmp = tmp + dz(k,l) * u(i,1,l)
1395 end do
1396 dudt(i,1,k) = tmp
1397 end do
1398 end do
1399
1400 do i = 1, lx * lx * lx
1401 du(i,1,1) = jacinv(i,1,1) &
1402 * ( vx(i,1,1) &
1403 * ( drdx(i,1,1) * dudr(i,1,1) &
1404 + dsdx(i,1,1) * duds(i,1,1) &
1405 + dtdx(i,1,1) * dudt(i,1,1) ) &
1406 + vy(i,1,1) &
1407 * ( drdy(i,1,1) * dudr(i,1,1) &
1408 + dsdy(i,1,1) * duds(i,1,1) &
1409 + dtdy(i,1,1) * dudt(i,1,1) ) &
1410 + vz(i,1,1) &
1411 * ( drdz(i,1,1) * dudr(i,1,1) &
1412 + dsdz(i,1,1) * duds(i,1,1) &
1413 + dtdz(i,1,1) * dudt(i,1,1) ) )
1414 end do
1415
1416 end subroutine cpu_conv1_lx_single
1417
1418 subroutine cpu_conv1_lx14_single(du, u, vx, vy, vz, dx, dy, dz, &
1419 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1420 integer, parameter :: lx = 14
1421 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
1422 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
1423 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
1424 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
1425 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
1426 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
1427 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1428 real(kind=rp), dimension(lx, lx, lx) :: dudr
1429 real(kind=rp), dimension(lx, lx, lx) :: duds
1430 real(kind=rp), dimension(lx, lx, lx) :: dudt
1431 integer :: i, j, k
1432
1433 do j = 1, lx * lx
1434 do i = 1, lx
1435 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1436 + dx(i,2) * u(2,j,1) &
1437 + dx(i,3) * u(3,j,1) &
1438 + dx(i,4) * u(4,j,1) &
1439 + dx(i,5) * u(5,j,1) &
1440 + dx(i,6) * u(6,j,1) &
1441 + dx(i,7) * u(7,j,1) &
1442 + dx(i,8) * u(8,j,1) &
1443 + dx(i,9) * u(9,j,1) &
1444 + dx(i,10) * u(10,j,1) &
1445 + dx(i,11) * u(11,j,1) &
1446 + dx(i,12) * u(12,j,1) &
1447 + dx(i,13) * u(13,j,1) &
1448 + dx(i,14) * u(14,j,1)
1449 end do
1450 end do
1451
1452 do k = 1, lx
1453 do j = 1, lx
1454 do i = 1, lx
1455 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1456 + dy(j,2) * u(i,2,k) &
1457 + dy(j,3) * u(i,3,k) &
1458 + dy(j,4) * u(i,4,k) &
1459 + dy(j,5) * u(i,5,k) &
1460 + dy(j,6) * u(i,6,k) &
1461 + dy(j,7) * u(i,7,k) &
1462 + dy(j,8) * u(i,8,k) &
1463 + dy(j,9) * u(i,9,k) &
1464 + dy(j,10) * u(i,10,k) &
1465 + dy(j,11) * u(i,11,k) &
1466 + dy(j,12) * u(i,12,k) &
1467 + dy(j,13) * u(i,13,k) &
1468 + dy(j,14) * u(i,14,k)
1469 end do
1470 end do
1471 end do
1472
1473 do k = 1, lx
1474 do i = 1, lx*lx
1475 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1476 + dz(k,2) * u(i,1,2) &
1477 + dz(k,3) * u(i,1,3) &
1478 + dz(k,4) * u(i,1,4) &
1479 + dz(k,5) * u(i,1,5) &
1480 + dz(k,6) * u(i,1,6) &
1481 + dz(k,7) * u(i,1,7) &
1482 + dz(k,8) * u(i,1,8) &
1483 + dz(k,9) * u(i,1,9) &
1484 + dz(k,10) * u(i,1,10) &
1485 + dz(k,11) * u(i,1,11) &
1486 + dz(k,12) * u(i,1,12) &
1487 + dz(k,13) * u(i,1,13) &
1488 + dz(k,14) * u(i,1,14)
1489 end do
1490 end do
1491
1492 do i = 1, lx * lx * lx
1493 du(i,1,1) = jacinv(i,1,1) &
1494 * ( vx(i,1,1) &
1495 * ( drdx(i,1,1) * dudr(i,1,1) &
1496 + dsdx(i,1,1) * duds(i,1,1) &
1497 + dtdx(i,1,1) * dudt(i,1,1) ) &
1498 + vy(i,1,1) &
1499 * ( drdy(i,1,1) * dudr(i,1,1) &
1500 + dsdy(i,1,1) * duds(i,1,1) &
1501 + dtdy(i,1,1) * dudt(i,1,1) ) &
1502 + vz(i,1,1) &
1503 * ( drdz(i,1,1) * dudr(i,1,1) &
1504 + dsdz(i,1,1) * duds(i,1,1) &
1505 + dtdz(i,1,1) * dudt(i,1,1) ) )
1506 end do
1507
1508 end subroutine cpu_conv1_lx14_single
1509
1510 subroutine cpu_conv1_lx13_single(du, u, vx, vy, vz, dx, dy, dz, &
1511 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1512 integer, parameter :: lx = 13
1513 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
1514 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
1515 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
1516 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
1517 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
1518 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
1519 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1520 real(kind=rp), dimension(lx, lx, lx) :: dudr
1521 real(kind=rp), dimension(lx, lx, lx) :: duds
1522 real(kind=rp), dimension(lx, lx, lx) :: dudt
1523 integer :: i, j, k
1524
1525 do j = 1, lx * lx
1526 do i = 1, lx
1527 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1528 + dx(i,2) * u(2,j,1) &
1529 + dx(i,3) * u(3,j,1) &
1530 + dx(i,4) * u(4,j,1) &
1531 + dx(i,5) * u(5,j,1) &
1532 + dx(i,6) * u(6,j,1) &
1533 + dx(i,7) * u(7,j,1) &
1534 + dx(i,8) * u(8,j,1) &
1535 + dx(i,9) * u(9,j,1) &
1536 + dx(i,10) * u(10,j,1) &
1537 + dx(i,11) * u(11,j,1) &
1538 + dx(i,12) * u(12,j,1) &
1539 + dx(i,13) * u(13,j,1)
1540 end do
1541 end do
1542
1543 do k = 1, lx
1544 do j = 1, lx
1545 do i = 1, lx
1546 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1547 + dy(j,2) * u(i,2,k) &
1548 + dy(j,3) * u(i,3,k) &
1549 + dy(j,4) * u(i,4,k) &
1550 + dy(j,5) * u(i,5,k) &
1551 + dy(j,6) * u(i,6,k) &
1552 + dy(j,7) * u(i,7,k) &
1553 + dy(j,8) * u(i,8,k) &
1554 + dy(j,9) * u(i,9,k) &
1555 + dy(j,10) * u(i,10,k) &
1556 + dy(j,11) * u(i,11,k) &
1557 + dy(j,12) * u(i,12,k) &
1558 + dy(j,13) * u(i,13,k)
1559 end do
1560 end do
1561 end do
1562
1563 do k = 1, lx
1564 do i = 1, lx*lx
1565 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1566 + dz(k,2) * u(i,1,2) &
1567 + dz(k,3) * u(i,1,3) &
1568 + dz(k,4) * u(i,1,4) &
1569 + dz(k,5) * u(i,1,5) &
1570 + dz(k,6) * u(i,1,6) &
1571 + dz(k,7) * u(i,1,7) &
1572 + dz(k,8) * u(i,1,8) &
1573 + dz(k,9) * u(i,1,9) &
1574 + dz(k,10) * u(i,1,10) &
1575 + dz(k,11) * u(i,1,11) &
1576 + dz(k,12) * u(i,1,12) &
1577 + dz(k,13) * u(i,1,13)
1578 end do
1579 end do
1580
1581 do i = 1, lx * lx * lx
1582 du(i,1,1) = jacinv(i,1,1) &
1583 * ( vx(i,1,1) &
1584 * ( drdx(i,1,1) * dudr(i,1,1) &
1585 + dsdx(i,1,1) * duds(i,1,1) &
1586 + dtdx(i,1,1) * dudt(i,1,1) ) &
1587 + vy(i,1,1) &
1588 * ( drdy(i,1,1) * dudr(i,1,1) &
1589 + dsdy(i,1,1) * duds(i,1,1) &
1590 + dtdy(i,1,1) * dudt(i,1,1) ) &
1591 + vz(i,1,1) &
1592 * ( drdz(i,1,1) * dudr(i,1,1) &
1593 + dsdz(i,1,1) * duds(i,1,1) &
1594 + dtdz(i,1,1) * dudt(i,1,1) ) )
1595 end do
1596
1597 end subroutine cpu_conv1_lx13_single
1598
1599
1600 subroutine cpu_conv1_lx12_single(du, u, vx, vy, vz, dx, dy, dz, &
1601 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1602 integer, parameter :: lx = 12
1603 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
1604 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
1605 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
1606 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
1607 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
1608 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
1609 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1610 real(kind=rp), dimension(lx, lx, lx) :: dudr
1611 real(kind=rp), dimension(lx, lx, lx) :: duds
1612 real(kind=rp), dimension(lx, lx, lx) :: dudt
1613 integer :: i, j, k
1614
1615 do j = 1, lx * lx
1616 do i = 1, lx
1617 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1618 + dx(i,2) * u(2,j,1) &
1619 + dx(i,3) * u(3,j,1) &
1620 + dx(i,4) * u(4,j,1) &
1621 + dx(i,5) * u(5,j,1) &
1622 + dx(i,6) * u(6,j,1) &
1623 + dx(i,7) * u(7,j,1) &
1624 + dx(i,8) * u(8,j,1) &
1625 + dx(i,9) * u(9,j,1) &
1626 + dx(i,10) * u(10,j,1) &
1627 + dx(i,11) * u(11,j,1) &
1628 + dx(i,12) * u(12,j,1)
1629 end do
1630 end do
1631
1632 do k = 1, lx
1633 do j = 1, lx
1634 do i = 1, lx
1635 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1636 + dy(j,2) * u(i,2,k) &
1637 + dy(j,3) * u(i,3,k) &
1638 + dy(j,4) * u(i,4,k) &
1639 + dy(j,5) * u(i,5,k) &
1640 + dy(j,6) * u(i,6,k) &
1641 + dy(j,7) * u(i,7,k) &
1642 + dy(j,8) * u(i,8,k) &
1643 + dy(j,9) * u(i,9,k) &
1644 + dy(j,10) * u(i,10,k) &
1645 + dy(j,11) * u(i,11,k) &
1646 + dy(j,12) * u(i,12,k)
1647 end do
1648 end do
1649 end do
1650
1651 do k = 1, lx
1652 do i = 1, lx*lx
1653 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1654 + dz(k,2) * u(i,1,2) &
1655 + dz(k,3) * u(i,1,3) &
1656 + dz(k,4) * u(i,1,4) &
1657 + dz(k,5) * u(i,1,5) &
1658 + dz(k,6) * u(i,1,6) &
1659 + dz(k,7) * u(i,1,7) &
1660 + dz(k,8) * u(i,1,8) &
1661 + dz(k,9) * u(i,1,9) &
1662 + dz(k,10) * u(i,1,10) &
1663 + dz(k,11) * u(i,1,11) &
1664 + dz(k,12) * u(i,1,12)
1665 end do
1666 end do
1667
1668 do i = 1, lx * lx * lx
1669 du(i,1,1) = jacinv(i,1,1) &
1670 * ( vx(i,1,1) &
1671 * ( drdx(i,1,1) * dudr(i,1,1) &
1672 + dsdx(i,1,1) * duds(i,1,1) &
1673 + dtdx(i,1,1) * dudt(i,1,1) ) &
1674 + vy(i,1,1) &
1675 * ( drdy(i,1,1) * dudr(i,1,1) &
1676 + dsdy(i,1,1) * duds(i,1,1) &
1677 + dtdy(i,1,1) * dudt(i,1,1) ) &
1678 + vz(i,1,1) &
1679 * ( drdz(i,1,1) * dudr(i,1,1) &
1680 + dsdz(i,1,1) * duds(i,1,1) &
1681 + dtdz(i,1,1) * dudt(i,1,1) ) )
1682 end do
1683
1684 end subroutine cpu_conv1_lx12_single
1685
1686 subroutine cpu_conv1_lx11_single(du, u, vx, vy, vz, dx, dy, dz, &
1687 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1688 integer, parameter :: lx = 11
1689 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
1690 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
1691 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
1692 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
1693 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
1694 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
1695 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1696 real(kind=rp), dimension(lx, lx, lx) :: dudr
1697 real(kind=rp), dimension(lx, lx, lx) :: duds
1698 real(kind=rp), dimension(lx, lx, lx) :: dudt
1699 integer :: i, j, k
1700
1701 do j = 1, lx * lx
1702 do i = 1, lx
1703 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1704 + dx(i,2) * u(2,j,1) &
1705 + dx(i,3) * u(3,j,1) &
1706 + dx(i,4) * u(4,j,1) &
1707 + dx(i,5) * u(5,j,1) &
1708 + dx(i,6) * u(6,j,1) &
1709 + dx(i,7) * u(7,j,1) &
1710 + dx(i,8) * u(8,j,1) &
1711 + dx(i,9) * u(9,j,1) &
1712 + dx(i,10) * u(10,j,1) &
1713 + dx(i,11) * u(11,j,1)
1714 end do
1715 end do
1716
1717 do k = 1, lx
1718 do j = 1, lx
1719 do i = 1, lx
1720 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1721 + dy(j,2) * u(i,2,k) &
1722 + dy(j,3) * u(i,3,k) &
1723 + dy(j,4) * u(i,4,k) &
1724 + dy(j,5) * u(i,5,k) &
1725 + dy(j,6) * u(i,6,k) &
1726 + dy(j,7) * u(i,7,k) &
1727 + dy(j,8) * u(i,8,k) &
1728 + dy(j,9) * u(i,9,k) &
1729 + dy(j,10) * u(i,10,k) &
1730 + dy(j,11) * u(i,11,k)
1731 end do
1732 end do
1733 end do
1734
1735 do k = 1, lx
1736 do i = 1, lx*lx
1737 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1738 + dz(k,2) * u(i,1,2) &
1739 + dz(k,3) * u(i,1,3) &
1740 + dz(k,4) * u(i,1,4) &
1741 + dz(k,5) * u(i,1,5) &
1742 + dz(k,6) * u(i,1,6) &
1743 + dz(k,7) * u(i,1,7) &
1744 + dz(k,8) * u(i,1,8) &
1745 + dz(k,9) * u(i,1,9) &
1746 + dz(k,10) * u(i,1,10) &
1747 + dz(k,11) * u(i,1,11)
1748 end do
1749 end do
1750
1751 do i = 1, lx * lx * lx
1752 du(i,1,1) = jacinv(i,1,1) &
1753 * ( vx(i,1,1) &
1754 * ( drdx(i,1,1) * dudr(i,1,1) &
1755 + dsdx(i,1,1) * duds(i,1,1) &
1756 + dtdx(i,1,1) * dudt(i,1,1) ) &
1757 + vy(i,1,1) &
1758 * ( drdy(i,1,1) * dudr(i,1,1) &
1759 + dsdy(i,1,1) * duds(i,1,1) &
1760 + dtdy(i,1,1) * dudt(i,1,1) ) &
1761 + vz(i,1,1) &
1762 * ( drdz(i,1,1) * dudr(i,1,1) &
1763 + dsdz(i,1,1) * duds(i,1,1) &
1764 + dtdz(i,1,1) * dudt(i,1,1) ) )
1765 end do
1766
1767 end subroutine cpu_conv1_lx11_single
1768
1769 subroutine cpu_conv1_lx10_single(du, u, vx, vy, vz, dx, dy, dz, &
1770 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1771 integer, parameter :: lx = 10
1772 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
1773 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
1774 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
1775 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
1776 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
1777 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
1778 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1779 real(kind=rp), dimension(lx, lx, lx) :: dudr
1780 real(kind=rp), dimension(lx, lx, lx) :: duds
1781 real(kind=rp), dimension(lx, lx, lx) :: dudt
1782 integer :: i, j, k
1783
1784 do j = 1, lx * lx
1785 do i = 1, lx
1786 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1787 + dx(i,2) * u(2,j,1) &
1788 + dx(i,3) * u(3,j,1) &
1789 + dx(i,4) * u(4,j,1) &
1790 + dx(i,5) * u(5,j,1) &
1791 + dx(i,6) * u(6,j,1) &
1792 + dx(i,7) * u(7,j,1) &
1793 + dx(i,8) * u(8,j,1) &
1794 + dx(i,9) * u(9,j,1) &
1795 + dx(i,10) * u(10,j,1)
1796 end do
1797 end do
1798
1799 do k = 1, lx
1800 do j = 1, lx
1801 do i = 1, lx
1802 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1803 + dy(j,2) * u(i,2,k) &
1804 + dy(j,3) * u(i,3,k) &
1805 + dy(j,4) * u(i,4,k) &
1806 + dy(j,5) * u(i,5,k) &
1807 + dy(j,6) * u(i,6,k) &
1808 + dy(j,7) * u(i,7,k) &
1809 + dy(j,8) * u(i,8,k) &
1810 + dy(j,9) * u(i,9,k) &
1811 + dy(j,10) * u(i,10,k)
1812 end do
1813 end do
1814 end do
1815
1816 do k = 1, lx
1817 do i = 1, lx*lx
1818 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1819 + dz(k,2) * u(i,1,2) &
1820 + dz(k,3) * u(i,1,3) &
1821 + dz(k,4) * u(i,1,4) &
1822 + dz(k,5) * u(i,1,5) &
1823 + dz(k,6) * u(i,1,6) &
1824 + dz(k,7) * u(i,1,7) &
1825 + dz(k,8) * u(i,1,8) &
1826 + dz(k,9) * u(i,1,9) &
1827 + dz(k,10) * u(i,1,10)
1828 end do
1829 end do
1830
1831 do i = 1, lx * lx * lx
1832 du(i,1,1) = jacinv(i,1,1) &
1833 * ( vx(i,1,1) &
1834 * ( drdx(i,1,1) * dudr(i,1,1) &
1835 + dsdx(i,1,1) * duds(i,1,1) &
1836 + dtdx(i,1,1) * dudt(i,1,1) ) &
1837 + vy(i,1,1) &
1838 * ( drdy(i,1,1) * dudr(i,1,1) &
1839 + dsdy(i,1,1) * duds(i,1,1) &
1840 + dtdy(i,1,1) * dudt(i,1,1) ) &
1841 + vz(i,1,1) &
1842 * ( drdz(i,1,1) * dudr(i,1,1) &
1843 + dsdz(i,1,1) * duds(i,1,1) &
1844 + dtdz(i,1,1) * dudt(i,1,1) ) )
1845 end do
1846
1847 end subroutine cpu_conv1_lx10_single
1848
1849 subroutine cpu_conv1_lx9_single(du, u, vx, vy, vz, dx, dy, dz, &
1850 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1851 integer, parameter :: lx = 9
1852 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
1853 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
1854 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
1855 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
1856 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
1857 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
1858 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1859 real(kind=rp), dimension(lx, lx, lx) :: dudr
1860 real(kind=rp), dimension(lx, lx, lx) :: duds
1861 real(kind=rp), dimension(lx, lx, lx) :: dudt
1862 integer :: i, j, k
1863
1864 do j = 1, lx * lx
1865 do i = 1, lx
1866 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1867 + dx(i,2) * u(2,j,1) &
1868 + dx(i,3) * u(3,j,1) &
1869 + dx(i,4) * u(4,j,1) &
1870 + dx(i,5) * u(5,j,1) &
1871 + dx(i,6) * u(6,j,1) &
1872 + dx(i,7) * u(7,j,1) &
1873 + dx(i,8) * u(8,j,1) &
1874 + dx(i,9) * u(9,j,1)
1875 end do
1876 end do
1877
1878 do k = 1, lx
1879 do j = 1, lx
1880 do i = 1, lx
1881 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1882 + dy(j,2) * u(i,2,k) &
1883 + dy(j,3) * u(i,3,k) &
1884 + dy(j,4) * u(i,4,k) &
1885 + dy(j,5) * u(i,5,k) &
1886 + dy(j,6) * u(i,6,k) &
1887 + dy(j,7) * u(i,7,k) &
1888 + dy(j,8) * u(i,8,k) &
1889 + dy(j,9) * u(i,9,k)
1890 end do
1891 end do
1892 end do
1893
1894 do k = 1, lx
1895 do i = 1, lx*lx
1896 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1897 + dz(k,2) * u(i,1,2) &
1898 + dz(k,3) * u(i,1,3) &
1899 + dz(k,4) * u(i,1,4) &
1900 + dz(k,5) * u(i,1,5) &
1901 + dz(k,6) * u(i,1,6) &
1902 + dz(k,7) * u(i,1,7) &
1903 + dz(k,8) * u(i,1,8) &
1904 + dz(k,9) * u(i,1,9)
1905 end do
1906 end do
1907
1908 do i = 1, lx * lx * lx
1909 du(i,1,1) = jacinv(i,1,1) &
1910 * ( vx(i,1,1) &
1911 * ( drdx(i,1,1) * dudr(i,1,1) &
1912 + dsdx(i,1,1) * duds(i,1,1) &
1913 + dtdx(i,1,1) * dudt(i,1,1) ) &
1914 + vy(i,1,1) &
1915 * ( drdy(i,1,1) * dudr(i,1,1) &
1916 + dsdy(i,1,1) * duds(i,1,1) &
1917 + dtdy(i,1,1) * dudt(i,1,1) ) &
1918 + vz(i,1,1) &
1919 * ( drdz(i,1,1) * dudr(i,1,1) &
1920 + dsdz(i,1,1) * duds(i,1,1) &
1921 + dtdz(i,1,1) * dudt(i,1,1) ) )
1922 end do
1923
1924 end subroutine cpu_conv1_lx9_single
1925
1926 subroutine cpu_conv1_lx8_single(du, u, vx, vy, vz, dx, dy, dz, &
1927 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
1928 integer, parameter :: lx = 8
1929 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
1930 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
1931 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
1932 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
1933 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
1934 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
1935 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1936 real(kind=rp), dimension(lx, lx, lx) :: dudr
1937 real(kind=rp), dimension(lx, lx, lx) :: duds
1938 real(kind=rp), dimension(lx, lx, lx) :: dudt
1939 integer :: i, j, k
1940
1941 do j = 1, lx * lx
1942 do i = 1, lx
1943 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
1944 + dx(i,2) * u(2,j,1) &
1945 + dx(i,3) * u(3,j,1) &
1946 + dx(i,4) * u(4,j,1) &
1947 + dx(i,5) * u(5,j,1) &
1948 + dx(i,6) * u(6,j,1) &
1949 + dx(i,7) * u(7,j,1) &
1950 + dx(i,8) * u(8,j,1)
1951 end do
1952 end do
1953
1954 do k = 1, lx
1955 do j = 1, lx
1956 do i = 1, lx
1957 duds(i,j,k) = dy(j,1) * u(i,1,k) &
1958 + dy(j,2) * u(i,2,k) &
1959 + dy(j,3) * u(i,3,k) &
1960 + dy(j,4) * u(i,4,k) &
1961 + dy(j,5) * u(i,5,k) &
1962 + dy(j,6) * u(i,6,k) &
1963 + dy(j,7) * u(i,7,k) &
1964 + dy(j,8) * u(i,8,k)
1965 end do
1966 end do
1967 end do
1968
1969 do k = 1, lx
1970 do i = 1, lx*lx
1971 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
1972 + dz(k,2) * u(i,1,2) &
1973 + dz(k,3) * u(i,1,3) &
1974 + dz(k,4) * u(i,1,4) &
1975 + dz(k,5) * u(i,1,5) &
1976 + dz(k,6) * u(i,1,6) &
1977 + dz(k,7) * u(i,1,7) &
1978 + dz(k,8) * u(i,1,8)
1979 end do
1980 end do
1981
1982 do i = 1, lx * lx * lx
1983 du(i,1,1) = jacinv(i,1,1) &
1984 * ( vx(i,1,1) &
1985 * ( drdx(i,1,1) * dudr(i,1,1) &
1986 + dsdx(i,1,1) * duds(i,1,1) &
1987 + dtdx(i,1,1) * dudt(i,1,1) ) &
1988 + vy(i,1,1) &
1989 * ( drdy(i,1,1) * dudr(i,1,1) &
1990 + dsdy(i,1,1) * duds(i,1,1) &
1991 + dtdy(i,1,1) * dudt(i,1,1) ) &
1992 + vz(i,1,1) &
1993 * ( drdz(i,1,1) * dudr(i,1,1) &
1994 + dsdz(i,1,1) * duds(i,1,1) &
1995 + dtdz(i,1,1) * dudt(i,1,1) ) )
1996 end do
1997
1998 end subroutine cpu_conv1_lx8_single
1999
2000 subroutine cpu_conv1_lx7_single(du, u, vx, vy, vz, dx, dy, dz, &
2001 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2002 integer, parameter :: lx = 7
2003 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
2004 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
2005 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2006 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2007 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2008 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
2009 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2010 real(kind=rp), dimension(lx, lx, lx) :: dudr
2011 real(kind=rp), dimension(lx, lx, lx) :: duds
2012 real(kind=rp), dimension(lx, lx, lx) :: dudt
2013 integer :: i, j, k
2014
2015 do j = 1, lx * lx
2016 do i = 1, lx
2017 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2018 + dx(i,2) * u(2,j,1) &
2019 + dx(i,3) * u(3,j,1) &
2020 + dx(i,4) * u(4,j,1) &
2021 + dx(i,5) * u(5,j,1) &
2022 + dx(i,6) * u(6,j,1) &
2023 + dx(i,7) * u(7,j,1)
2024 end do
2025 end do
2026
2027 do k = 1, lx
2028 do j = 1, lx
2029 do i = 1, lx
2030 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2031 + dy(j,2) * u(i,2,k) &
2032 + dy(j,3) * u(i,3,k) &
2033 + dy(j,4) * u(i,4,k) &
2034 + dy(j,5) * u(i,5,k) &
2035 + dy(j,6) * u(i,6,k) &
2036 + dy(j,7) * u(i,7,k)
2037 end do
2038 end do
2039 end do
2040
2041 do k = 1, lx
2042 do i = 1, lx*lx
2043 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2044 + dz(k,2) * u(i,1,2) &
2045 + dz(k,3) * u(i,1,3) &
2046 + dz(k,4) * u(i,1,4) &
2047 + dz(k,5) * u(i,1,5) &
2048 + dz(k,6) * u(i,1,6) &
2049 + dz(k,7) * u(i,1,7)
2050 end do
2051 end do
2052
2053 do i = 1, lx * lx * lx
2054 du(i,1,1) = jacinv(i,1,1) &
2055 * ( vx(i,1,1) &
2056 * ( drdx(i,1,1) * dudr(i,1,1) &
2057 + dsdx(i,1,1) * duds(i,1,1) &
2058 + dtdx(i,1,1) * dudt(i,1,1) ) &
2059 + vy(i,1,1) &
2060 * ( drdy(i,1,1) * dudr(i,1,1) &
2061 + dsdy(i,1,1) * duds(i,1,1) &
2062 + dtdy(i,1,1) * dudt(i,1,1) ) &
2063 + vz(i,1,1) &
2064 * ( drdz(i,1,1) * dudr(i,1,1) &
2065 + dsdz(i,1,1) * duds(i,1,1) &
2066 + dtdz(i,1,1) * dudt(i,1,1) ) )
2067 end do
2068
2069 end subroutine cpu_conv1_lx7_single
2070
2071 subroutine cpu_conv1_lx6_single(du, u, vx, vy, vz, dx, dy, dz, &
2072 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2073 integer, parameter :: lx = 6
2074 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
2075 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
2076 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2077 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2078 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2079 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
2080 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2081 real(kind=rp), dimension(lx, lx, lx) :: dudr
2082 real(kind=rp), dimension(lx, lx, lx) :: duds
2083 real(kind=rp), dimension(lx, lx, lx) :: dudt
2084 integer :: i, j, k
2085
2086 do j = 1, lx * lx
2087 do i = 1, lx
2088 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2089 + dx(i,2) * u(2,j,1) &
2090 + dx(i,3) * u(3,j,1) &
2091 + dx(i,4) * u(4,j,1) &
2092 + dx(i,5) * u(5,j,1) &
2093 + dx(i,6) * u(6,j,1)
2094 end do
2095 end do
2096
2097 do k = 1, lx
2098 do j = 1, lx
2099 do i = 1, lx
2100 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2101 + dy(j,2) * u(i,2,k) &
2102 + dy(j,3) * u(i,3,k) &
2103 + dy(j,4) * u(i,4,k) &
2104 + dy(j,5) * u(i,5,k) &
2105 + dy(j,6) * u(i,6,k)
2106 end do
2107 end do
2108 end do
2109
2110 do k = 1, lx
2111 do i = 1, lx*lx
2112 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2113 + dz(k,2) * u(i,1,2) &
2114 + dz(k,3) * u(i,1,3) &
2115 + dz(k,4) * u(i,1,4) &
2116 + dz(k,5) * u(i,1,5) &
2117 + dz(k,6) * u(i,1,6)
2118 end do
2119 end do
2120
2121 do i = 1, lx * lx * lx
2122 du(i,1,1) = jacinv(i,1,1) &
2123 * ( vx(i,1,1) &
2124 * ( drdx(i,1,1) * dudr(i,1,1) &
2125 + dsdx(i,1,1) * duds(i,1,1) &
2126 + dtdx(i,1,1) * dudt(i,1,1) ) &
2127 + vy(i,1,1) &
2128 * ( drdy(i,1,1) * dudr(i,1,1) &
2129 + dsdy(i,1,1) * duds(i,1,1) &
2130 + dtdy(i,1,1) * dudt(i,1,1) ) &
2131 + vz(i,1,1) &
2132 * ( drdz(i,1,1) * dudr(i,1,1) &
2133 + dsdz(i,1,1) * duds(i,1,1) &
2134 + dtdz(i,1,1) * dudt(i,1,1) ) )
2135 end do
2136
2137 end subroutine cpu_conv1_lx6_single
2138
2139 subroutine cpu_conv1_lx5_single(du, u, vx, vy, vz, dx, dy, dz, &
2140 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2141 integer, parameter :: lx = 5
2142 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
2143 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
2144 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2145 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2146 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2147 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
2148 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2149 real(kind=rp), dimension(lx, lx, lx) :: dudr
2150 real(kind=rp), dimension(lx, lx, lx) :: duds
2151 real(kind=rp), dimension(lx, lx, lx) :: dudt
2152 integer :: i, j, k
2153
2154 do j = 1, lx * lx
2155 do i = 1, lx
2156 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2157 + dx(i,2) * u(2,j,1) &
2158 + dx(i,3) * u(3,j,1) &
2159 + dx(i,4) * u(4,j,1) &
2160 + dx(i,5) * u(5,j,1)
2161 end do
2162 end do
2163
2164 do k = 1, lx
2165 do j = 1, lx
2166 do i = 1, lx
2167 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2168 + dy(j,2) * u(i,2,k) &
2169 + dy(j,3) * u(i,3,k) &
2170 + dy(j,4) * u(i,4,k) &
2171 + dy(j,5) * u(i,5,k)
2172 end do
2173 end do
2174 end do
2175
2176 do k = 1, lx
2177 do i = 1, lx*lx
2178 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2179 + dz(k,2) * u(i,1,2) &
2180 + dz(k,3) * u(i,1,3) &
2181 + dz(k,4) * u(i,1,4) &
2182 + dz(k,5) * u(i,1,5)
2183 end do
2184 end do
2185
2186 do i = 1, lx * lx * lx
2187 du(i,1,1) = jacinv(i,1,1) &
2188 * ( vx(i,1,1) &
2189 * ( drdx(i,1,1) * dudr(i,1,1) &
2190 + dsdx(i,1,1) * duds(i,1,1) &
2191 + dtdx(i,1,1) * dudt(i,1,1) ) &
2192 + vy(i,1,1) &
2193 * ( drdy(i,1,1) * dudr(i,1,1) &
2194 + dsdy(i,1,1) * duds(i,1,1) &
2195 + dtdy(i,1,1) * dudt(i,1,1) ) &
2196 + vz(i,1,1) &
2197 * ( drdz(i,1,1) * dudr(i,1,1) &
2198 + dsdz(i,1,1) * duds(i,1,1) &
2199 + dtdz(i,1,1) * dudt(i,1,1) ) )
2200 end do
2201
2202 end subroutine cpu_conv1_lx5_single
2203
2204 subroutine cpu_conv1_lx4_single(du, u, vx, vy, vz, dx, dy, dz, &
2205 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2206 integer, parameter :: lx = 4
2207 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
2208 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
2209 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2210 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2211 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2212 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
2213 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2214 real(kind=rp), dimension(lx, lx, lx) :: dudr
2215 real(kind=rp), dimension(lx, lx, lx) :: duds
2216 real(kind=rp), dimension(lx, lx, lx) :: dudt
2217 integer :: i, j, k
2218
2219 do j = 1, lx * lx
2220 do i = 1, lx
2221 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2222 + dx(i,2) * u(2,j,1) &
2223 + dx(i,3) * u(3,j,1) &
2224 + dx(i,4) * u(4,j,1)
2225 end do
2226 end do
2227
2228 do k = 1, lx
2229 do j = 1, lx
2230 do i = 1, lx
2231 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2232 + dy(j,2) * u(i,2,k) &
2233 + dy(j,3) * u(i,3,k) &
2234 + dy(j,4) * u(i,4,k)
2235 end do
2236 end do
2237 end do
2238
2239 do k = 1, lx
2240 do i = 1, lx*lx
2241 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2242 + dz(k,2) * u(i,1,2) &
2243 + dz(k,3) * u(i,1,3) &
2244 + dz(k,4) * u(i,1,4)
2245 end do
2246 end do
2247
2248 do i = 1, lx * lx * lx
2249 du(i,1,1) = jacinv(i,1,1) &
2250 * ( vx(i,1,1) &
2251 * ( drdx(i,1,1) * dudr(i,1,1) &
2252 + dsdx(i,1,1) * duds(i,1,1) &
2253 + dtdx(i,1,1) * dudt(i,1,1) ) &
2254 + vy(i,1,1) &
2255 * ( drdy(i,1,1) * dudr(i,1,1) &
2256 + dsdy(i,1,1) * duds(i,1,1) &
2257 + dtdy(i,1,1) * dudt(i,1,1) ) &
2258 + vz(i,1,1) &
2259 * ( drdz(i,1,1) * dudr(i,1,1) &
2260 + dsdz(i,1,1) * duds(i,1,1) &
2261 + dtdz(i,1,1) * dudt(i,1,1) ) )
2262 end do
2263
2264 end subroutine cpu_conv1_lx4_single
2265
2266 subroutine cpu_conv1_lx3_single(du, u, vx, vy, vz, dx, dy, dz, &
2267 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2268 integer, parameter :: lx = 3
2269 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
2270 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
2271 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2272 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2273 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2274 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
2275 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2276 real(kind=rp), dimension(lx, lx, lx) :: dudr
2277 real(kind=rp), dimension(lx, lx, lx) :: duds
2278 real(kind=rp), dimension(lx, lx, lx) :: dudt
2279 integer :: i, j, k
2280
2281 do j = 1, lx * lx
2282 do i = 1, lx
2283 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2284 + dx(i,2) * u(2,j,1) &
2285 + dx(i,3) * u(3,j,1)
2286 end do
2287 end do
2288
2289 do k = 1, lx
2290 do j = 1, lx
2291 do i = 1, lx
2292 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2293 + dy(j,2) * u(i,2,k) &
2294 + dy(j,3) * u(i,3,k)
2295 end do
2296 end do
2297 end do
2298
2299 do k = 1, lx
2300 do i = 1, lx*lx
2301 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2302 + dz(k,2) * u(i,1,2) &
2303 + dz(k,3) * u(i,1,3)
2304 end do
2305 end do
2306
2307 do i = 1, lx * lx * lx
2308 du(i,1,1) = jacinv(i,1,1) &
2309 * ( vx(i,1,1) &
2310 * ( drdx(i,1,1) * dudr(i,1,1) &
2311 + dsdx(i,1,1) * duds(i,1,1) &
2312 + dtdx(i,1,1) * dudt(i,1,1) ) &
2313 + vy(i,1,1) &
2314 * ( drdy(i,1,1) * dudr(i,1,1) &
2315 + dsdy(i,1,1) * duds(i,1,1) &
2316 + dtdy(i,1,1) * dudt(i,1,1) ) &
2317 + vz(i,1,1) &
2318 * ( drdz(i,1,1) * dudr(i,1,1) &
2319 + dsdz(i,1,1) * duds(i,1,1) &
2320 + dtdz(i,1,1) * dudt(i,1,1) ) )
2321 end do
2322
2323 end subroutine cpu_conv1_lx3_single
2324
2325 subroutine cpu_conv1_lx2_single(du, u, vx, vy, vz, dx, dy, dz, &
2326 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, jacinv)
2327 integer, parameter :: lx = 2
2328 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: du
2329 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u, vx, vy, vz
2330 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2331 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2332 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2333 real(kind=rp), dimension(lx, lx, lx), intent(in) :: jacinv
2334 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2335 real(kind=rp), dimension(lx, lx, lx) :: dudr
2336 real(kind=rp), dimension(lx, lx, lx) :: duds
2337 real(kind=rp), dimension(lx, lx, lx) :: dudt
2338 integer :: i, j, k
2339
2340 do j = 1, lx * lx
2341 do i = 1, lx
2342 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2343 + dx(i,2) * u(2,j,1)
2344 end do
2345 end do
2346
2347 do k = 1, lx
2348 do j = 1, lx
2349 do i = 1, lx
2350 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2351 + dy(j,2) * u(i,2,k)
2352 end do
2353 end do
2354 end do
2355
2356 do k = 1, lx
2357 do i = 1, lx*lx
2358 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2359 + dz(k,2) * u(i,1,2)
2360 end do
2361 end do
2362
2363 do i = 1, lx * lx * lx
2364 du(i,1,1) = jacinv(i,1,1) &
2365 * ( vx(i,1,1) &
2366 * ( drdx(i,1,1) * dudr(i,1,1) &
2367 + dsdx(i,1,1) * duds(i,1,1) &
2368 + dtdx(i,1,1) * dudt(i,1,1) ) &
2369 + vy(i,1,1) &
2370 * ( drdy(i,1,1) * dudr(i,1,1) &
2371 + dsdy(i,1,1) * duds(i,1,1) &
2372 + dtdy(i,1,1) * dudt(i,1,1) ) &
2373 + vz(i,1,1) &
2374 * ( drdz(i,1,1) * dudr(i,1,1) &
2375 + dsdz(i,1,1) * duds(i,1,1) &
2376 + dtdz(i,1,1) * dudt(i,1,1) ) )
2377 end do
2378
2379 end subroutine cpu_conv1_lx2_single
2380
2381end submodule cpu_conv1
Operators CPU backend.
Definition opr_cpu.f90:34