Neko 1.99.3
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-2026, 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 !$omp parallel
74 select case (xh%lx)
75 case (14)
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)
81 case (13)
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)
87 case (12)
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)
93 case (11)
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)
99 case (10)
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)
105 case (9)
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)
111 case (8)
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)
117 case (7)
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)
123 case (6)
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)
129 case (5)
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)
135 case (4)
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)
141 case (3)
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)
147 case (2)
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)
153 case default
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)
159 end select
160 !$omp end parallel
161 end associate
162
163 end subroutine opr_cpu_conv1_many
164
165 subroutine cpu_conv1_lx(du, u, vx, vy, vz, dx, dy, dz, &
166 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
167 jacinv, nelv, lx)
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
179 real(kind=rp) :: tmp
180 integer :: e, i, j, k, l
181 !$omp do
182 do e = 1, nelv
183 do j = 1, lx * lx
184 do i = 1, lx
185 tmp = 0.0_rp
186 do k = 1, lx
187 tmp = tmp + dx(i,k) * u(k,j,1,e)
188 end do
189 dudr(i,j,1) = tmp
190 end do
191 end do
192
193 do k = 1, lx
194 do j = 1, lx
195 do i = 1, lx
196 tmp = 0.0_rp
197 do l = 1, lx
198 tmp = tmp + dy(j,l) * u(i,l,k,e)
199 end do
200 duds(i,j,k) = tmp
201 end do
202 end do
203 end do
204
205 do k = 1, lx
206 do i = 1, lx*lx
207 tmp = 0.0_rp
208 do l = 1, lx
209 tmp = tmp + dz(k,l) * u(i,1,l,e)
210 end do
211 dudt(i,1,k) = tmp
212 end do
213 end do
214
215 do i = 1, lx * lx * lx
216 du(i,1,1,e) = jacinv(i,1,1,e) &
217 * ( vx(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) ) &
221 + vy(i,1,1,e) &
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) ) &
225 + vz(i,1,1,e) &
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) ) )
229 end do
230 end do
231 !$omp end do
232 end subroutine cpu_conv1_lx
233
234 subroutine cpu_conv1_lx14(du, u, vx, vy, vz, dx, dy, dz, &
235 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
236 jacinv, nelv)
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
250 !$omp do
251 do e = 1, nelv
252 do j = 1, lx * lx
253 do i = 1, lx
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)
268 end do
269 end do
270
271 do k = 1, lx
272 do j = 1, lx
273 do i = 1, lx
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)
288 end do
289 end do
290 end do
291
292 do k = 1, lx
293 do i = 1, lx*lx
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)
308 end do
309 end do
310
311 do i = 1, lx * lx * lx
312 du(i,1,1,e) = jacinv(i,1,1,e) &
313 * ( vx(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) ) &
317 + vy(i,1,1,e) &
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) ) &
321 + vz(i,1,1,e) &
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) ) )
325 end do
326 end do
327 !$omp end do
328 end subroutine cpu_conv1_lx14
329
330 subroutine cpu_conv1_lx13(du, u, vx, vy, vz, dx, dy, dz, &
331 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
332 jacinv, nelv)
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
346 !$omp do
347 do e = 1, nelv
348 do j = 1, lx * lx
349 do i = 1, lx
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)
363 end do
364 end do
365
366 do k = 1, lx
367 do j = 1, lx
368 do i = 1, lx
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)
382 end do
383 end do
384 end do
385
386 do k = 1, lx
387 do i = 1, lx*lx
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)
401 end do
402 end do
403
404 do i = 1, lx * lx * lx
405 du(i,1,1,e) = jacinv(i,1,1,e) &
406 * ( vx(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) ) &
410 + vy(i,1,1,e) &
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) ) &
414 + vz(i,1,1,e) &
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) ) )
418 end do
419 end do
420 !$omp end do
421 end subroutine cpu_conv1_lx13
422
423 subroutine cpu_conv1_lx12(du, u, vx, vy, vz, dx, dy, dz, &
424 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
425 jacinv, nelv)
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
439 !$omp do
440 do e = 1, nelv
441 do j = 1, lx * lx
442 do i = 1, lx
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)
455 end do
456 end do
457
458 do k = 1, lx
459 do j = 1, lx
460 do i = 1, lx
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)
473 end do
474 end do
475 end do
476
477 do k = 1, lx
478 do i = 1, lx*lx
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)
491 end do
492 end do
493
494 do i = 1, lx * lx * lx
495 du(i,1,1,e) = jacinv(i,1,1,e) &
496 * ( vx(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) ) &
500 + vy(i,1,1,e) &
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) ) &
504 + vz(i,1,1,e) &
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) ) )
508 end do
509 end do
510 !$omp end do
511 end subroutine cpu_conv1_lx12
512
513 subroutine cpu_conv1_lx11(du, u, vx, vy, vz, dx, dy, dz, &
514 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
515 jacinv, nelv)
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
529 !$omp do
530 do e = 1, nelv
531 do j = 1, lx * lx
532 do i = 1, lx
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)
544 end do
545 end do
546
547 do k = 1, lx
548 do j = 1, lx
549 do i = 1, lx
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)
561 end do
562 end do
563 end do
564
565 do k = 1, lx
566 do i = 1, lx*lx
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)
578 end do
579 end do
580
581 do i = 1, lx * lx * lx
582 du(i,1,1,e) = jacinv(i,1,1,e) &
583 * ( vx(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) ) &
587 + vy(i,1,1,e) &
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) ) &
591 + vz(i,1,1,e) &
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) ) )
595 end do
596 end do
597 !$omp end do
598 end subroutine cpu_conv1_lx11
599
600 subroutine cpu_conv1_lx10(du, u, vx, vy, vz, dx, dy, dz, &
601 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
602 jacinv, nelv)
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
616 !$omp do
617 do e = 1, nelv
618 do j = 1, lx * lx
619 do i = 1, lx
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)
630 end do
631 end do
632
633 do k = 1, lx
634 do j = 1, lx
635 do i = 1, lx
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)
646 end do
647 end do
648 end do
649
650 do k = 1, lx
651 do i = 1, lx*lx
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)
662 end do
663 end do
664
665 do i = 1, lx * lx * lx
666 du(i,1,1,e) = jacinv(i,1,1,e) &
667 * ( vx(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) ) &
671 + vy(i,1,1,e) &
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) ) &
675 + vz(i,1,1,e) &
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) ) )
679 end do
680 end do
681 !$omp end do
682 end subroutine cpu_conv1_lx10
683
684 subroutine cpu_conv1_lx9(du, u, vx, vy, vz, dx, dy, dz, &
685 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
686 jacinv, nelv)
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
700 !$omp do
701 do e = 1, nelv
702 do j = 1, lx * lx
703 do i = 1, lx
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)
713 end do
714 end do
715
716 do k = 1, lx
717 do j = 1, lx
718 do i = 1, lx
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)
728 end do
729 end do
730 end do
731
732 do k = 1, lx
733 do i = 1, lx*lx
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)
743 end do
744 end do
745
746 do i = 1, lx * lx * lx
747 du(i,1,1,e) = jacinv(i,1,1,e) &
748 * ( vx(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) ) &
752 + vy(i,1,1,e) &
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) ) &
756 + vz(i,1,1,e) &
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) ) )
760 end do
761 end do
762 !$omp end do
763 end subroutine cpu_conv1_lx9
764
765 subroutine cpu_conv1_lx8(du, u, vx, vy, vz, dx, dy, dz, &
766 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
767 jacinv, nelv)
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
781 !$omp do
782 do e = 1, nelv
783 do j = 1, lx * lx
784 do i = 1, lx
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)
793 end do
794 end do
795
796 do k = 1, lx
797 do j = 1, lx
798 do i = 1, lx
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)
807 end do
808 end do
809 end do
810
811 do k = 1, lx
812 do i = 1, lx*lx
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)
821 end do
822 end do
823
824 do i = 1, lx * lx * lx
825 du(i,1,1,e) = jacinv(i,1,1,e) &
826 * ( vx(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) ) &
830 + vy(i,1,1,e) &
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) ) &
834 + vz(i,1,1,e) &
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) ) )
838 end do
839 end do
840 !$omp end do
841 end subroutine cpu_conv1_lx8
842
843 subroutine cpu_conv1_lx7(du, u, vx, vy, vz, dx, dy, dz, &
844 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
845 jacinv, nelv)
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
859 !$omp do
860 do e = 1, nelv
861 do j = 1, lx * lx
862 do i = 1, lx
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)
870 end do
871 end do
872
873 do k = 1, lx
874 do j = 1, lx
875 do i = 1, lx
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)
883 end do
884 end do
885 end do
886
887 do k = 1, lx
888 do i = 1, lx*lx
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)
896 end do
897 end do
898
899 do i = 1, lx * lx * lx
900 du(i,1,1,e) = jacinv(i,1,1,e) &
901 * ( vx(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) ) &
905 + vy(i,1,1,e) &
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) ) &
909 + vz(i,1,1,e) &
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) ) )
913 end do
914 end do
915 !$omp end do
916 end subroutine cpu_conv1_lx7
917
918 subroutine cpu_conv1_lx6(du, u, vx, vy, vz, dx, dy, dz, &
919 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
920 jacinv, nelv)
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
934 !$omp do
935 do e = 1, nelv
936 do j = 1, lx * lx
937 do i = 1, lx
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)
944 end do
945 end do
946
947 do k = 1, lx
948 do j = 1, lx
949 do i = 1, lx
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)
956 end do
957 end do
958 end do
959
960 do k = 1, lx
961 do i = 1, lx*lx
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)
968 end do
969 end do
970
971 do i = 1, lx * lx * lx
972 du(i,1,1,e) = jacinv(i,1,1,e) &
973 * ( vx(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) ) &
977 + vy(i,1,1,e) &
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) ) &
981 + vz(i,1,1,e) &
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) ) )
985 end do
986 end do
987 !$omp end do
988 end subroutine cpu_conv1_lx6
989
990 subroutine cpu_conv1_lx5(du, u, vx, vy, vz, dx, dy, dz, &
991 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
992 jacinv, nelv)
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
1006 !$omp do
1007 do e = 1, nelv
1008 do j = 1, lx * lx
1009 do i = 1, lx
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)
1015 end do
1016 end do
1017
1018 do k = 1, lx
1019 do j = 1, lx
1020 do i = 1, lx
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)
1026 end do
1027 end do
1028 end do
1029
1030 do k = 1, lx
1031 do i = 1, lx*lx
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)
1037 end do
1038 end do
1039
1040 do i = 1, lx * lx * lx
1041 du(i,1,1,e) = jacinv(i,1,1,e) &
1042 * ( vx(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) ) &
1046 + vy(i,1,1,e) &
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) ) &
1050 + vz(i,1,1,e) &
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) ) )
1054 end do
1055 end do
1056 !$omp end do
1057 end subroutine cpu_conv1_lx5
1058
1059 subroutine cpu_conv1_lx4(du, u, vx, vy, vz, dx, dy, dz, &
1060 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1061 jacinv, nelv)
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
1075 !$omp do
1076 do e = 1, nelv
1077 do j = 1, lx * lx
1078 do i = 1, lx
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)
1083 end do
1084 end do
1085
1086 do k = 1, lx
1087 do j = 1, lx
1088 do i = 1, lx
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)
1093 end do
1094 end do
1095 end do
1096
1097 do k = 1, lx
1098 do i = 1, lx*lx
1099 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
1100 + dz(k,2) * u(i,1,2,e) &
1101 + dz(k,3) * u(i,1,3,e) &
1102 + dz(k,4) * u(i,1,4,e)
1103 end do
1104 end do
1105
1106 do i = 1, lx * lx * lx
1107 du(i,1,1,e) = jacinv(i,1,1,e) &
1108 * ( vx(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) ) &
1112 + vy(i,1,1,e) &
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) ) &
1116 + vz(i,1,1,e) &
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) ) )
1120 end do
1121 end do
1122 !$omp end do
1123 end subroutine cpu_conv1_lx4
1124
1125 subroutine cpu_conv1_lx3(du, u, vx, vy, vz, dx, dy, dz, &
1126 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1127 jacinv, nelv)
1128 integer, parameter :: lx = 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
1141 !$omp do
1142 do e = 1, nelv
1143 do j = 1, lx * lx
1144 do i = 1, lx
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)
1148 end do
1149 end do
1150
1151 do k = 1, lx
1152 do j = 1, lx
1153 do i = 1, lx
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)
1157 end do
1158 end do
1159 end do
1160
1161 do k = 1, lx
1162 do i = 1, lx*lx
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)
1166 end do
1167 end do
1168
1169 do i = 1, lx * lx * lx
1170 du(i,1,1,e) = jacinv(i,1,1,e) &
1171 * ( vx(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) ) &
1175 + vy(i,1,1,e) &
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) ) &
1179 + vz(i,1,1,e) &
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) ) )
1183 end do
1184 end do
1185 !$omp end do
1186 end subroutine cpu_conv1_lx3
1187
1188 subroutine cpu_conv1_lx2(du, u, vx, vy, vz, dx, dy, dz, &
1189 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1190 jacinv, nelv)
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
1204 !$omp do
1205 do e = 1, nelv
1206 do j = 1, lx * lx
1207 do i = 1, lx
1208 dudr(i,j,1) = dx(i,1) * u(1,j,1,e) &
1209 + dx(i,2) * u(2,j,1,e)
1210 end do
1211 end do
1212
1213 do k = 1, lx
1214 do j = 1, lx
1215 do i = 1, lx
1216 duds(i,j,k) = dy(j,1) * u(i,1,k,e) &
1217 + dy(j,2) * u(i,2,k,e)
1218 end do
1219 end do
1220 end do
1221
1222 do k = 1, lx
1223 do i = 1, lx*lx
1224 dudt(i,1,k) = dz(k,1) * u(i,1,1,e) &
1225 + dz(k,2) * u(i,1,2,e)
1226 end do
1227 end do
1228
1229 do i = 1, lx * lx * lx
1230 du(i,1,1,e) = jacinv(i,1,1,e) &
1231 * ( vx(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) ) &
1235 + vy(i,1,1,e) &
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) ) &
1239 + vz(i,1,1,e) &
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) ) )
1243 end do
1244 end do
1245 !$omp end do
1246 end subroutine cpu_conv1_lx2
1247
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)
1258
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)
1263 select case (xh%lx)
1264 case (14)
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),&
1269 jacinv(1,1,1,e))
1270 case (13)
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),&
1275 jacinv(1,1,1,e))
1276 case (12)
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),&
1281 jacinv(1,1,1,e))
1282 case (11)
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),&
1287 jacinv(1,1,1,e))
1288 case (10)
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),&
1293 jacinv(1,1,1,e))
1294 case (9)
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),&
1299 jacinv(1,1,1,e))
1300 case (8)
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),&
1305 jacinv(1,1,1,e))
1306 case (7)
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),&
1311 jacinv(1,1,1,e))
1312 case (6)
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),&
1317 jacinv(1,1,1,e))
1318 case (5)
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),&
1323 jacinv(1,1,1,e))
1324 case (4)
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),&
1329 jacinv(1,1,1,e))
1330 case (3)
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),&
1335 jacinv(1,1,1,e))
1336 case (2)
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),&
1341 jacinv(1,1,1,e))
1342 case default
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)
1348 end select
1349 end associate
1350
1351 end subroutine opr_cpu_conv1_single
1352
1353 !OCL SERIAL
1354 subroutine cpu_conv1_lx_single(du, u, vx, vy, vz, dx, dy, dz, &
1355 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, &
1356 jacinv, lx)
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
1370
1371 do j = 1, lx * lx
1372 do i = 1, lx
1373 tmp = 0.0_rp
1374 do k = 1, lx
1375 tmp = tmp + dx(i,k) * u(k,j,1)
1376 end do
1377 dudr(i,j,1) = tmp
1378 end do
1379 end do
1380
1381 do k = 1, lx
1382 do j = 1, lx
1383 do i = 1, lx
1384 tmp = 0.0_rp
1385 do l = 1, lx
1386 tmp = tmp + dy(j,l) * u(i,l,k)
1387 end do
1388 duds(i,j,k) = tmp
1389 end do
1390 end do
1391 end do
1392
1393 do k = 1, lx
1394 do i = 1, lx*lx
1395 tmp = 0.0_rp
1396 do l = 1, lx
1397 tmp = tmp + dz(k,l) * u(i,1,l)
1398 end do
1399 dudt(i,1,k) = tmp
1400 end do
1401 end do
1402
1403 do i = 1, lx * lx * lx
1404 du(i,1,1) = jacinv(i,1,1) &
1405 * ( vx(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) ) &
1409 + vy(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) ) &
1413 + vz(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) ) )
1417 end do
1418
1419 end subroutine cpu_conv1_lx_single
1420
1421 !OCL SERIAL
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
1435 integer :: i, j, k
1436
1437 do j = 1, lx * lx
1438 do i = 1, lx
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)
1453 end do
1454 end do
1455
1456 do k = 1, lx
1457 do j = 1, lx
1458 do i = 1, lx
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)
1473 end do
1474 end do
1475 end do
1476
1477 do k = 1, lx
1478 do i = 1, lx*lx
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)
1493 end do
1494 end do
1495
1496 do i = 1, lx * lx * lx
1497 du(i,1,1) = jacinv(i,1,1) &
1498 * ( vx(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) ) &
1502 + vy(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) ) &
1506 + vz(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) ) )
1510 end do
1511
1512 end subroutine cpu_conv1_lx14_single
1513
1514 !OCL SERIAL
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
1528 integer :: i, j, k
1529
1530 do j = 1, lx * lx
1531 do i = 1, lx
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)
1545 end do
1546 end do
1547
1548 do k = 1, lx
1549 do j = 1, lx
1550 do i = 1, lx
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)
1564 end do
1565 end do
1566 end do
1567
1568 do k = 1, lx
1569 do i = 1, lx*lx
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)
1583 end do
1584 end do
1585
1586 do i = 1, lx * lx * lx
1587 du(i,1,1) = jacinv(i,1,1) &
1588 * ( vx(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) ) &
1592 + vy(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) ) &
1596 + vz(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) ) )
1600 end do
1601
1602 end subroutine cpu_conv1_lx13_single
1603
1604 !OCL SERIAL
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
1618 integer :: i, j, k
1619
1620 do j = 1, lx * lx
1621 do i = 1, lx
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)
1634 end do
1635 end do
1636
1637 do k = 1, lx
1638 do j = 1, lx
1639 do i = 1, lx
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)
1652 end do
1653 end do
1654 end do
1655
1656 do k = 1, lx
1657 do i = 1, lx*lx
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)
1670 end do
1671 end do
1672
1673 do i = 1, lx * lx * lx
1674 du(i,1,1) = jacinv(i,1,1) &
1675 * ( vx(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) ) &
1679 + vy(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) ) &
1683 + vz(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) ) )
1687 end do
1688
1689 end subroutine cpu_conv1_lx12_single
1690
1691 !OCL SERIAL
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
1705 integer :: i, j, k
1706
1707 do j = 1, lx * lx
1708 do i = 1, lx
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)
1720 end do
1721 end do
1722
1723 do k = 1, lx
1724 do j = 1, lx
1725 do i = 1, lx
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)
1737 end do
1738 end do
1739 end do
1740
1741 do k = 1, lx
1742 do i = 1, lx*lx
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)
1754 end do
1755 end do
1756
1757 do i = 1, lx * lx * lx
1758 du(i,1,1) = jacinv(i,1,1) &
1759 * ( vx(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) ) &
1763 + vy(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) ) &
1767 + vz(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) ) )
1771 end do
1772
1773 end subroutine cpu_conv1_lx11_single
1774
1775 !OCL SERIAL
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
1789 integer :: i, j, k
1790
1791 do j = 1, lx * lx
1792 do i = 1, lx
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)
1803 end do
1804 end do
1805
1806 do k = 1, lx
1807 do j = 1, lx
1808 do i = 1, lx
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)
1819 end do
1820 end do
1821 end do
1822
1823 do k = 1, lx
1824 do i = 1, lx*lx
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)
1835 end do
1836 end do
1837
1838 do i = 1, lx * lx * lx
1839 du(i,1,1) = jacinv(i,1,1) &
1840 * ( vx(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) ) &
1844 + vy(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) ) &
1848 + vz(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) ) )
1852 end do
1853
1854 end subroutine cpu_conv1_lx10_single
1855
1856 !OCL SERIAL
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
1870 integer :: i, j, k
1871
1872 do j = 1, lx * lx
1873 do i = 1, lx
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)
1883 end do
1884 end do
1885
1886 do k = 1, lx
1887 do j = 1, lx
1888 do i = 1, lx
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)
1898 end do
1899 end do
1900 end do
1901
1902 do k = 1, lx
1903 do i = 1, lx*lx
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)
1913 end do
1914 end do
1915
1916 do i = 1, lx * lx * lx
1917 du(i,1,1) = jacinv(i,1,1) &
1918 * ( vx(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) ) &
1922 + vy(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) ) &
1926 + vz(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) ) )
1930 end do
1931
1932 end subroutine cpu_conv1_lx9_single
1933
1934 !OCL SERIAL
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
1948 integer :: i, j, k
1949
1950 do j = 1, lx * lx
1951 do i = 1, lx
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)
1960 end do
1961 end do
1962
1963 do k = 1, lx
1964 do j = 1, lx
1965 do i = 1, lx
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)
1974 end do
1975 end do
1976 end do
1977
1978 do k = 1, lx
1979 do i = 1, lx*lx
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)
1988 end do
1989 end do
1990
1991 do i = 1, lx * lx * lx
1992 du(i,1,1) = jacinv(i,1,1) &
1993 * ( vx(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) ) &
1997 + vy(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) ) &
2001 + vz(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) ) )
2005 end do
2006
2007 end subroutine cpu_conv1_lx8_single
2008
2009 !OCL SERIAL
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
2023 integer :: i, j, k
2024
2025 do j = 1, lx * lx
2026 do i = 1, lx
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)
2034 end do
2035 end do
2036
2037 do k = 1, lx
2038 do j = 1, lx
2039 do i = 1, lx
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)
2047 end do
2048 end do
2049 end do
2050
2051 do k = 1, lx
2052 do i = 1, lx*lx
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)
2060 end do
2061 end do
2062
2063 do i = 1, lx * lx * lx
2064 du(i,1,1) = jacinv(i,1,1) &
2065 * ( vx(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) ) &
2069 + vy(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) ) &
2073 + vz(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) ) )
2077 end do
2078
2079 end subroutine cpu_conv1_lx7_single
2080
2081 !OCL SERIAL
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
2095 integer :: i, j, k
2096
2097 do j = 1, lx * lx
2098 do i = 1, lx
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)
2105 end do
2106 end do
2107
2108 do k = 1, lx
2109 do j = 1, lx
2110 do i = 1, lx
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)
2117 end do
2118 end do
2119 end do
2120
2121 do k = 1, lx
2122 do i = 1, lx*lx
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)
2129 end do
2130 end do
2131
2132 do i = 1, lx * lx * lx
2133 du(i,1,1) = jacinv(i,1,1) &
2134 * ( vx(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) ) &
2138 + vy(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) ) &
2142 + vz(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) ) )
2146 end do
2147
2148 end subroutine cpu_conv1_lx6_single
2149
2150 !OCL SERIAL
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
2164 integer :: i, j, k
2165
2166 do j = 1, lx * lx
2167 do i = 1, lx
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)
2173 end do
2174 end do
2175
2176 do k = 1, lx
2177 do j = 1, lx
2178 do i = 1, lx
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)
2184 end do
2185 end do
2186 end do
2187
2188 do k = 1, lx
2189 do i = 1, lx*lx
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)
2195 end do
2196 end do
2197
2198 do i = 1, lx * lx * lx
2199 du(i,1,1) = jacinv(i,1,1) &
2200 * ( vx(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) ) &
2204 + vy(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) ) &
2208 + vz(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) ) )
2212 end do
2213
2214 end subroutine cpu_conv1_lx5_single
2215
2216 !OCL SERIAL
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
2230 integer :: i, j, k
2231
2232 do j = 1, lx * lx
2233 do i = 1, lx
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)
2238 end do
2239 end do
2240
2241 do k = 1, lx
2242 do j = 1, lx
2243 do i = 1, lx
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)
2248 end do
2249 end do
2250 end do
2251
2252 do k = 1, lx
2253 do i = 1, lx*lx
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)
2258 end do
2259 end do
2260
2261 do i = 1, lx * lx * lx
2262 du(i,1,1) = jacinv(i,1,1) &
2263 * ( vx(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) ) &
2267 + vy(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) ) &
2271 + vz(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) ) )
2275 end do
2276
2277 end subroutine cpu_conv1_lx4_single
2278
2279 !OCL SERIAL
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
2293 integer :: i, j, k
2294
2295 do j = 1, lx * lx
2296 do i = 1, lx
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)
2300 end do
2301 end do
2302
2303 do k = 1, lx
2304 do j = 1, lx
2305 do i = 1, lx
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)
2309 end do
2310 end do
2311 end do
2312
2313 do k = 1, lx
2314 do i = 1, lx*lx
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)
2318 end do
2319 end do
2320
2321 do i = 1, lx * lx * lx
2322 du(i,1,1) = jacinv(i,1,1) &
2323 * ( vx(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) ) &
2327 + vy(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) ) &
2331 + vz(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) ) )
2335 end do
2336
2337 end subroutine cpu_conv1_lx3_single
2338
2339 !OCL SERIAL
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
2353 integer :: i, j, k
2354
2355 do j = 1, lx * lx
2356 do i = 1, lx
2357 dudr(i,j,1) = dx(i,1) * u(1,j,1) &
2358 + dx(i,2) * u(2,j,1)
2359 end do
2360 end do
2361
2362 do k = 1, lx
2363 do j = 1, lx
2364 do i = 1, lx
2365 duds(i,j,k) = dy(j,1) * u(i,1,k) &
2366 + dy(j,2) * u(i,2,k)
2367 end do
2368 end do
2369 end do
2370
2371 do k = 1, lx
2372 do i = 1, lx*lx
2373 dudt(i,1,k) = dz(k,1) * u(i,1,1) &
2374 + dz(k,2) * u(i,1,2)
2375 end do
2376 end do
2377
2378 do i = 1, lx * lx * lx
2379 du(i,1,1) = jacinv(i,1,1) &
2380 * ( vx(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) ) &
2384 + vy(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) ) &
2388 + vz(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) ) )
2392 end do
2393
2394 end subroutine cpu_conv1_lx2_single
2395
2396end submodule cpu_conv1
Operators CPU backend.
Definition opr_cpu.f90:34