Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
cpu_opgrad.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_opgrad
35 implicit none
36
37contains
38
39 module subroutine opr_cpu_opgrad(ux, uy, uz, u, coef, e_start, e_end)
40 type(coef_t), intent(in) :: coef
41 integer, intent(in) :: e_start, e_end
42 real(kind=rp), dimension(coef%Xh%lxyz, e_end - e_start + 1), &
43 intent(inout) :: ux
44 real(kind=rp), dimension(coef%Xh%lxyz, e_end - e_start + 1), &
45 intent(inout) :: uy
46 real(kind=rp), dimension(coef%Xh%lxyz, e_end - e_start + 1), &
47 intent(inout) :: uz
48 real(kind=rp), dimension(coef%Xh%lxyz, e_end - e_start + 1), &
49 intent(in) :: u
50 integer :: e_len
51
52 e_len = e_end-e_start+1
53
54 if (e_len .eq. 1) then
55 call opr_cpu_opgrad_single(ux, uy, uz, u, coef, e_start)
56 else
57 call opr_cpu_opgrad_many(ux, uy, uz, u, coef, e_start, e_len)
58 end if
59 end subroutine opr_cpu_opgrad
60
61 subroutine opr_cpu_opgrad_many(ux, uy, uz, u, coef, e_start, e_len)
62 type(coef_t), intent(in) :: coef
63 integer, intent(in) :: e_start, e_len
64 real(kind=rp), dimension(coef%Xh%lxyz, e_len), intent(inout) :: ux
65 real(kind=rp), dimension(coef%Xh%lxyz, e_len), intent(inout) :: uy
66 real(kind=rp), dimension(coef%Xh%lxyz, e_len), intent(inout) :: uz
67 real(kind=rp), dimension(coef%Xh%lxyz, e_len), intent(in) :: u
68
69 associate(xh => coef%Xh, msh => coef%msh, &
70 drdx => coef%drdx, drdy => coef%drdy, drdz => coef%drdz, &
71 dsdx => coef%dsdx, dsdy => coef%dsdy, dsdz => coef%dsdz, &
72 dtdx => coef%dtdx, dtdy => coef%dtdy, dtdz => coef%dtdz)
73
74 !$omp parallel
75 select case (xh%lx)
76 case (18)
77 call cpu_opgrad_lx18(ux, uy, uz, u, &
78 xh%dx, xh%dy, xh%dz, &
79 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
80 dtdx(1, 1, 1, e_start), &
81 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
82 dtdy(1, 1, 1, e_start), &
83 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
84 dtdz(1, 1, 1, e_start), &
85 xh%w3, e_len)
86 case (17)
87 call cpu_opgrad_lx17(ux, uy, uz, u, &
88 xh%dx, xh%dy, xh%dz, &
89 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
90 dtdx(1, 1, 1, e_start), &
91 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
92 dtdy(1, 1, 1, e_start), &
93 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
94 dtdz(1, 1, 1, e_start), &
95 xh%w3, e_len)
96 case (16)
97 call cpu_opgrad_lx16(ux, uy, uz, u, &
98 xh%dx, xh%dy, xh%dz, &
99 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
100 dtdx(1, 1, 1, e_start), &
101 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
102 dtdy(1, 1, 1, e_start), &
103 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
104 dtdz(1, 1, 1, e_start), &
105 xh%w3, e_len)
106 case (15)
107 call cpu_opgrad_lx15(ux, uy, uz, u, &
108 xh%dx, xh%dy, xh%dz, &
109 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
110 dtdx(1, 1, 1, e_start), &
111 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
112 dtdy(1, 1, 1, e_start), &
113 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
114 dtdz(1, 1, 1, e_start), &
115 xh%w3, e_len)
116 case (14)
117 call cpu_opgrad_lx14(ux, uy, uz, u, &
118 xh%dx, xh%dy, xh%dz, &
119 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
120 dtdx(1, 1, 1, e_start), &
121 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
122 dtdy(1, 1, 1, e_start), &
123 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
124 dtdz(1, 1, 1, e_start), &
125 xh%w3, e_len)
126 case (13)
127 call cpu_opgrad_lx13(ux, uy, uz, u, &
128 xh%dx, xh%dy, xh%dz, &
129 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
130 dtdx(1, 1, 1, e_start), &
131 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
132 dtdy(1, 1, 1, e_start), &
133 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
134 dtdz(1, 1, 1, e_start), &
135 xh%w3, e_len)
136 case (12)
137 call cpu_opgrad_lx12(ux, uy, uz, u, &
138 xh%dx, xh%dy, xh%dz, &
139 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
140 dtdx(1, 1, 1, e_start), &
141 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
142 dtdy(1, 1, 1, e_start), &
143 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
144 dtdz(1, 1, 1, e_start), &
145 xh%w3, e_len)
146 case (11)
147 call cpu_opgrad_lx11(ux, uy, uz, u, &
148 xh%dx, xh%dy, xh%dz, &
149 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
150 dtdx(1, 1, 1, e_start), &
151 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
152 dtdy(1, 1, 1, e_start), &
153 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
154 dtdz(1, 1, 1, e_start), &
155 xh%w3, e_len)
156 case (10)
157 call cpu_opgrad_lx10(ux, uy, uz, u, &
158 xh%dx, xh%dy, xh%dz, &
159 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
160 dtdx(1, 1, 1, e_start), &
161 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
162 dtdy(1, 1, 1, e_start), &
163 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
164 dtdz(1, 1, 1, e_start), &
165 xh%w3, e_len)
166
167 case (9)
168 call cpu_opgrad_lx9(ux, uy, uz, u, &
169 xh%dx, xh%dy, xh%dz, &
170 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
171 dtdx(1, 1, 1, e_start), &
172 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
173 dtdy(1, 1, 1, e_start), &
174 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
175 dtdz(1, 1, 1, e_start), &
176 xh%w3, e_len)
177 case (8)
178 call cpu_opgrad_lx8(ux, uy, uz, u, &
179 xh%dx, xh%dy, xh%dz, &
180 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
181 dtdx(1, 1, 1, e_start), &
182 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
183 dtdy(1, 1, 1, e_start), &
184 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
185 dtdz(1, 1, 1, e_start), &
186 xh%w3, e_len)
187 case (7)
188 call cpu_opgrad_lx7(ux, uy, uz, u, &
189 xh%dx, xh%dy, xh%dz, &
190 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
191 dtdx(1, 1, 1, e_start), &
192 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
193 dtdy(1, 1, 1, e_start), &
194 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
195 dtdz(1, 1, 1, e_start), &
196 xh%w3, e_len)
197 case (6)
198 call cpu_opgrad_lx6(ux, uy, uz, u, &
199 xh%dx, xh%dy, xh%dz, &
200 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
201 dtdx(1, 1, 1, e_start), &
202 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
203 dtdy(1, 1, 1, e_start), &
204 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
205 dtdz(1, 1, 1, e_start), &
206 xh%w3, e_len)
207 case (5)
208 call cpu_opgrad_lx5(ux, uy, uz, u, &
209 xh%dx, xh%dy, xh%dz, &
210 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
211 dtdx(1, 1, 1, e_start), &
212 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
213 dtdy(1, 1, 1, e_start), &
214 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
215 dtdz(1, 1, 1, e_start), &
216 xh%w3, e_len)
217 case (4)
218 call cpu_opgrad_lx4(ux, uy, uz, u, &
219 xh%dx, xh%dy, xh%dz, &
220 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
221 dtdx(1, 1, 1, e_start), &
222 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
223 dtdy(1, 1, 1, e_start), &
224 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
225 dtdz(1, 1, 1, e_start), &
226 xh%w3, e_len)
227 case (3)
228 call cpu_opgrad_lx3(ux, uy, uz, u, &
229 xh%dx, xh%dy, xh%dz, &
230 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
231 dtdx(1, 1, 1, e_start), &
232 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
233 dtdy(1, 1, 1, e_start), &
234 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
235 dtdz(1, 1, 1, e_start), &
236 xh%w3, e_len)
237 case (2)
238 call cpu_opgrad_lx2(ux, uy, uz, u, &
239 xh%dx, xh%dy, xh%dz, &
240 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
241 dtdx(1, 1, 1, e_start), &
242 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
243 dtdy(1, 1, 1, e_start), &
244 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
245 dtdz(1, 1, 1, e_start), &
246 xh%w3, e_len)
247 case default
248 call cpu_opgrad_lx(ux, uy, uz, u, &
249 xh%dx, xh%dy, xh%dz, &
250 drdx(1, 1, 1, e_start), dsdx(1, 1, 1, e_start), &
251 dtdx(1, 1, 1, e_start), &
252 drdy(1, 1, 1, e_start), dsdy(1, 1, 1, e_start), &
253 dtdy(1, 1, 1, e_start), &
254 drdz(1, 1, 1, e_start), dsdz(1, 1, 1, e_start), &
255 dtdz(1, 1, 1, e_start), &
256 xh%w3, e_len, xh%lx)
257 end select
258 !$omp end parallel
259 end associate
260
261 end subroutine opr_cpu_opgrad_many
262
263 subroutine cpu_opgrad_lx(ux, uy, uz, u, dx, dy, dz, &
264 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n, lx)
265 integer, intent(in) :: n, lx
266 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
267 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
268 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
269 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
270 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
271 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
272 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
273 real(kind=rp) :: ur(lx, lx, lx)
274 real(kind=rp) :: us(lx, lx, lx)
275 real(kind=rp) :: ut(lx, lx, lx)
276 real(kind=rp) :: tmp
277 integer :: e, i, j, k, l
278
279 !$omp do
280 do e = 1, n
281 do j = 1, lx * lx
282 do i = 1, lx
283 tmp = 0.0_rp
284 do k = 1, lx
285 tmp = tmp + dx(i,k) * u(k,j,1,e)
286 end do
287 ur(i,j,1) = tmp
288 end do
289 end do
290
291 do k = 1, lx
292 do j = 1, lx
293 do i = 1, lx
294 tmp = 0.0_rp
295 do l = 1, lx
296 tmp = tmp + dy(j,l) * u(i,l,k,e)
297 end do
298 us(i,j,k) = tmp
299 end do
300 end do
301 end do
302
303 do k = 1, lx
304 do i = 1, lx*lx
305 tmp = 0.0_rp
306 do l = 1, lx
307 tmp = tmp + dz(k,l) * u(i,1,l,e)
308 end do
309 ut(i,1,k) = tmp
310 end do
311 end do
312
313 do i = 1, lx * lx * lx
314 ux(i,1,1,e) = w3(i,1,1) &
315 * ( drdx(i,1,1,e) * ur(i,1,1) &
316 + dsdx(i,1,1,e) * us(i,1,1) &
317 + dtdx(i,1,1,e) * ut(i,1,1) )
318 uy(i,1,1,e) = w3(i,1,1) &
319 * ( dsdy(i,1,1,e) * us(i,1,1) &
320 + drdy(i,1,1,e) * ur(i,1,1) &
321 + dtdy(i,1,1,e) * ut(i,1,1) )
322 uz(i,1,1,e) = w3(i,1,1) &
323 * ( dtdz(i,1,1,e) * ut(i,1,1) &
324 + drdz(i,1,1,e) * ur(i,1,1) &
325 + dsdz(i,1,1,e) * us(i,1,1) )
326 end do
327 end do
328 !$omp end do
329 end subroutine cpu_opgrad_lx
330
331 subroutine cpu_opgrad_lx18(ux, uy, uz, u, dx, dy, dz, &
332 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
333 integer, parameter :: lx = 18
334 integer, intent(in) :: n
335 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
336 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
337 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
338 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
339 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
340 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
341 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
342 real(kind=rp) :: ur(lx, lx, lx)
343 real(kind=rp) :: us(lx, lx, lx)
344 real(kind=rp) :: ut(lx, lx, lx)
345 integer :: e, i, j, k
346
347 !$omp do
348 do e = 1, n
349 do j = 1, lx * lx
350 do i = 1, lx
351 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
352 + dx(i,2) * u(2,j,1,e) &
353 + dx(i,3) * u(3,j,1,e) &
354 + dx(i,4) * u(4,j,1,e) &
355 + dx(i,5) * u(5,j,1,e) &
356 + dx(i,6) * u(6,j,1,e) &
357 + dx(i,7) * u(7,j,1,e) &
358 + dx(i,8) * u(8,j,1,e) &
359 + dx(i,9) * u(9,j,1,e) &
360 + dx(i,10) * u(10,j,1,e) &
361 + dx(i,11) * u(11,j,1,e) &
362 + dx(i,12) * u(12,j,1,e) &
363 + dx(i,13) * u(13,j,1,e) &
364 + dx(i,14) * u(14,j,1,e) &
365 + dx(i,15) * u(15,j,1,e) &
366 + dx(i,16) * u(16,j,1,e) &
367 + dx(i,17) * u(17,j,1,e) &
368 + dx(i,18) * u(18,j,1,e)
369 end do
370 end do
371
372 do k = 1, lx
373 do j = 1, lx
374 do i = 1, lx
375 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
376 + dy(j,2) * u(i,2,k,e) &
377 + dy(j,3) * u(i,3,k,e) &
378 + dy(j,4) * u(i,4,k,e) &
379 + dy(j,5) * u(i,5,k,e) &
380 + dy(j,6) * u(i,6,k,e) &
381 + dy(j,7) * u(i,7,k,e) &
382 + dy(j,8) * u(i,8,k,e) &
383 + dy(j,9) * u(i,9,k,e) &
384 + dy(j,10) * u(i,10,k,e) &
385 + dy(j,11) * u(i,11,k,e) &
386 + dy(j,12) * u(i,12,k,e) &
387 + dy(j,13) * u(i,13,k,e) &
388 + dy(j,14) * u(i,14,k,e) &
389 + dy(j,15) * u(i,15,k,e) &
390 + dy(j,16) * u(i,16,k,e) &
391 + dy(j,17) * u(i,17,k,e) &
392 + dy(j,18) * u(i,18,k,e)
393 end do
394 end do
395 end do
396
397 do k = 1, lx
398 do i = 1, lx*lx
399 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
400 + dz(k,2) * u(i,1,2,e) &
401 + dz(k,3) * u(i,1,3,e) &
402 + dz(k,4) * u(i,1,4,e) &
403 + dz(k,5) * u(i,1,5,e) &
404 + dz(k,6) * u(i,1,6,e) &
405 + dz(k,7) * u(i,1,7,e) &
406 + dz(k,8) * u(i,1,8,e) &
407 + dz(k,9) * u(i,1,9,e) &
408 + dz(k,10) * u(i,1,10,e) &
409 + dz(k,11) * u(i,1,11,e) &
410 + dz(k,12) * u(i,1,12,e) &
411 + dz(k,13) * u(i,1,13,e) &
412 + dz(k,14) * u(i,1,14,e) &
413 + dz(k,15) * u(i,1,15,e) &
414 + dz(k,16) * u(i,1,16,e) &
415 + dz(k,17) * u(i,1,17,e) &
416 + dz(k,18) * u(i,1,18,e)
417 end do
418 end do
419
420 do i = 1, lx * lx * lx
421 ux(i,1,1,e) = w3(i,1,1) &
422 * ( drdx(i,1,1,e) * ur(i,1,1) &
423 + dsdx(i,1,1,e) * us(i,1,1) &
424 + dtdx(i,1,1,e) * ut(i,1,1) )
425 uy(i,1,1,e) = w3(i,1,1) &
426 * ( dsdy(i,1,1,e) * us(i,1,1) &
427 + drdy(i,1,1,e) * ur(i,1,1) &
428 + dtdy(i,1,1,e) * ut(i,1,1) )
429 uz(i,1,1,e) = w3(i,1,1) &
430 * ( dtdz(i,1,1,e) * ut(i,1,1) &
431 + drdz(i,1,1,e) * ur(i,1,1) &
432 + dsdz(i,1,1,e) * us(i,1,1) )
433 end do
434 end do
435 !$omp end do
436 end subroutine cpu_opgrad_lx18
437
438 subroutine cpu_opgrad_lx17(ux, uy, uz, u, dx, dy, dz, &
439 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
440 integer, parameter :: lx = 17
441 integer, intent(in) :: n
442 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
443 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
444 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
445 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
446 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
447 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
448 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
449 real(kind=rp) :: ur(lx, lx, lx)
450 real(kind=rp) :: us(lx, lx, lx)
451 real(kind=rp) :: ut(lx, lx, lx)
452 integer :: e, i, j, k
453
454 !$omp do
455 do e = 1, n
456 do j = 1, lx * lx
457 do i = 1, lx
458 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
459 + dx(i,2) * u(2,j,1,e) &
460 + dx(i,3) * u(3,j,1,e) &
461 + dx(i,4) * u(4,j,1,e) &
462 + dx(i,5) * u(5,j,1,e) &
463 + dx(i,6) * u(6,j,1,e) &
464 + dx(i,7) * u(7,j,1,e) &
465 + dx(i,8) * u(8,j,1,e) &
466 + dx(i,9) * u(9,j,1,e) &
467 + dx(i,10) * u(10,j,1,e) &
468 + dx(i,11) * u(11,j,1,e) &
469 + dx(i,12) * u(12,j,1,e) &
470 + dx(i,13) * u(13,j,1,e) &
471 + dx(i,14) * u(14,j,1,e) &
472 + dx(i,15) * u(15,j,1,e) &
473 + dx(i,16) * u(16,j,1,e) &
474 + dx(i,17) * u(17,j,1,e)
475 end do
476 end do
477
478 do k = 1, lx
479 do j = 1, lx
480 do i = 1, lx
481 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
482 + dy(j,2) * u(i,2,k,e) &
483 + dy(j,3) * u(i,3,k,e) &
484 + dy(j,4) * u(i,4,k,e) &
485 + dy(j,5) * u(i,5,k,e) &
486 + dy(j,6) * u(i,6,k,e) &
487 + dy(j,7) * u(i,7,k,e) &
488 + dy(j,8) * u(i,8,k,e) &
489 + dy(j,9) * u(i,9,k,e) &
490 + dy(j,10) * u(i,10,k,e) &
491 + dy(j,11) * u(i,11,k,e) &
492 + dy(j,12) * u(i,12,k,e) &
493 + dy(j,13) * u(i,13,k,e) &
494 + dy(j,14) * u(i,14,k,e) &
495 + dy(j,15) * u(i,15,k,e) &
496 + dy(j,16) * u(i,16,k,e) &
497 + dy(j,17) * u(i,17,k,e)
498 end do
499 end do
500 end do
501
502 do k = 1, lx
503 do i = 1, lx*lx
504 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
505 + dz(k,2) * u(i,1,2,e) &
506 + dz(k,3) * u(i,1,3,e) &
507 + dz(k,4) * u(i,1,4,e) &
508 + dz(k,5) * u(i,1,5,e) &
509 + dz(k,6) * u(i,1,6,e) &
510 + dz(k,7) * u(i,1,7,e) &
511 + dz(k,8) * u(i,1,8,e) &
512 + dz(k,9) * u(i,1,9,e) &
513 + dz(k,10) * u(i,1,10,e) &
514 + dz(k,11) * u(i,1,11,e) &
515 + dz(k,12) * u(i,1,12,e) &
516 + dz(k,13) * u(i,1,13,e) &
517 + dz(k,14) * u(i,1,14,e) &
518 + dz(k,15) * u(i,1,15,e) &
519 + dz(k,16) * u(i,1,16,e) &
520 + dz(k,17) * u(i,1,17,e)
521 end do
522 end do
523
524 do i = 1, lx * lx * lx
525 ux(i,1,1,e) = w3(i,1,1) &
526 * ( drdx(i,1,1,e) * ur(i,1,1) &
527 + dsdx(i,1,1,e) * us(i,1,1) &
528 + dtdx(i,1,1,e) * ut(i,1,1) )
529 uy(i,1,1,e) = w3(i,1,1) &
530 * ( dsdy(i,1,1,e) * us(i,1,1) &
531 + drdy(i,1,1,e) * ur(i,1,1) &
532 + dtdy(i,1,1,e) * ut(i,1,1) )
533 uz(i,1,1,e) = w3(i,1,1) &
534 * ( dtdz(i,1,1,e) * ut(i,1,1) &
535 + drdz(i,1,1,e) * ur(i,1,1) &
536 + dsdz(i,1,1,e) * us(i,1,1) )
537 end do
538 end do
539 !$omp end do
540 end subroutine cpu_opgrad_lx17
541
542 subroutine cpu_opgrad_lx16(ux, uy, uz, u, dx, dy, dz, &
543 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
544 integer, parameter :: lx = 16
545 integer, intent(in) :: n
546 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
547 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
548 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
549 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
550 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
551 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
552 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
553 real(kind=rp) :: ur(lx, lx, lx)
554 real(kind=rp) :: us(lx, lx, lx)
555 real(kind=rp) :: ut(lx, lx, lx)
556 integer :: e, i, j, k
557
558 !$omp do
559 do e = 1, n
560 do j = 1, lx * lx
561 do i = 1, lx
562 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
563 + dx(i,2) * u(2,j,1,e) &
564 + dx(i,3) * u(3,j,1,e) &
565 + dx(i,4) * u(4,j,1,e) &
566 + dx(i,5) * u(5,j,1,e) &
567 + dx(i,6) * u(6,j,1,e) &
568 + dx(i,7) * u(7,j,1,e) &
569 + dx(i,8) * u(8,j,1,e) &
570 + dx(i,9) * u(9,j,1,e) &
571 + dx(i,10) * u(10,j,1,e) &
572 + dx(i,11) * u(11,j,1,e) &
573 + dx(i,12) * u(12,j,1,e) &
574 + dx(i,13) * u(13,j,1,e) &
575 + dx(i,14) * u(14,j,1,e) &
576 + dx(i,15) * u(15,j,1,e) &
577 + dx(i,16) * u(16,j,1,e)
578 end do
579 end do
580
581 do k = 1, lx
582 do j = 1, lx
583 do i = 1, lx
584 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
585 + dy(j,2) * u(i,2,k,e) &
586 + dy(j,3) * u(i,3,k,e) &
587 + dy(j,4) * u(i,4,k,e) &
588 + dy(j,5) * u(i,5,k,e) &
589 + dy(j,6) * u(i,6,k,e) &
590 + dy(j,7) * u(i,7,k,e) &
591 + dy(j,8) * u(i,8,k,e) &
592 + dy(j,9) * u(i,9,k,e) &
593 + dy(j,10) * u(i,10,k,e) &
594 + dy(j,11) * u(i,11,k,e) &
595 + dy(j,12) * u(i,12,k,e) &
596 + dy(j,13) * u(i,13,k,e) &
597 + dy(j,14) * u(i,14,k,e) &
598 + dy(j,15) * u(i,15,k,e) &
599 + dy(j,16) * u(i,16,k,e)
600 end do
601 end do
602 end do
603
604 do k = 1, lx
605 do i = 1, lx*lx
606 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
607 + dz(k,2) * u(i,1,2,e) &
608 + dz(k,3) * u(i,1,3,e) &
609 + dz(k,4) * u(i,1,4,e) &
610 + dz(k,5) * u(i,1,5,e) &
611 + dz(k,6) * u(i,1,6,e) &
612 + dz(k,7) * u(i,1,7,e) &
613 + dz(k,8) * u(i,1,8,e) &
614 + dz(k,9) * u(i,1,9,e) &
615 + dz(k,10) * u(i,1,10,e) &
616 + dz(k,11) * u(i,1,11,e) &
617 + dz(k,12) * u(i,1,12,e) &
618 + dz(k,13) * u(i,1,13,e) &
619 + dz(k,14) * u(i,1,14,e) &
620 + dz(k,15) * u(i,1,15,e) &
621 + dz(k,16) * u(i,1,16,e)
622 end do
623 end do
624
625 do i = 1, lx * lx * lx
626 ux(i,1,1,e) = w3(i,1,1) &
627 * ( drdx(i,1,1,e) * ur(i,1,1) &
628 + dsdx(i,1,1,e) * us(i,1,1) &
629 + dtdx(i,1,1,e) * ut(i,1,1) )
630 uy(i,1,1,e) = w3(i,1,1) &
631 * ( dsdy(i,1,1,e) * us(i,1,1) &
632 + drdy(i,1,1,e) * ur(i,1,1) &
633 + dtdy(i,1,1,e) * ut(i,1,1) )
634 uz(i,1,1,e) = w3(i,1,1) &
635 * ( dtdz(i,1,1,e) * ut(i,1,1) &
636 + drdz(i,1,1,e) * ur(i,1,1) &
637 + dsdz(i,1,1,e) * us(i,1,1) )
638 end do
639 end do
640 !$omp end do
641 end subroutine cpu_opgrad_lx16
642
643 subroutine cpu_opgrad_lx15(ux, uy, uz, u, dx, dy, dz, &
644 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
645 integer, parameter :: lx = 15
646 integer, intent(in) :: n
647 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
648 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
649 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
650 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
651 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
652 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
653 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
654 real(kind=rp) :: ur(lx, lx, lx)
655 real(kind=rp) :: us(lx, lx, lx)
656 real(kind=rp) :: ut(lx, lx, lx)
657 integer :: e, i, j, k
658
659 !$omp do
660 do e = 1, n
661 do j = 1, lx * lx
662 do i = 1, lx
663 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
664 + dx(i,2) * u(2,j,1,e) &
665 + dx(i,3) * u(3,j,1,e) &
666 + dx(i,4) * u(4,j,1,e) &
667 + dx(i,5) * u(5,j,1,e) &
668 + dx(i,6) * u(6,j,1,e) &
669 + dx(i,7) * u(7,j,1,e) &
670 + dx(i,8) * u(8,j,1,e) &
671 + dx(i,9) * u(9,j,1,e) &
672 + dx(i,10) * u(10,j,1,e) &
673 + dx(i,11) * u(11,j,1,e) &
674 + dx(i,12) * u(12,j,1,e) &
675 + dx(i,13) * u(13,j,1,e) &
676 + dx(i,14) * u(14,j,1,e) &
677 + dx(i,15) * u(15,j,1,e)
678 end do
679 end do
680
681 do k = 1, lx
682 do j = 1, lx
683 do i = 1, lx
684 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
685 + dy(j,2) * u(i,2,k,e) &
686 + dy(j,3) * u(i,3,k,e) &
687 + dy(j,4) * u(i,4,k,e) &
688 + dy(j,5) * u(i,5,k,e) &
689 + dy(j,6) * u(i,6,k,e) &
690 + dy(j,7) * u(i,7,k,e) &
691 + dy(j,8) * u(i,8,k,e) &
692 + dy(j,9) * u(i,9,k,e) &
693 + dy(j,10) * u(i,10,k,e) &
694 + dy(j,11) * u(i,11,k,e) &
695 + dy(j,12) * u(i,12,k,e) &
696 + dy(j,13) * u(i,13,k,e) &
697 + dy(j,14) * u(i,14,k,e) &
698 + dy(j,15) * u(i,15,k,e)
699 end do
700 end do
701 end do
702
703 do k = 1, lx
704 do i = 1, lx*lx
705 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
706 + dz(k,2) * u(i,1,2,e) &
707 + dz(k,3) * u(i,1,3,e) &
708 + dz(k,4) * u(i,1,4,e) &
709 + dz(k,5) * u(i,1,5,e) &
710 + dz(k,6) * u(i,1,6,e) &
711 + dz(k,7) * u(i,1,7,e) &
712 + dz(k,8) * u(i,1,8,e) &
713 + dz(k,9) * u(i,1,9,e) &
714 + dz(k,10) * u(i,1,10,e) &
715 + dz(k,11) * u(i,1,11,e) &
716 + dz(k,12) * u(i,1,12,e) &
717 + dz(k,13) * u(i,1,13,e) &
718 + dz(k,14) * u(i,1,14,e) &
719 + dz(k,15) * u(i,1,15,e)
720 end do
721 end do
722
723 do i = 1, lx * lx * lx
724 ux(i,1,1,e) = w3(i,1,1) &
725 * ( drdx(i,1,1,e) * ur(i,1,1) &
726 + dsdx(i,1,1,e) * us(i,1,1) &
727 + dtdx(i,1,1,e) * ut(i,1,1) )
728 uy(i,1,1,e) = w3(i,1,1) &
729 * ( dsdy(i,1,1,e) * us(i,1,1) &
730 + drdy(i,1,1,e) * ur(i,1,1) &
731 + dtdy(i,1,1,e) * ut(i,1,1) )
732 uz(i,1,1,e) = w3(i,1,1) &
733 * ( dtdz(i,1,1,e) * ut(i,1,1) &
734 + drdz(i,1,1,e) * ur(i,1,1) &
735 + dsdz(i,1,1,e) * us(i,1,1) )
736 end do
737 end do
738 !$omp end do
739 end subroutine cpu_opgrad_lx15
740
741 subroutine cpu_opgrad_lx14(ux, uy, uz, u, dx, dy, dz, &
742 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
743 integer, parameter :: lx = 14
744 integer, intent(in) :: n
745 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
746 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
747 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
748 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
749 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
750 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
751 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
752 real(kind=rp) :: ur(lx, lx, lx)
753 real(kind=rp) :: us(lx, lx, lx)
754 real(kind=rp) :: ut(lx, lx, lx)
755 integer :: e, i, j, k
756
757 !$omp do
758 do e = 1, n
759 do j = 1, lx * lx
760 do i = 1, lx
761 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
762 + dx(i,2) * u(2,j,1,e) &
763 + dx(i,3) * u(3,j,1,e) &
764 + dx(i,4) * u(4,j,1,e) &
765 + dx(i,5) * u(5,j,1,e) &
766 + dx(i,6) * u(6,j,1,e) &
767 + dx(i,7) * u(7,j,1,e) &
768 + dx(i,8) * u(8,j,1,e) &
769 + dx(i,9) * u(9,j,1,e) &
770 + dx(i,10) * u(10,j,1,e) &
771 + dx(i,11) * u(11,j,1,e) &
772 + dx(i,12) * u(12,j,1,e) &
773 + dx(i,13) * u(13,j,1,e) &
774 + dx(i,14) * u(14,j,1,e)
775 end do
776 end do
777
778 do k = 1, lx
779 do j = 1, lx
780 do i = 1, lx
781 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
782 + dy(j,2) * u(i,2,k,e) &
783 + dy(j,3) * u(i,3,k,e) &
784 + dy(j,4) * u(i,4,k,e) &
785 + dy(j,5) * u(i,5,k,e) &
786 + dy(j,6) * u(i,6,k,e) &
787 + dy(j,7) * u(i,7,k,e) &
788 + dy(j,8) * u(i,8,k,e) &
789 + dy(j,9) * u(i,9,k,e) &
790 + dy(j,10) * u(i,10,k,e) &
791 + dy(j,11) * u(i,11,k,e) &
792 + dy(j,12) * u(i,12,k,e) &
793 + dy(j,13) * u(i,13,k,e) &
794 + dy(j,14) * u(i,14,k,e)
795 end do
796 end do
797 end do
798
799 do k = 1, lx
800 do i = 1, lx*lx
801 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
802 + dz(k,2) * u(i,1,2,e) &
803 + dz(k,3) * u(i,1,3,e) &
804 + dz(k,4) * u(i,1,4,e) &
805 + dz(k,5) * u(i,1,5,e) &
806 + dz(k,6) * u(i,1,6,e) &
807 + dz(k,7) * u(i,1,7,e) &
808 + dz(k,8) * u(i,1,8,e) &
809 + dz(k,9) * u(i,1,9,e) &
810 + dz(k,10) * u(i,1,10,e) &
811 + dz(k,11) * u(i,1,11,e) &
812 + dz(k,12) * u(i,1,12,e) &
813 + dz(k,13) * u(i,1,13,e) &
814 + dz(k,14) * u(i,1,14,e)
815 end do
816 end do
817
818 do i = 1, lx * lx * lx
819 ux(i,1,1,e) = w3(i,1,1) &
820 * ( drdx(i,1,1,e) * ur(i,1,1) &
821 + dsdx(i,1,1,e) * us(i,1,1) &
822 + dtdx(i,1,1,e) * ut(i,1,1) )
823 uy(i,1,1,e) = w3(i,1,1) &
824 * ( dsdy(i,1,1,e) * us(i,1,1) &
825 + drdy(i,1,1,e) * ur(i,1,1) &
826 + dtdy(i,1,1,e) * ut(i,1,1) )
827 uz(i,1,1,e) = w3(i,1,1) &
828 * ( dtdz(i,1,1,e) * ut(i,1,1) &
829 + drdz(i,1,1,e) * ur(i,1,1) &
830 + dsdz(i,1,1,e) * us(i,1,1) )
831 end do
832 end do
833 !$omp end do
834 end subroutine cpu_opgrad_lx14
835
836 subroutine cpu_opgrad_lx13(ux, uy, uz, u, dx, dy, dz, &
837 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
838 integer, parameter :: lx = 13
839 integer, intent(in) :: n
840 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
841 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
842 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
843 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
844 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
845 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
846 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
847 real(kind=rp) :: ur(lx, lx, lx)
848 real(kind=rp) :: us(lx, lx, lx)
849 real(kind=rp) :: ut(lx, lx, lx)
850 integer :: e, i, j, k
851
852 !$omp do
853 do e = 1, n
854 do j = 1, lx * lx
855 do i = 1, lx
856 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
857 + dx(i,2) * u(2,j,1,e) &
858 + dx(i,3) * u(3,j,1,e) &
859 + dx(i,4) * u(4,j,1,e) &
860 + dx(i,5) * u(5,j,1,e) &
861 + dx(i,6) * u(6,j,1,e) &
862 + dx(i,7) * u(7,j,1,e) &
863 + dx(i,8) * u(8,j,1,e) &
864 + dx(i,9) * u(9,j,1,e) &
865 + dx(i,10) * u(10,j,1,e) &
866 + dx(i,11) * u(11,j,1,e) &
867 + dx(i,12) * u(12,j,1,e) &
868 + dx(i,13) * u(13,j,1,e)
869 end do
870 end do
871
872 do k = 1, lx
873 do j = 1, lx
874 do i = 1, lx
875 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
876 + dy(j,2) * u(i,2,k,e) &
877 + dy(j,3) * u(i,3,k,e) &
878 + dy(j,4) * u(i,4,k,e) &
879 + dy(j,5) * u(i,5,k,e) &
880 + dy(j,6) * u(i,6,k,e) &
881 + dy(j,7) * u(i,7,k,e) &
882 + dy(j,8) * u(i,8,k,e) &
883 + dy(j,9) * u(i,9,k,e) &
884 + dy(j,10) * u(i,10,k,e) &
885 + dy(j,11) * u(i,11,k,e) &
886 + dy(j,12) * u(i,12,k,e) &
887 + dy(j,13) * u(i,13,k,e)
888 end do
889 end do
890 end do
891
892 do k = 1, lx
893 do i = 1, lx*lx
894 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
895 + dz(k,2) * u(i,1,2,e) &
896 + dz(k,3) * u(i,1,3,e) &
897 + dz(k,4) * u(i,1,4,e) &
898 + dz(k,5) * u(i,1,5,e) &
899 + dz(k,6) * u(i,1,6,e) &
900 + dz(k,7) * u(i,1,7,e) &
901 + dz(k,8) * u(i,1,8,e) &
902 + dz(k,9) * u(i,1,9,e) &
903 + dz(k,10) * u(i,1,10,e) &
904 + dz(k,11) * u(i,1,11,e) &
905 + dz(k,12) * u(i,1,12,e) &
906 + dz(k,13) * u(i,1,13,e)
907 end do
908 end do
909
910 do i = 1, lx * lx * lx
911 ux(i,1,1,e) = w3(i,1,1) &
912 * ( drdx(i,1,1,e) * ur(i,1,1) &
913 + dsdx(i,1,1,e) * us(i,1,1) &
914 + dtdx(i,1,1,e) * ut(i,1,1) )
915 uy(i,1,1,e) = w3(i,1,1) &
916 * ( dsdy(i,1,1,e) * us(i,1,1) &
917 + drdy(i,1,1,e) * ur(i,1,1) &
918 + dtdy(i,1,1,e) * ut(i,1,1) )
919 uz(i,1,1,e) = w3(i,1,1) &
920 * ( dtdz(i,1,1,e) * ut(i,1,1) &
921 + drdz(i,1,1,e) * ur(i,1,1) &
922 + dsdz(i,1,1,e) * us(i,1,1) )
923 end do
924 end do
925 !$omp end do
926 end subroutine cpu_opgrad_lx13
927
928 subroutine cpu_opgrad_lx12(ux, uy, uz, u, dx, dy, dz, &
929 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
930 integer, parameter :: lx = 12
931 integer, intent(in) :: n
932 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
933 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
934 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
935 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
936 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
937 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
938 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
939 real(kind=rp) :: ur(lx, lx, lx)
940 real(kind=rp) :: us(lx, lx, lx)
941 real(kind=rp) :: ut(lx, lx, lx)
942 integer :: e, i, j, k
943
944 !$omp do
945 do e = 1, n
946 do j = 1, lx * lx
947 do i = 1, lx
948 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
949 + dx(i,2) * u(2,j,1,e) &
950 + dx(i,3) * u(3,j,1,e) &
951 + dx(i,4) * u(4,j,1,e) &
952 + dx(i,5) * u(5,j,1,e) &
953 + dx(i,6) * u(6,j,1,e) &
954 + dx(i,7) * u(7,j,1,e) &
955 + dx(i,8) * u(8,j,1,e) &
956 + dx(i,9) * u(9,j,1,e) &
957 + dx(i,10) * u(10,j,1,e) &
958 + dx(i,11) * u(11,j,1,e) &
959 + dx(i,12) * u(12,j,1,e)
960 end do
961 end do
962
963 do k = 1, lx
964 do j = 1, lx
965 do i = 1, lx
966 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
967 + dy(j,2) * u(i,2,k,e) &
968 + dy(j,3) * u(i,3,k,e) &
969 + dy(j,4) * u(i,4,k,e) &
970 + dy(j,5) * u(i,5,k,e) &
971 + dy(j,6) * u(i,6,k,e) &
972 + dy(j,7) * u(i,7,k,e) &
973 + dy(j,8) * u(i,8,k,e) &
974 + dy(j,9) * u(i,9,k,e) &
975 + dy(j,10) * u(i,10,k,e) &
976 + dy(j,11) * u(i,11,k,e) &
977 + dy(j,12) * u(i,12,k,e)
978 end do
979 end do
980 end do
981
982 do k = 1, lx
983 do i = 1, lx*lx
984 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
985 + dz(k,2) * u(i,1,2,e) &
986 + dz(k,3) * u(i,1,3,e) &
987 + dz(k,4) * u(i,1,4,e) &
988 + dz(k,5) * u(i,1,5,e) &
989 + dz(k,6) * u(i,1,6,e) &
990 + dz(k,7) * u(i,1,7,e) &
991 + dz(k,8) * u(i,1,8,e) &
992 + dz(k,9) * u(i,1,9,e) &
993 + dz(k,10) * u(i,1,10,e) &
994 + dz(k,11) * u(i,1,11,e) &
995 + dz(k,12) * u(i,1,12,e)
996 end do
997 end do
998
999 do i = 1, lx * lx * lx
1000 ux(i,1,1,e) = w3(i,1,1) &
1001 * ( drdx(i,1,1,e) * ur(i,1,1) &
1002 + dsdx(i,1,1,e) * us(i,1,1) &
1003 + dtdx(i,1,1,e) * ut(i,1,1) )
1004 uy(i,1,1,e) = w3(i,1,1) &
1005 * ( dsdy(i,1,1,e) * us(i,1,1) &
1006 + drdy(i,1,1,e) * ur(i,1,1) &
1007 + dtdy(i,1,1,e) * ut(i,1,1) )
1008 uz(i,1,1,e) = w3(i,1,1) &
1009 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1010 + drdz(i,1,1,e) * ur(i,1,1) &
1011 + dsdz(i,1,1,e) * us(i,1,1) )
1012 end do
1013 end do
1014 !$omp end do
1015 end subroutine cpu_opgrad_lx12
1016
1017 subroutine cpu_opgrad_lx11(ux, uy, uz, u, dx, dy, dz, &
1018 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1019 integer, parameter :: lx = 11
1020 integer, intent(in) :: n
1021 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1022 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1023 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1024 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1025 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1026 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1027 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1028 real(kind=rp) :: ur(lx, lx, lx)
1029 real(kind=rp) :: us(lx, lx, lx)
1030 real(kind=rp) :: ut(lx, lx, lx)
1031 integer :: e, i, j, k
1032
1033 !$omp do
1034 do e = 1, n
1035 do j = 1, lx * lx
1036 do i = 1, lx
1037 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1038 + dx(i,2) * u(2,j,1,e) &
1039 + dx(i,3) * u(3,j,1,e) &
1040 + dx(i,4) * u(4,j,1,e) &
1041 + dx(i,5) * u(5,j,1,e) &
1042 + dx(i,6) * u(6,j,1,e) &
1043 + dx(i,7) * u(7,j,1,e) &
1044 + dx(i,8) * u(8,j,1,e) &
1045 + dx(i,9) * u(9,j,1,e) &
1046 + dx(i,10) * u(10,j,1,e) &
1047 + dx(i,11) * u(11,j,1,e)
1048 end do
1049 end do
1050
1051 do k = 1, lx
1052 do j = 1, lx
1053 do i = 1, lx
1054 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1055 + dy(j,2) * u(i,2,k,e) &
1056 + dy(j,3) * u(i,3,k,e) &
1057 + dy(j,4) * u(i,4,k,e) &
1058 + dy(j,5) * u(i,5,k,e) &
1059 + dy(j,6) * u(i,6,k,e) &
1060 + dy(j,7) * u(i,7,k,e) &
1061 + dy(j,8) * u(i,8,k,e) &
1062 + dy(j,9) * u(i,9,k,e) &
1063 + dy(j,10) * u(i,10,k,e) &
1064 + dy(j,11) * u(i,11,k,e)
1065 end do
1066 end do
1067 end do
1068
1069 do k = 1, lx
1070 do i = 1, lx*lx
1071 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1072 + dz(k,2) * u(i,1,2,e) &
1073 + dz(k,3) * u(i,1,3,e) &
1074 + dz(k,4) * u(i,1,4,e) &
1075 + dz(k,5) * u(i,1,5,e) &
1076 + dz(k,6) * u(i,1,6,e) &
1077 + dz(k,7) * u(i,1,7,e) &
1078 + dz(k,8) * u(i,1,8,e) &
1079 + dz(k,9) * u(i,1,9,e) &
1080 + dz(k,10) * u(i,1,10,e) &
1081 + dz(k,11) * u(i,1,11,e)
1082 end do
1083 end do
1084
1085 do i = 1, lx * lx * lx
1086 ux(i,1,1,e) = w3(i,1,1) &
1087 * ( drdx(i,1,1,e) * ur(i,1,1) &
1088 + dsdx(i,1,1,e) * us(i,1,1) &
1089 + dtdx(i,1,1,e) * ut(i,1,1) )
1090 uy(i,1,1,e) = w3(i,1,1) &
1091 * ( dsdy(i,1,1,e) * us(i,1,1) &
1092 + drdy(i,1,1,e) * ur(i,1,1) &
1093 + dtdy(i,1,1,e) * ut(i,1,1) )
1094 uz(i,1,1,e) = w3(i,1,1) &
1095 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1096 + drdz(i,1,1,e) * ur(i,1,1) &
1097 + dsdz(i,1,1,e) * us(i,1,1) )
1098 end do
1099 end do
1100 !$omp end do
1101 end subroutine cpu_opgrad_lx11
1102
1103 subroutine cpu_opgrad_lx10(ux, uy, uz, u, dx, dy, dz, &
1104 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1105 integer, parameter :: lx = 10
1106 integer, intent(in) :: n
1107 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1108 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1109 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1110 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1111 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1112 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1113 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1114 real(kind=rp) :: ur(lx, lx, lx)
1115 real(kind=rp) :: us(lx, lx, lx)
1116 real(kind=rp) :: ut(lx, lx, lx)
1117 integer :: e, i, j, k
1118
1119 !$omp do
1120 do e = 1, n
1121 do j = 1, lx * lx
1122 do i = 1, lx
1123 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1124 + dx(i,2) * u(2,j,1,e) &
1125 + dx(i,3) * u(3,j,1,e) &
1126 + dx(i,4) * u(4,j,1,e) &
1127 + dx(i,5) * u(5,j,1,e) &
1128 + dx(i,6) * u(6,j,1,e) &
1129 + dx(i,7) * u(7,j,1,e) &
1130 + dx(i,8) * u(8,j,1,e) &
1131 + dx(i,9) * u(9,j,1,e) &
1132 + dx(i,10) * u(10,j,1,e)
1133 end do
1134 end do
1135
1136 do k = 1, lx
1137 do j = 1, lx
1138 do i = 1, lx
1139 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1140 + dy(j,2) * u(i,2,k,e) &
1141 + dy(j,3) * u(i,3,k,e) &
1142 + dy(j,4) * u(i,4,k,e) &
1143 + dy(j,5) * u(i,5,k,e) &
1144 + dy(j,6) * u(i,6,k,e) &
1145 + dy(j,7) * u(i,7,k,e) &
1146 + dy(j,8) * u(i,8,k,e) &
1147 + dy(j,9) * u(i,9,k,e) &
1148 + dy(j,10) * u(i,10,k,e)
1149 end do
1150 end do
1151 end do
1152
1153 do k = 1, lx
1154 do i = 1, lx*lx
1155 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1156 + dz(k,2) * u(i,1,2,e) &
1157 + dz(k,3) * u(i,1,3,e) &
1158 + dz(k,4) * u(i,1,4,e) &
1159 + dz(k,5) * u(i,1,5,e) &
1160 + dz(k,6) * u(i,1,6,e) &
1161 + dz(k,7) * u(i,1,7,e) &
1162 + dz(k,8) * u(i,1,8,e) &
1163 + dz(k,9) * u(i,1,9,e) &
1164 + dz(k,10) * u(i,1,10,e)
1165 end do
1166 end do
1167
1168 do i = 1, lx * lx * lx
1169 ux(i,1,1,e) = w3(i,1,1) &
1170 * ( drdx(i,1,1,e) * ur(i,1,1) &
1171 + dsdx(i,1,1,e) * us(i,1,1) &
1172 + dtdx(i,1,1,e) * ut(i,1,1) )
1173 uy(i,1,1,e) = w3(i,1,1) &
1174 * ( dsdy(i,1,1,e) * us(i,1,1) &
1175 + drdy(i,1,1,e) * ur(i,1,1) &
1176 + dtdy(i,1,1,e) * ut(i,1,1) )
1177 uz(i,1,1,e) = w3(i,1,1) &
1178 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1179 + drdz(i,1,1,e) * ur(i,1,1) &
1180 + dsdz(i,1,1,e) * us(i,1,1) )
1181 end do
1182 end do
1183 !$omp end do
1184 end subroutine cpu_opgrad_lx10
1185
1186 subroutine cpu_opgrad_lx9(ux, uy, uz, u, dx, dy, dz, &
1187 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1188 integer, parameter :: lx = 9
1189 integer, intent(in) :: n
1190 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1191 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1192 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1193 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1194 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1195 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1196 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1197 real(kind=rp) :: ur(lx, lx, lx)
1198 real(kind=rp) :: us(lx, lx, lx)
1199 real(kind=rp) :: ut(lx, lx, lx)
1200 integer :: e, i, j, k
1201
1202 !$omp do
1203 do e = 1, n
1204 do j = 1, lx * lx
1205 do i = 1, lx
1206 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1207 + dx(i,2) * u(2,j,1,e) &
1208 + dx(i,3) * u(3,j,1,e) &
1209 + dx(i,4) * u(4,j,1,e) &
1210 + dx(i,5) * u(5,j,1,e) &
1211 + dx(i,6) * u(6,j,1,e) &
1212 + dx(i,7) * u(7,j,1,e) &
1213 + dx(i,8) * u(8,j,1,e) &
1214 + dx(i,9) * u(9,j,1,e)
1215 end do
1216 end do
1217
1218 do k = 1, lx
1219 do j = 1, lx
1220 do i = 1, lx
1221 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1222 + dy(j,2) * u(i,2,k,e) &
1223 + dy(j,3) * u(i,3,k,e) &
1224 + dy(j,4) * u(i,4,k,e) &
1225 + dy(j,5) * u(i,5,k,e) &
1226 + dy(j,6) * u(i,6,k,e) &
1227 + dy(j,7) * u(i,7,k,e) &
1228 + dy(j,8) * u(i,8,k,e) &
1229 + dy(j,9) * u(i,9,k,e)
1230 end do
1231 end do
1232 end do
1233
1234 do k = 1, lx
1235 do i = 1, lx*lx
1236 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1237 + dz(k,2) * u(i,1,2,e) &
1238 + dz(k,3) * u(i,1,3,e) &
1239 + dz(k,4) * u(i,1,4,e) &
1240 + dz(k,5) * u(i,1,5,e) &
1241 + dz(k,6) * u(i,1,6,e) &
1242 + dz(k,7) * u(i,1,7,e) &
1243 + dz(k,8) * u(i,1,8,e) &
1244 + dz(k,9) * u(i,1,9,e)
1245 end do
1246 end do
1247
1248 do i = 1, lx * lx * lx
1249 ux(i,1,1,e) = w3(i,1,1) &
1250 * ( drdx(i,1,1,e) * ur(i,1,1) &
1251 + dsdx(i,1,1,e) * us(i,1,1) &
1252 + dtdx(i,1,1,e) * ut(i,1,1) )
1253 uy(i,1,1,e) = w3(i,1,1) &
1254 * ( dsdy(i,1,1,e) * us(i,1,1) &
1255 + drdy(i,1,1,e) * ur(i,1,1) &
1256 + dtdy(i,1,1,e) * ut(i,1,1) )
1257 uz(i,1,1,e) = w3(i,1,1) &
1258 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1259 + drdz(i,1,1,e) * ur(i,1,1) &
1260 + dsdz(i,1,1,e) * us(i,1,1) )
1261 end do
1262 end do
1263 !$omp end do
1264 end subroutine cpu_opgrad_lx9
1265
1266 subroutine cpu_opgrad_lx8(ux, uy, uz, u, dx, dy, dz, &
1267 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1268 integer, parameter :: lx = 8
1269 integer, intent(in) :: n
1270 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1271 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1272 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1273 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1274 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1275 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1276 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1277 real(kind=rp) :: ur(lx, lx, lx)
1278 real(kind=rp) :: us(lx, lx, lx)
1279 real(kind=rp) :: ut(lx, lx, lx)
1280 integer :: e, i, j, k
1281
1282 !$omp do
1283 do e = 1, n
1284 do j = 1, lx * lx
1285 do i = 1, lx
1286 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1287 + dx(i,2) * u(2,j,1,e) &
1288 + dx(i,3) * u(3,j,1,e) &
1289 + dx(i,4) * u(4,j,1,e) &
1290 + dx(i,5) * u(5,j,1,e) &
1291 + dx(i,6) * u(6,j,1,e) &
1292 + dx(i,7) * u(7,j,1,e) &
1293 + dx(i,8) * u(8,j,1,e)
1294 end do
1295 end do
1296
1297 do k = 1, lx
1298 do j = 1, lx
1299 do i = 1, lx
1300 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1301 + dy(j,2) * u(i,2,k,e) &
1302 + dy(j,3) * u(i,3,k,e) &
1303 + dy(j,4) * u(i,4,k,e) &
1304 + dy(j,5) * u(i,5,k,e) &
1305 + dy(j,6) * u(i,6,k,e) &
1306 + dy(j,7) * u(i,7,k,e) &
1307 + dy(j,8) * u(i,8,k,e)
1308 end do
1309 end do
1310 end do
1311
1312 do k = 1, lx
1313 do i = 1, lx*lx
1314 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1315 + dz(k,2) * u(i,1,2,e) &
1316 + dz(k,3) * u(i,1,3,e) &
1317 + dz(k,4) * u(i,1,4,e) &
1318 + dz(k,5) * u(i,1,5,e) &
1319 + dz(k,6) * u(i,1,6,e) &
1320 + dz(k,7) * u(i,1,7,e) &
1321 + dz(k,8) * u(i,1,8,e)
1322 end do
1323 end do
1324
1325 do i = 1, lx * lx * lx
1326 ux(i,1,1,e) = w3(i,1,1) &
1327 * ( drdx(i,1,1,e) * ur(i,1,1) &
1328 + dsdx(i,1,1,e) * us(i,1,1) &
1329 + dtdx(i,1,1,e) * ut(i,1,1) )
1330 uy(i,1,1,e) = w3(i,1,1) &
1331 * ( dsdy(i,1,1,e) * us(i,1,1) &
1332 + drdy(i,1,1,e) * ur(i,1,1) &
1333 + dtdy(i,1,1,e) * ut(i,1,1) )
1334 uz(i,1,1,e) = w3(i,1,1) &
1335 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1336 + drdz(i,1,1,e) * ur(i,1,1) &
1337 + dsdz(i,1,1,e) * us(i,1,1) )
1338 end do
1339 end do
1340 !$omp end do
1341 end subroutine cpu_opgrad_lx8
1342
1343 subroutine cpu_opgrad_lx7(ux, uy, uz, u, dx, dy, dz, &
1344 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1345 integer, parameter :: lx = 7
1346 integer, intent(in) :: n
1347 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1348 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1349 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1350 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1351 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1352 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1353 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1354 real(kind=rp) :: ur(lx, lx, lx)
1355 real(kind=rp) :: us(lx, lx, lx)
1356 real(kind=rp) :: ut(lx, lx, lx)
1357 integer :: e, i, j, k
1358
1359 !$omp do
1360 do e = 1, n
1361 do j = 1, lx * lx
1362 do i = 1, lx
1363 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1364 + dx(i,2) * u(2,j,1,e) &
1365 + dx(i,3) * u(3,j,1,e) &
1366 + dx(i,4) * u(4,j,1,e) &
1367 + dx(i,5) * u(5,j,1,e) &
1368 + dx(i,6) * u(6,j,1,e) &
1369 + dx(i,7) * u(7,j,1,e)
1370 end do
1371 end do
1372
1373 do k = 1, lx
1374 do j = 1, lx
1375 do i = 1, lx
1376 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1377 + dy(j,2) * u(i,2,k,e) &
1378 + dy(j,3) * u(i,3,k,e) &
1379 + dy(j,4) * u(i,4,k,e) &
1380 + dy(j,5) * u(i,5,k,e) &
1381 + dy(j,6) * u(i,6,k,e) &
1382 + dy(j,7) * u(i,7,k,e)
1383 end do
1384 end do
1385 end do
1386
1387 do k = 1, lx
1388 do i = 1, lx*lx
1389 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1390 + dz(k,2) * u(i,1,2,e) &
1391 + dz(k,3) * u(i,1,3,e) &
1392 + dz(k,4) * u(i,1,4,e) &
1393 + dz(k,5) * u(i,1,5,e) &
1394 + dz(k,6) * u(i,1,6,e) &
1395 + dz(k,7) * u(i,1,7,e)
1396 end do
1397 end do
1398
1399 do i = 1, lx * lx * lx
1400 ux(i,1,1,e) = w3(i,1,1) &
1401 * ( drdx(i,1,1,e) * ur(i,1,1) &
1402 + dsdx(i,1,1,e) * us(i,1,1) &
1403 + dtdx(i,1,1,e) * ut(i,1,1) )
1404 uy(i,1,1,e) = w3(i,1,1) &
1405 * ( dsdy(i,1,1,e) * us(i,1,1) &
1406 + drdy(i,1,1,e) * ur(i,1,1) &
1407 + dtdy(i,1,1,e) * ut(i,1,1) )
1408 uz(i,1,1,e) = w3(i,1,1) &
1409 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1410 + drdz(i,1,1,e) * ur(i,1,1) &
1411 + dsdz(i,1,1,e) * us(i,1,1) )
1412 end do
1413 end do
1414 !$omp end do
1415 end subroutine cpu_opgrad_lx7
1416
1417 subroutine cpu_opgrad_lx6(ux, uy, uz, u, dx, dy, dz, &
1418 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1419 integer, parameter :: lx = 6
1420 integer, intent(in) :: n
1421 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1422 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1423 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1424 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1425 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1426 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1427 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1428 real(kind=rp) :: ur(lx, lx, lx)
1429 real(kind=rp) :: us(lx, lx, lx)
1430 real(kind=rp) :: ut(lx, lx, lx)
1431 integer :: e, i, j, k
1432
1433 !$omp do
1434 do e = 1, n
1435 do j = 1, lx * lx
1436 do i = 1, lx
1437 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1438 + dx(i,2) * u(2,j,1,e) &
1439 + dx(i,3) * u(3,j,1,e) &
1440 + dx(i,4) * u(4,j,1,e) &
1441 + dx(i,5) * u(5,j,1,e) &
1442 + dx(i,6) * u(6,j,1,e)
1443 end do
1444 end do
1445
1446 do k = 1, lx
1447 do j = 1, lx
1448 do i = 1, lx
1449 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1450 + dy(j,2) * u(i,2,k,e) &
1451 + dy(j,3) * u(i,3,k,e) &
1452 + dy(j,4) * u(i,4,k,e) &
1453 + dy(j,5) * u(i,5,k,e) &
1454 + dy(j,6) * u(i,6,k,e)
1455 end do
1456 end do
1457 end do
1458
1459 do k = 1, lx
1460 do i = 1, lx*lx
1461 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1462 + dz(k,2) * u(i,1,2,e) &
1463 + dz(k,3) * u(i,1,3,e) &
1464 + dz(k,4) * u(i,1,4,e) &
1465 + dz(k,5) * u(i,1,5,e) &
1466 + dz(k,6) * u(i,1,6,e)
1467 end do
1468 end do
1469
1470 do i = 1, lx * lx * lx
1471 ux(i,1,1,e) = w3(i,1,1) &
1472 * ( drdx(i,1,1,e) * ur(i,1,1) &
1473 + dsdx(i,1,1,e) * us(i,1,1) &
1474 + dtdx(i,1,1,e) * ut(i,1,1) )
1475 uy(i,1,1,e) = w3(i,1,1) &
1476 * ( dsdy(i,1,1,e) * us(i,1,1) &
1477 + drdy(i,1,1,e) * ur(i,1,1) &
1478 + dtdy(i,1,1,e) * ut(i,1,1) )
1479 uz(i,1,1,e) = w3(i,1,1) &
1480 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1481 + drdz(i,1,1,e) * ur(i,1,1) &
1482 + dsdz(i,1,1,e) * us(i,1,1) )
1483 end do
1484 end do
1485 !$omp end do
1486 end subroutine cpu_opgrad_lx6
1487
1488 subroutine cpu_opgrad_lx5(ux, uy, uz, u, dx, dy, dz, &
1489 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1490 integer, parameter :: lx = 5
1491 integer, intent(in) :: n
1492 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1493 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1494 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1495 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1496 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1497 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1498 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1499 real(kind=rp) :: ur(lx, lx, lx)
1500 real(kind=rp) :: us(lx, lx, lx)
1501 real(kind=rp) :: ut(lx, lx, lx)
1502 integer :: e, i, j, k
1503
1504 !$omp do
1505 do e = 1, n
1506 do j = 1, lx * lx
1507 do i = 1, lx
1508 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1509 + dx(i,2) * u(2,j,1,e) &
1510 + dx(i,3) * u(3,j,1,e) &
1511 + dx(i,4) * u(4,j,1,e) &
1512 + dx(i,5) * u(5,j,1,e)
1513 end do
1514 end do
1515
1516 do k = 1, lx
1517 do j = 1, lx
1518 do i = 1, lx
1519 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1520 + dy(j,2) * u(i,2,k,e) &
1521 + dy(j,3) * u(i,3,k,e) &
1522 + dy(j,4) * u(i,4,k,e) &
1523 + dy(j,5) * u(i,5,k,e)
1524 end do
1525 end do
1526 end do
1527
1528 do k = 1, lx
1529 do i = 1, lx*lx
1530 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1531 + dz(k,2) * u(i,1,2,e) &
1532 + dz(k,3) * u(i,1,3,e) &
1533 + dz(k,4) * u(i,1,4,e) &
1534 + dz(k,5) * u(i,1,5,e)
1535 end do
1536 end do
1537
1538 do i = 1, lx * lx * lx
1539 ux(i,1,1,e) = w3(i,1,1) &
1540 * ( drdx(i,1,1,e) * ur(i,1,1) &
1541 + dsdx(i,1,1,e) * us(i,1,1) &
1542 + dtdx(i,1,1,e) * ut(i,1,1) )
1543 uy(i,1,1,e) = w3(i,1,1) &
1544 * ( dsdy(i,1,1,e) * us(i,1,1) &
1545 + drdy(i,1,1,e) * ur(i,1,1) &
1546 + dtdy(i,1,1,e) * ut(i,1,1) )
1547 uz(i,1,1,e) = w3(i,1,1) &
1548 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1549 + drdz(i,1,1,e) * ur(i,1,1) &
1550 + dsdz(i,1,1,e) * us(i,1,1) )
1551 end do
1552 end do
1553 !$omp end do
1554 end subroutine cpu_opgrad_lx5
1555
1556 subroutine cpu_opgrad_lx4(ux, uy, uz, u, dx, dy, dz, &
1557 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1558 integer, parameter :: lx = 4
1559 integer, intent(in) :: n
1560 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1561 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1562 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1563 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1564 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1565 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1566 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1567 real(kind=rp) :: ur(lx, lx, lx)
1568 real(kind=rp) :: us(lx, lx, lx)
1569 real(kind=rp) :: ut(lx, lx, lx)
1570 integer :: e, i, j, k
1571
1572 !$omp do
1573 do e = 1, n
1574 do j = 1, lx * lx
1575 do i = 1, lx
1576 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1577 + dx(i,2) * u(2,j,1,e) &
1578 + dx(i,3) * u(3,j,1,e) &
1579 + dx(i,4) * u(4,j,1,e)
1580 end do
1581 end do
1582
1583 do k = 1, lx
1584 do j = 1, lx
1585 do i = 1, lx
1586 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1587 + dy(j,2) * u(i,2,k,e) &
1588 + dy(j,3) * u(i,3,k,e) &
1589 + dy(j,4) * u(i,4,k,e)
1590 end do
1591 end do
1592 end do
1593
1594 do k = 1, lx
1595 do i = 1, lx*lx
1596 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1597 + dz(k,2) * u(i,1,2,e) &
1598 + dz(k,3) * u(i,1,3,e) &
1599 + dz(k,4) * u(i,1,4,e)
1600 end do
1601 end do
1602
1603 do i = 1, lx * lx * lx
1604 ux(i,1,1,e) = w3(i,1,1) &
1605 * ( drdx(i,1,1,e) * ur(i,1,1) &
1606 + dsdx(i,1,1,e) * us(i,1,1) &
1607 + dtdx(i,1,1,e) * ut(i,1,1) )
1608 uy(i,1,1,e) = w3(i,1,1) &
1609 * ( dsdy(i,1,1,e) * us(i,1,1) &
1610 + drdy(i,1,1,e) * ur(i,1,1) &
1611 + dtdy(i,1,1,e) * ut(i,1,1) )
1612 uz(i,1,1,e) = w3(i,1,1) &
1613 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1614 + drdz(i,1,1,e) * ur(i,1,1) &
1615 + dsdz(i,1,1,e) * us(i,1,1) )
1616 end do
1617 end do
1618 !$omp end do
1619 end subroutine cpu_opgrad_lx4
1620
1621 subroutine cpu_opgrad_lx3(ux, uy, uz, u, dx, dy, dz, &
1622 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1623 integer, parameter :: lx = 3
1624 integer, intent(in) :: n
1625 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1626 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1627 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1628 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1629 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1630 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1631 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1632 real(kind=rp) :: ur(lx, lx, lx)
1633 real(kind=rp) :: us(lx, lx, lx)
1634 real(kind=rp) :: ut(lx, lx, lx)
1635 integer :: e, i, j, k
1636
1637 !$omp do
1638 do e = 1, n
1639 do j = 1, lx * lx
1640 do i = 1, lx
1641 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1642 + dx(i,2) * u(2,j,1,e) &
1643 + dx(i,3) * u(3,j,1,e)
1644 end do
1645 end do
1646
1647 do k = 1, lx
1648 do j = 1, lx
1649 do i = 1, lx
1650 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1651 + dy(j,2) * u(i,2,k,e) &
1652 + dy(j,3) * u(i,3,k,e)
1653 end do
1654 end do
1655 end do
1656
1657 do k = 1, lx
1658 do i = 1, lx*lx
1659 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1660 + dz(k,2) * u(i,1,2,e) &
1661 + dz(k,3) * u(i,1,3,e)
1662 end do
1663 end do
1664
1665 do i = 1, lx * lx * lx
1666 ux(i,1,1,e) = w3(i,1,1) &
1667 * ( drdx(i,1,1,e) * ur(i,1,1) &
1668 + dsdx(i,1,1,e) * us(i,1,1) &
1669 + dtdx(i,1,1,e) * ut(i,1,1) )
1670 uy(i,1,1,e) = w3(i,1,1) &
1671 * ( dsdy(i,1,1,e) * us(i,1,1) &
1672 + drdy(i,1,1,e) * ur(i,1,1) &
1673 + dtdy(i,1,1,e) * ut(i,1,1) )
1674 uz(i,1,1,e) = w3(i,1,1) &
1675 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1676 + drdz(i,1,1,e) * ur(i,1,1) &
1677 + dsdz(i,1,1,e) * us(i,1,1) )
1678 end do
1679 end do
1680 !$omp end do
1681 end subroutine cpu_opgrad_lx3
1682
1683 subroutine cpu_opgrad_lx2(ux, uy, uz, u, dx, dy, dz, &
1684 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1685 integer, parameter :: lx = 2
1686 integer, intent(in) :: n
1687 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1688 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1689 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1690 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1691 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1692 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1693 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1694 real(kind=rp) :: ur(lx, lx, lx)
1695 real(kind=rp) :: us(lx, lx, lx)
1696 real(kind=rp) :: ut(lx, lx, lx)
1697 integer :: e, i, j, k
1698
1699 !$omp do
1700 do e = 1, n
1701 do j = 1, lx * lx
1702 do i = 1, lx
1703 ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1704 + dx(i,2) * u(2,j,1,e)
1705 end do
1706 end do
1707
1708 do k = 1, lx
1709 do j = 1, lx
1710 do i = 1, lx
1711 us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1712 + dy(j,2) * u(i,2,k,e)
1713 end do
1714 end do
1715 end do
1716
1717 do k = 1, lx
1718 do i = 1, lx*lx
1719 ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1720 + dz(k,2) * u(i,1,2,e)
1721 end do
1722 end do
1723
1724 do i = 1, lx * lx * lx
1725 ux(i,1,1,e) = w3(i,1,1) &
1726 * ( drdx(i,1,1,e) * ur(i,1,1) &
1727 + dsdx(i,1,1,e) * us(i,1,1) &
1728 + dtdx(i,1,1,e) * ut(i,1,1) )
1729 uy(i,1,1,e) = w3(i,1,1) &
1730 * ( dsdy(i,1,1,e) * us(i,1,1) &
1731 + drdy(i,1,1,e) * ur(i,1,1) &
1732 + dtdy(i,1,1,e) * ut(i,1,1) )
1733 uz(i,1,1,e) = w3(i,1,1) &
1734 * ( dtdz(i,1,1,e) * ut(i,1,1) &
1735 + drdz(i,1,1,e) * ur(i,1,1) &
1736 + dsdz(i,1,1,e) * us(i,1,1) )
1737 end do
1738 end do
1739 !$omp end do
1740 end subroutine cpu_opgrad_lx2
1741
1742
1743 subroutine opr_cpu_opgrad_single(ux, uy, uz, u, coef, e)
1744 integer, parameter :: e_len = 1
1745 type(coef_t), intent(in) :: coef
1746 integer, intent(in) :: e
1747 real(kind=rp), dimension(coef%Xh%lxyz, e_len), intent(inout) :: ux
1748 real(kind=rp), dimension(coef%Xh%lxyz, e_len), intent(inout) :: uy
1749 real(kind=rp), dimension(coef%Xh%lxyz, e_len), intent(inout) :: uz
1750 real(kind=rp), dimension(coef%Xh%lxyz, e_len), intent(in) :: u
1751
1752 associate(xh => coef%Xh, msh => coef%msh, &
1753 drdx => coef%drdx, drdy => coef%drdy, drdz => coef%drdz, &
1754 dsdx => coef%dsdx, dsdy => coef%dsdy, dsdz => coef%dsdz, &
1755 dtdx => coef%dtdx, dtdy => coef%dtdy, dtdz => coef%dtdz)
1756
1757 select case (xh%lx)
1758 case (18)
1759 call cpu_opgrad_lx18_single(ux, uy, uz, u, &
1760 xh%dx, xh%dy, xh%dz, &
1761 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1762 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1763 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1764 xh%w3)
1765 case (17)
1766 call cpu_opgrad_lx17_single(ux, uy, uz, u, &
1767 xh%dx, xh%dy, xh%dz, &
1768 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1769 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1770 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1771 xh%w3)
1772 case (16)
1773 call cpu_opgrad_lx16_single(ux, uy, uz, u, &
1774 xh%dx, xh%dy, xh%dz, &
1775 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1776 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1777 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1778 xh%w3)
1779 case (15)
1780 call cpu_opgrad_lx15_single(ux, uy, uz, u, &
1781 xh%dx, xh%dy, xh%dz, &
1782 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1783 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1784 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1785 xh%w3)
1786 case (14)
1787 call cpu_opgrad_lx14_single(ux, uy, uz, u, &
1788 xh%dx, xh%dy, xh%dz, &
1789 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1790 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1791 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1792 xh%w3)
1793 case (13)
1794 call cpu_opgrad_lx13_single(ux, uy, uz, u, &
1795 xh%dx, xh%dy, xh%dz, &
1796 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1797 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1798 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1799 xh%w3)
1800 case (12)
1801 call cpu_opgrad_lx12_single(ux, uy, uz, u, &
1802 xh%dx, xh%dy, xh%dz, &
1803 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1804 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1805 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1806 xh%w3)
1807 case (11)
1808 call cpu_opgrad_lx11_single(ux, uy, uz, u, &
1809 xh%dx, xh%dy, xh%dz, &
1810 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1811 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1812 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1813 xh%w3)
1814 case (10)
1815 call cpu_opgrad_lx10_single(ux, uy, uz, u, &
1816 xh%dx, xh%dy, xh%dz, &
1817 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1818 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1819 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1820 xh%w3)
1821
1822 case (9)
1823 call cpu_opgrad_lx9_single(ux, uy, uz, u, &
1824 xh%dx, xh%dy, xh%dz, &
1825 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1826 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1827 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1828 xh%w3)
1829 case (8)
1830 call cpu_opgrad_lx8_single(ux, uy, uz, u, &
1831 xh%dx, xh%dy, xh%dz, &
1832 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1833 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1834 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1835 xh%w3)
1836 case (7)
1837 call cpu_opgrad_lx7_single(ux, uy, uz, u, &
1838 xh%dx, xh%dy, xh%dz, &
1839 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1840 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1841 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1842 xh%w3)
1843 case (6)
1844 call cpu_opgrad_lx6_single(ux, uy, uz, u, &
1845 xh%dx, xh%dy, xh%dz, &
1846 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1847 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1848 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1849 xh%w3)
1850 case (5)
1851 call cpu_opgrad_lx5_single(ux, uy, uz, u, &
1852 xh%dx, xh%dy, xh%dz, &
1853 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1854 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1855 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1856 xh%w3)
1857 case (4)
1858 call cpu_opgrad_lx4_single(ux, uy, uz, u, &
1859 xh%dx, xh%dy, xh%dz, &
1860 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1861 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1862 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1863 xh%w3)
1864 case (3)
1865 call cpu_opgrad_lx3_single(ux, uy, uz, u, &
1866 xh%dx, xh%dy, xh%dz, &
1867 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1868 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1869 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1870 xh%w3)
1871 case (2)
1872 call cpu_opgrad_lx2_single(ux, uy, uz, u, &
1873 xh%dx, xh%dy, xh%dz, &
1874 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1875 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1876 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1877 xh%w3)
1878 case default
1879 call cpu_opgrad_lx_single(ux, uy, uz, u, &
1880 xh%dx, xh%dy, xh%dz, &
1881 drdx(1,1,1,e), dsdx(1,1,1,e), dtdx(1,1,1,e), &
1882 drdy(1,1,1,e), dsdy(1,1,1,e), dtdy(1,1,1,e), &
1883 drdz(1,1,1,e), dsdz(1,1,1,e), dtdz(1,1,1,e), &
1884 xh%w3, xh%lx)
1885 end select
1886 end associate
1887
1888 end subroutine opr_cpu_opgrad_single
1889
1890 !OCL SERIAL
1891 subroutine cpu_opgrad_lx_single(ux, uy, uz, u, dx, dy, dz, &
1892 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, lx)
1893 integer, intent(in) :: lx
1894 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
1895 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
1896 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1897 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
1898 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
1899 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
1900 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1901 real(kind=rp) :: ur(lx, lx, lx)
1902 real(kind=rp) :: us(lx, lx, lx)
1903 real(kind=rp) :: ut(lx, lx, lx)
1904 real(kind=rp) :: tmp
1905 integer :: i, j, k, l
1906
1907 do j = 1, lx * lx
1908 do i = 1, lx
1909 tmp = 0.0_rp
1910 do k = 1, lx
1911 tmp = tmp + dx(i,k) * u(k,j,1)
1912 end do
1913 ur(i,j,1) = tmp
1914 end do
1915 end do
1916
1917 do k = 1, lx
1918 do j = 1, lx
1919 do i = 1, lx
1920 tmp = 0.0_rp
1921 do l = 1, lx
1922 tmp = tmp + dy(j,l) * u(i,l,k)
1923 end do
1924 us(i,j,k) = tmp
1925 end do
1926 end do
1927 end do
1928
1929 do k = 1, lx
1930 do i = 1, lx*lx
1931 tmp = 0.0_rp
1932 do l = 1, lx
1933 tmp = tmp + dz(k,l) * u(i,1,l)
1934 end do
1935 ut(i,1,k) = tmp
1936 end do
1937 end do
1938
1939 do i = 1, lx * lx * lx
1940 ux(i,1,1) = w3(i,1,1) &
1941 * ( drdx(i,1,1) * ur(i,1,1) &
1942 + dsdx(i,1,1) * us(i,1,1) &
1943 + dtdx(i,1,1) * ut(i,1,1) )
1944 uy(i,1,1) = w3(i,1,1) &
1945 * ( dsdy(i,1,1) * us(i,1,1) &
1946 + drdy(i,1,1) * ur(i,1,1) &
1947 + dtdy(i,1,1) * ut(i,1,1) )
1948 uz(i,1,1) = w3(i,1,1) &
1949 * ( dtdz(i,1,1) * ut(i,1,1) &
1950 + drdz(i,1,1) * ur(i,1,1) &
1951 + dsdz(i,1,1) * us(i,1,1) )
1952 end do
1953
1954 end subroutine cpu_opgrad_lx_single
1955
1956 !OCL SERIAL
1957 subroutine cpu_opgrad_lx18_single(ux, uy, uz, u, dx, dy, dz, &
1958 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
1959 integer, parameter :: lx = 18
1960 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
1961 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
1962 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1963 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
1964 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
1965 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
1966 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1967 real(kind=rp) :: ur(lx, lx, lx)
1968 real(kind=rp) :: us(lx, lx, lx)
1969 real(kind=rp) :: ut(lx, lx, lx)
1970 integer :: i, j, k
1971
1972 do j = 1, lx * lx
1973 do i = 1, lx
1974 ur(i,j,1) = dx(i,1) * u(1,j,1) &
1975 + dx(i,2) * u(2,j,1) &
1976 + dx(i,3) * u(3,j,1) &
1977 + dx(i,4) * u(4,j,1) &
1978 + dx(i,5) * u(5,j,1) &
1979 + dx(i,6) * u(6,j,1) &
1980 + dx(i,7) * u(7,j,1) &
1981 + dx(i,8) * u(8,j,1) &
1982 + dx(i,9) * u(9,j,1) &
1983 + dx(i,10) * u(10,j,1) &
1984 + dx(i,11) * u(11,j,1) &
1985 + dx(i,12) * u(12,j,1) &
1986 + dx(i,13) * u(13,j,1) &
1987 + dx(i,14) * u(14,j,1) &
1988 + dx(i,15) * u(15,j,1) &
1989 + dx(i,16) * u(16,j,1) &
1990 + dx(i,17) * u(17,j,1) &
1991 + dx(i,18) * u(18,j,1)
1992 end do
1993 end do
1994
1995 do k = 1, lx
1996 do j = 1, lx
1997 do i = 1, lx
1998 us(i,j,k) = dy(j,1) * u(i,1,k) &
1999 + dy(j,2) * u(i,2,k) &
2000 + dy(j,3) * u(i,3,k) &
2001 + dy(j,4) * u(i,4,k) &
2002 + dy(j,5) * u(i,5,k) &
2003 + dy(j,6) * u(i,6,k) &
2004 + dy(j,7) * u(i,7,k) &
2005 + dy(j,8) * u(i,8,k) &
2006 + dy(j,9) * u(i,9,k) &
2007 + dy(j,10) * u(i,10,k) &
2008 + dy(j,11) * u(i,11,k) &
2009 + dy(j,12) * u(i,12,k) &
2010 + dy(j,13) * u(i,13,k) &
2011 + dy(j,14) * u(i,14,k) &
2012 + dy(j,15) * u(i,15,k) &
2013 + dy(j,16) * u(i,16,k) &
2014 + dy(j,17) * u(i,17,k) &
2015 + dy(j,18) * u(i,18,k)
2016 end do
2017 end do
2018 end do
2019
2020 do k = 1, lx
2021 do i = 1, lx*lx
2022 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2023 + dz(k,2) * u(i,1,2) &
2024 + dz(k,3) * u(i,1,3) &
2025 + dz(k,4) * u(i,1,4) &
2026 + dz(k,5) * u(i,1,5) &
2027 + dz(k,6) * u(i,1,6) &
2028 + dz(k,7) * u(i,1,7) &
2029 + dz(k,8) * u(i,1,8) &
2030 + dz(k,9) * u(i,1,9) &
2031 + dz(k,10) * u(i,1,10) &
2032 + dz(k,11) * u(i,1,11) &
2033 + dz(k,12) * u(i,1,12) &
2034 + dz(k,13) * u(i,1,13) &
2035 + dz(k,14) * u(i,1,14) &
2036 + dz(k,15) * u(i,1,15) &
2037 + dz(k,16) * u(i,1,16) &
2038 + dz(k,17) * u(i,1,17) &
2039 + dz(k,18) * u(i,1,18)
2040 end do
2041 end do
2042
2043 do i = 1, lx * lx * lx
2044 ux(i,1,1) = w3(i,1,1) &
2045 * ( drdx(i,1,1) * ur(i,1,1) &
2046 + dsdx(i,1,1) * us(i,1,1) &
2047 + dtdx(i,1,1) * ut(i,1,1) )
2048 uy(i,1,1) = w3(i,1,1) &
2049 * ( dsdy(i,1,1) * us(i,1,1) &
2050 + drdy(i,1,1) * ur(i,1,1) &
2051 + dtdy(i,1,1) * ut(i,1,1) )
2052 uz(i,1,1) = w3(i,1,1) &
2053 * ( dtdz(i,1,1) * ut(i,1,1) &
2054 + drdz(i,1,1) * ur(i,1,1) &
2055 + dsdz(i,1,1) * us(i,1,1) )
2056 end do
2057
2058 end subroutine cpu_opgrad_lx18_single
2059
2060 !OCL SERIAL
2061 subroutine cpu_opgrad_lx17_single(ux, uy, uz, u, dx, dy, dz, &
2062 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2063 integer, parameter :: lx = 17
2064 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2065 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2066 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2067 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2068 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2069 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2070 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2071 real(kind=rp) :: ur(lx, lx, lx)
2072 real(kind=rp) :: us(lx, lx, lx)
2073 real(kind=rp) :: ut(lx, lx, lx)
2074 integer :: i, j, k
2075
2076 do j = 1, lx * lx
2077 do i = 1, lx
2078 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2079 + dx(i,2) * u(2,j,1) &
2080 + dx(i,3) * u(3,j,1) &
2081 + dx(i,4) * u(4,j,1) &
2082 + dx(i,5) * u(5,j,1) &
2083 + dx(i,6) * u(6,j,1) &
2084 + dx(i,7) * u(7,j,1) &
2085 + dx(i,8) * u(8,j,1) &
2086 + dx(i,9) * u(9,j,1) &
2087 + dx(i,10) * u(10,j,1) &
2088 + dx(i,11) * u(11,j,1) &
2089 + dx(i,12) * u(12,j,1) &
2090 + dx(i,13) * u(13,j,1) &
2091 + dx(i,14) * u(14,j,1) &
2092 + dx(i,15) * u(15,j,1) &
2093 + dx(i,16) * u(16,j,1) &
2094 + dx(i,17) * u(17,j,1)
2095 end do
2096 end do
2097
2098 do k = 1, lx
2099 do j = 1, lx
2100 do i = 1, lx
2101 us(i,j,k) = dy(j,1) * u(i,1,k) &
2102 + dy(j,2) * u(i,2,k) &
2103 + dy(j,3) * u(i,3,k) &
2104 + dy(j,4) * u(i,4,k) &
2105 + dy(j,5) * u(i,5,k) &
2106 + dy(j,6) * u(i,6,k) &
2107 + dy(j,7) * u(i,7,k) &
2108 + dy(j,8) * u(i,8,k) &
2109 + dy(j,9) * u(i,9,k) &
2110 + dy(j,10) * u(i,10,k) &
2111 + dy(j,11) * u(i,11,k) &
2112 + dy(j,12) * u(i,12,k) &
2113 + dy(j,13) * u(i,13,k) &
2114 + dy(j,14) * u(i,14,k) &
2115 + dy(j,15) * u(i,15,k) &
2116 + dy(j,16) * u(i,16,k) &
2117 + dy(j,17) * u(i,17,k)
2118 end do
2119 end do
2120 end do
2121
2122 do k = 1, lx
2123 do i = 1, lx*lx
2124 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2125 + dz(k,2) * u(i,1,2) &
2126 + dz(k,3) * u(i,1,3) &
2127 + dz(k,4) * u(i,1,4) &
2128 + dz(k,5) * u(i,1,5) &
2129 + dz(k,6) * u(i,1,6) &
2130 + dz(k,7) * u(i,1,7) &
2131 + dz(k,8) * u(i,1,8) &
2132 + dz(k,9) * u(i,1,9) &
2133 + dz(k,10) * u(i,1,10) &
2134 + dz(k,11) * u(i,1,11) &
2135 + dz(k,12) * u(i,1,12) &
2136 + dz(k,13) * u(i,1,13) &
2137 + dz(k,14) * u(i,1,14) &
2138 + dz(k,15) * u(i,1,15) &
2139 + dz(k,16) * u(i,1,16) &
2140 + dz(k,17) * u(i,1,17)
2141 end do
2142 end do
2143
2144 do i = 1, lx * lx * lx
2145 ux(i,1,1) = w3(i,1,1) &
2146 * ( drdx(i,1,1) * ur(i,1,1) &
2147 + dsdx(i,1,1) * us(i,1,1) &
2148 + dtdx(i,1,1) * ut(i,1,1) )
2149 uy(i,1,1) = w3(i,1,1) &
2150 * ( dsdy(i,1,1) * us(i,1,1) &
2151 + drdy(i,1,1) * ur(i,1,1) &
2152 + dtdy(i,1,1) * ut(i,1,1) )
2153 uz(i,1,1) = w3(i,1,1) &
2154 * ( dtdz(i,1,1) * ut(i,1,1) &
2155 + drdz(i,1,1) * ur(i,1,1) &
2156 + dsdz(i,1,1) * us(i,1,1) )
2157 end do
2158 end subroutine cpu_opgrad_lx17_single
2159
2160 !OCL SERIAL
2161 subroutine cpu_opgrad_lx16_single(ux, uy, uz, u, dx, dy, dz, &
2162 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2163 integer, parameter :: lx = 16
2164 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2165 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2166 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2167 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2168 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2169 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2170 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2171 real(kind=rp) :: ur(lx, lx, lx)
2172 real(kind=rp) :: us(lx, lx, lx)
2173 real(kind=rp) :: ut(lx, lx, lx)
2174 integer :: i, j, k
2175
2176 do j = 1, lx * lx
2177 do i = 1, lx
2178 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2179 + dx(i,2) * u(2,j,1) &
2180 + dx(i,3) * u(3,j,1) &
2181 + dx(i,4) * u(4,j,1) &
2182 + dx(i,5) * u(5,j,1) &
2183 + dx(i,6) * u(6,j,1) &
2184 + dx(i,7) * u(7,j,1) &
2185 + dx(i,8) * u(8,j,1) &
2186 + dx(i,9) * u(9,j,1) &
2187 + dx(i,10) * u(10,j,1) &
2188 + dx(i,11) * u(11,j,1) &
2189 + dx(i,12) * u(12,j,1) &
2190 + dx(i,13) * u(13,j,1) &
2191 + dx(i,14) * u(14,j,1) &
2192 + dx(i,15) * u(15,j,1) &
2193 + dx(i,16) * u(16,j,1)
2194 end do
2195 end do
2196
2197 do k = 1, lx
2198 do j = 1, lx
2199 do i = 1, lx
2200 us(i,j,k) = dy(j,1) * u(i,1,k) &
2201 + dy(j,2) * u(i,2,k) &
2202 + dy(j,3) * u(i,3,k) &
2203 + dy(j,4) * u(i,4,k) &
2204 + dy(j,5) * u(i,5,k) &
2205 + dy(j,6) * u(i,6,k) &
2206 + dy(j,7) * u(i,7,k) &
2207 + dy(j,8) * u(i,8,k) &
2208 + dy(j,9) * u(i,9,k) &
2209 + dy(j,10) * u(i,10,k) &
2210 + dy(j,11) * u(i,11,k) &
2211 + dy(j,12) * u(i,12,k) &
2212 + dy(j,13) * u(i,13,k) &
2213 + dy(j,14) * u(i,14,k) &
2214 + dy(j,15) * u(i,15,k) &
2215 + dy(j,16) * u(i,16,k)
2216 end do
2217 end do
2218 end do
2219
2220 do k = 1, lx
2221 do i = 1, lx*lx
2222 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2223 + dz(k,2) * u(i,1,2) &
2224 + dz(k,3) * u(i,1,3) &
2225 + dz(k,4) * u(i,1,4) &
2226 + dz(k,5) * u(i,1,5) &
2227 + dz(k,6) * u(i,1,6) &
2228 + dz(k,7) * u(i,1,7) &
2229 + dz(k,8) * u(i,1,8) &
2230 + dz(k,9) * u(i,1,9) &
2231 + dz(k,10) * u(i,1,10) &
2232 + dz(k,11) * u(i,1,11) &
2233 + dz(k,12) * u(i,1,12) &
2234 + dz(k,13) * u(i,1,13) &
2235 + dz(k,14) * u(i,1,14) &
2236 + dz(k,15) * u(i,1,15) &
2237 + dz(k,16) * u(i,1,16)
2238 end do
2239 end do
2240
2241 do i = 1, lx * lx * lx
2242 ux(i,1,1) = w3(i,1,1) &
2243 * ( drdx(i,1,1) * ur(i,1,1) &
2244 + dsdx(i,1,1) * us(i,1,1) &
2245 + dtdx(i,1,1) * ut(i,1,1) )
2246 uy(i,1,1) = w3(i,1,1) &
2247 * ( dsdy(i,1,1) * us(i,1,1) &
2248 + drdy(i,1,1) * ur(i,1,1) &
2249 + dtdy(i,1,1) * ut(i,1,1) )
2250 uz(i,1,1) = w3(i,1,1) &
2251 * ( dtdz(i,1,1) * ut(i,1,1) &
2252 + drdz(i,1,1) * ur(i,1,1) &
2253 + dsdz(i,1,1) * us(i,1,1) )
2254 end do
2255 end subroutine cpu_opgrad_lx16_single
2256
2257 !OCL SERIAL
2258 subroutine cpu_opgrad_lx15_single(ux, uy, uz, u, dx, dy, dz, &
2259 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2260 integer, parameter :: lx = 15
2261 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2262 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2263 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2264 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2265 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2266 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2267 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2268 real(kind=rp) :: ur(lx, lx, lx)
2269 real(kind=rp) :: us(lx, lx, lx)
2270 real(kind=rp) :: ut(lx, lx, lx)
2271 integer :: i, j, k
2272
2273 do j = 1, lx * lx
2274 do i = 1, lx
2275 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2276 + dx(i,2) * u(2,j,1) &
2277 + dx(i,3) * u(3,j,1) &
2278 + dx(i,4) * u(4,j,1) &
2279 + dx(i,5) * u(5,j,1) &
2280 + dx(i,6) * u(6,j,1) &
2281 + dx(i,7) * u(7,j,1) &
2282 + dx(i,8) * u(8,j,1) &
2283 + dx(i,9) * u(9,j,1) &
2284 + dx(i,10) * u(10,j,1) &
2285 + dx(i,11) * u(11,j,1) &
2286 + dx(i,12) * u(12,j,1) &
2287 + dx(i,13) * u(13,j,1) &
2288 + dx(i,14) * u(14,j,1) &
2289 + dx(i,15) * u(15,j,1)
2290 end do
2291 end do
2292
2293 do k = 1, lx
2294 do j = 1, lx
2295 do i = 1, lx
2296 us(i,j,k) = dy(j,1) * u(i,1,k) &
2297 + dy(j,2) * u(i,2,k) &
2298 + dy(j,3) * u(i,3,k) &
2299 + dy(j,4) * u(i,4,k) &
2300 + dy(j,5) * u(i,5,k) &
2301 + dy(j,6) * u(i,6,k) &
2302 + dy(j,7) * u(i,7,k) &
2303 + dy(j,8) * u(i,8,k) &
2304 + dy(j,9) * u(i,9,k) &
2305 + dy(j,10) * u(i,10,k) &
2306 + dy(j,11) * u(i,11,k) &
2307 + dy(j,12) * u(i,12,k) &
2308 + dy(j,13) * u(i,13,k) &
2309 + dy(j,14) * u(i,14,k) &
2310 + dy(j,15) * u(i,15,k)
2311 end do
2312 end do
2313 end do
2314
2315 do k = 1, lx
2316 do i = 1, lx*lx
2317 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2318 + dz(k,2) * u(i,1,2) &
2319 + dz(k,3) * u(i,1,3) &
2320 + dz(k,4) * u(i,1,4) &
2321 + dz(k,5) * u(i,1,5) &
2322 + dz(k,6) * u(i,1,6) &
2323 + dz(k,7) * u(i,1,7) &
2324 + dz(k,8) * u(i,1,8) &
2325 + dz(k,9) * u(i,1,9) &
2326 + dz(k,10) * u(i,1,10) &
2327 + dz(k,11) * u(i,1,11) &
2328 + dz(k,12) * u(i,1,12) &
2329 + dz(k,13) * u(i,1,13) &
2330 + dz(k,14) * u(i,1,14) &
2331 + dz(k,15) * u(i,1,15)
2332 end do
2333 end do
2334
2335 do i = 1, lx * lx * lx
2336 ux(i,1,1) = w3(i,1,1) &
2337 * ( drdx(i,1,1) * ur(i,1,1) &
2338 + dsdx(i,1,1) * us(i,1,1) &
2339 + dtdx(i,1,1) * ut(i,1,1) )
2340 uy(i,1,1) = w3(i,1,1) &
2341 * ( dsdy(i,1,1) * us(i,1,1) &
2342 + drdy(i,1,1) * ur(i,1,1) &
2343 + dtdy(i,1,1) * ut(i,1,1) )
2344 uz(i,1,1) = w3(i,1,1) &
2345 * ( dtdz(i,1,1) * ut(i,1,1) &
2346 + drdz(i,1,1) * ur(i,1,1) &
2347 + dsdz(i,1,1) * us(i,1,1) )
2348 end do
2349 end subroutine cpu_opgrad_lx15_single
2350
2351 !OCL SERIAL
2352 subroutine cpu_opgrad_lx14_single(ux, uy, uz, u, dx, dy, dz, &
2353 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2354 integer, parameter :: lx = 14
2355 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2356 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2357 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2358 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2359 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2360 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2361 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2362 real(kind=rp) :: ur(lx, lx, lx)
2363 real(kind=rp) :: us(lx, lx, lx)
2364 real(kind=rp) :: ut(lx, lx, lx)
2365 integer :: i, j, k
2366
2367 do j = 1, lx * lx
2368 do i = 1, lx
2369 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2370 + dx(i,2) * u(2,j,1) &
2371 + dx(i,3) * u(3,j,1) &
2372 + dx(i,4) * u(4,j,1) &
2373 + dx(i,5) * u(5,j,1) &
2374 + dx(i,6) * u(6,j,1) &
2375 + dx(i,7) * u(7,j,1) &
2376 + dx(i,8) * u(8,j,1) &
2377 + dx(i,9) * u(9,j,1) &
2378 + dx(i,10) * u(10,j,1) &
2379 + dx(i,11) * u(11,j,1) &
2380 + dx(i,12) * u(12,j,1) &
2381 + dx(i,13) * u(13,j,1) &
2382 + dx(i,14) * u(14,j,1)
2383 end do
2384 end do
2385
2386 do k = 1, lx
2387 do j = 1, lx
2388 do i = 1, lx
2389 us(i,j,k) = dy(j,1) * u(i,1,k) &
2390 + dy(j,2) * u(i,2,k) &
2391 + dy(j,3) * u(i,3,k) &
2392 + dy(j,4) * u(i,4,k) &
2393 + dy(j,5) * u(i,5,k) &
2394 + dy(j,6) * u(i,6,k) &
2395 + dy(j,7) * u(i,7,k) &
2396 + dy(j,8) * u(i,8,k) &
2397 + dy(j,9) * u(i,9,k) &
2398 + dy(j,10) * u(i,10,k) &
2399 + dy(j,11) * u(i,11,k) &
2400 + dy(j,12) * u(i,12,k) &
2401 + dy(j,13) * u(i,13,k) &
2402 + dy(j,14) * u(i,14,k)
2403 end do
2404 end do
2405 end do
2406
2407 do k = 1, lx
2408 do i = 1, lx*lx
2409 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2410 + dz(k,2) * u(i,1,2) &
2411 + dz(k,3) * u(i,1,3) &
2412 + dz(k,4) * u(i,1,4) &
2413 + dz(k,5) * u(i,1,5) &
2414 + dz(k,6) * u(i,1,6) &
2415 + dz(k,7) * u(i,1,7) &
2416 + dz(k,8) * u(i,1,8) &
2417 + dz(k,9) * u(i,1,9) &
2418 + dz(k,10) * u(i,1,10) &
2419 + dz(k,11) * u(i,1,11) &
2420 + dz(k,12) * u(i,1,12) &
2421 + dz(k,13) * u(i,1,13) &
2422 + dz(k,14) * u(i,1,14)
2423 end do
2424 end do
2425
2426 do i = 1, lx * lx * lx
2427 ux(i,1,1) = w3(i,1,1) &
2428 * ( drdx(i,1,1) * ur(i,1,1) &
2429 + dsdx(i,1,1) * us(i,1,1) &
2430 + dtdx(i,1,1) * ut(i,1,1) )
2431 uy(i,1,1) = w3(i,1,1) &
2432 * ( dsdy(i,1,1) * us(i,1,1) &
2433 + drdy(i,1,1) * ur(i,1,1) &
2434 + dtdy(i,1,1) * ut(i,1,1) )
2435 uz(i,1,1) = w3(i,1,1) &
2436 * ( dtdz(i,1,1) * ut(i,1,1) &
2437 + drdz(i,1,1) * ur(i,1,1) &
2438 + dsdz(i,1,1) * us(i,1,1) )
2439 end do
2440 end subroutine cpu_opgrad_lx14_single
2441
2442 !OCL SERIAL
2443 subroutine cpu_opgrad_lx13_single(ux, uy, uz, u, dx, dy, dz, &
2444 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2445 integer, parameter :: lx = 13
2446 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2447 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2448 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2449 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2450 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2451 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2452 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2453 real(kind=rp) :: ur(lx, lx, lx)
2454 real(kind=rp) :: us(lx, lx, lx)
2455 real(kind=rp) :: ut(lx, lx, lx)
2456 integer :: i, j, k
2457
2458 do j = 1, lx * lx
2459 do i = 1, lx
2460 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2461 + dx(i,2) * u(2,j,1) &
2462 + dx(i,3) * u(3,j,1) &
2463 + dx(i,4) * u(4,j,1) &
2464 + dx(i,5) * u(5,j,1) &
2465 + dx(i,6) * u(6,j,1) &
2466 + dx(i,7) * u(7,j,1) &
2467 + dx(i,8) * u(8,j,1) &
2468 + dx(i,9) * u(9,j,1) &
2469 + dx(i,10) * u(10,j,1) &
2470 + dx(i,11) * u(11,j,1) &
2471 + dx(i,12) * u(12,j,1) &
2472 + dx(i,13) * u(13,j,1)
2473 end do
2474 end do
2475
2476 do k = 1, lx
2477 do j = 1, lx
2478 do i = 1, lx
2479 us(i,j,k) = dy(j,1) * u(i,1,k) &
2480 + dy(j,2) * u(i,2,k) &
2481 + dy(j,3) * u(i,3,k) &
2482 + dy(j,4) * u(i,4,k) &
2483 + dy(j,5) * u(i,5,k) &
2484 + dy(j,6) * u(i,6,k) &
2485 + dy(j,7) * u(i,7,k) &
2486 + dy(j,8) * u(i,8,k) &
2487 + dy(j,9) * u(i,9,k) &
2488 + dy(j,10) * u(i,10,k) &
2489 + dy(j,11) * u(i,11,k) &
2490 + dy(j,12) * u(i,12,k) &
2491 + dy(j,13) * u(i,13,k)
2492 end do
2493 end do
2494 end do
2495
2496 do k = 1, lx
2497 do i = 1, lx*lx
2498 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2499 + dz(k,2) * u(i,1,2) &
2500 + dz(k,3) * u(i,1,3) &
2501 + dz(k,4) * u(i,1,4) &
2502 + dz(k,5) * u(i,1,5) &
2503 + dz(k,6) * u(i,1,6) &
2504 + dz(k,7) * u(i,1,7) &
2505 + dz(k,8) * u(i,1,8) &
2506 + dz(k,9) * u(i,1,9) &
2507 + dz(k,10) * u(i,1,10) &
2508 + dz(k,11) * u(i,1,11) &
2509 + dz(k,12) * u(i,1,12) &
2510 + dz(k,13) * u(i,1,13)
2511 end do
2512 end do
2513
2514 do i = 1, lx * lx * lx
2515 ux(i,1,1) = w3(i,1,1) &
2516 * ( drdx(i,1,1) * ur(i,1,1) &
2517 + dsdx(i,1,1) * us(i,1,1) &
2518 + dtdx(i,1,1) * ut(i,1,1) )
2519 uy(i,1,1) = w3(i,1,1) &
2520 * ( dsdy(i,1,1) * us(i,1,1) &
2521 + drdy(i,1,1) * ur(i,1,1) &
2522 + dtdy(i,1,1) * ut(i,1,1) )
2523 uz(i,1,1) = w3(i,1,1) &
2524 * ( dtdz(i,1,1) * ut(i,1,1) &
2525 + drdz(i,1,1) * ur(i,1,1) &
2526 + dsdz(i,1,1) * us(i,1,1) )
2527 end do
2528 end subroutine cpu_opgrad_lx13_single
2529
2530 !OCL SERIAL
2531 subroutine cpu_opgrad_lx12_single(ux, uy, uz, u, dx, dy, dz, &
2532 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2533 integer, parameter :: lx = 12
2534 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2535 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2536 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2537 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2538 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2539 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2540 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2541 real(kind=rp) :: ur(lx, lx, lx)
2542 real(kind=rp) :: us(lx, lx, lx)
2543 real(kind=rp) :: ut(lx, lx, lx)
2544 integer :: i, j, k
2545
2546 do j = 1, lx * lx
2547 do i = 1, lx
2548 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2549 + dx(i,2) * u(2,j,1) &
2550 + dx(i,3) * u(3,j,1) &
2551 + dx(i,4) * u(4,j,1) &
2552 + dx(i,5) * u(5,j,1) &
2553 + dx(i,6) * u(6,j,1) &
2554 + dx(i,7) * u(7,j,1) &
2555 + dx(i,8) * u(8,j,1) &
2556 + dx(i,9) * u(9,j,1) &
2557 + dx(i,10) * u(10,j,1) &
2558 + dx(i,11) * u(11,j,1) &
2559 + dx(i,12) * u(12,j,1)
2560 end do
2561 end do
2562
2563 do k = 1, lx
2564 do j = 1, lx
2565 do i = 1, lx
2566 us(i,j,k) = dy(j,1) * u(i,1,k) &
2567 + dy(j,2) * u(i,2,k) &
2568 + dy(j,3) * u(i,3,k) &
2569 + dy(j,4) * u(i,4,k) &
2570 + dy(j,5) * u(i,5,k) &
2571 + dy(j,6) * u(i,6,k) &
2572 + dy(j,7) * u(i,7,k) &
2573 + dy(j,8) * u(i,8,k) &
2574 + dy(j,9) * u(i,9,k) &
2575 + dy(j,10) * u(i,10,k) &
2576 + dy(j,11) * u(i,11,k) &
2577 + dy(j,12) * u(i,12,k)
2578 end do
2579 end do
2580 end do
2581
2582 do k = 1, lx
2583 do i = 1, lx*lx
2584 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2585 + dz(k,2) * u(i,1,2) &
2586 + dz(k,3) * u(i,1,3) &
2587 + dz(k,4) * u(i,1,4) &
2588 + dz(k,5) * u(i,1,5) &
2589 + dz(k,6) * u(i,1,6) &
2590 + dz(k,7) * u(i,1,7) &
2591 + dz(k,8) * u(i,1,8) &
2592 + dz(k,9) * u(i,1,9) &
2593 + dz(k,10) * u(i,1,10) &
2594 + dz(k,11) * u(i,1,11) &
2595 + dz(k,12) * u(i,1,12)
2596 end do
2597 end do
2598
2599 do i = 1, lx * lx * lx
2600 ux(i,1,1) = w3(i,1,1) &
2601 * ( drdx(i,1,1) * ur(i,1,1) &
2602 + dsdx(i,1,1) * us(i,1,1) &
2603 + dtdx(i,1,1) * ut(i,1,1) )
2604 uy(i,1,1) = w3(i,1,1) &
2605 * ( dsdy(i,1,1) * us(i,1,1) &
2606 + drdy(i,1,1) * ur(i,1,1) &
2607 + dtdy(i,1,1) * ut(i,1,1) )
2608 uz(i,1,1) = w3(i,1,1) &
2609 * ( dtdz(i,1,1) * ut(i,1,1) &
2610 + drdz(i,1,1) * ur(i,1,1) &
2611 + dsdz(i,1,1) * us(i,1,1) )
2612 end do
2613 end subroutine cpu_opgrad_lx12_single
2614
2615 !OCL SERIAL
2616 subroutine cpu_opgrad_lx11_single(ux, uy, uz, u, dx, dy, dz, &
2617 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2618 integer, parameter :: lx = 11
2619 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2620 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2621 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2622 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2623 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2624 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2625 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2626 real(kind=rp) :: ur(lx, lx, lx)
2627 real(kind=rp) :: us(lx, lx, lx)
2628 real(kind=rp) :: ut(lx, lx, lx)
2629 integer :: i, j, k
2630
2631 do j = 1, lx * lx
2632 do i = 1, lx
2633 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2634 + dx(i,2) * u(2,j,1) &
2635 + dx(i,3) * u(3,j,1) &
2636 + dx(i,4) * u(4,j,1) &
2637 + dx(i,5) * u(5,j,1) &
2638 + dx(i,6) * u(6,j,1) &
2639 + dx(i,7) * u(7,j,1) &
2640 + dx(i,8) * u(8,j,1) &
2641 + dx(i,9) * u(9,j,1) &
2642 + dx(i,10) * u(10,j,1) &
2643 + dx(i,11) * u(11,j,1)
2644 end do
2645 end do
2646
2647 do k = 1, lx
2648 do j = 1, lx
2649 do i = 1, lx
2650 us(i,j,k) = dy(j,1) * u(i,1,k) &
2651 + dy(j,2) * u(i,2,k) &
2652 + dy(j,3) * u(i,3,k) &
2653 + dy(j,4) * u(i,4,k) &
2654 + dy(j,5) * u(i,5,k) &
2655 + dy(j,6) * u(i,6,k) &
2656 + dy(j,7) * u(i,7,k) &
2657 + dy(j,8) * u(i,8,k) &
2658 + dy(j,9) * u(i,9,k) &
2659 + dy(j,10) * u(i,10,k) &
2660 + dy(j,11) * u(i,11,k)
2661 end do
2662 end do
2663 end do
2664
2665 do k = 1, lx
2666 do i = 1, lx*lx
2667 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2668 + dz(k,2) * u(i,1,2) &
2669 + dz(k,3) * u(i,1,3) &
2670 + dz(k,4) * u(i,1,4) &
2671 + dz(k,5) * u(i,1,5) &
2672 + dz(k,6) * u(i,1,6) &
2673 + dz(k,7) * u(i,1,7) &
2674 + dz(k,8) * u(i,1,8) &
2675 + dz(k,9) * u(i,1,9) &
2676 + dz(k,10) * u(i,1,10) &
2677 + dz(k,11) * u(i,1,11)
2678 end do
2679 end do
2680
2681 do i = 1, lx * lx * lx
2682 ux(i,1,1) = w3(i,1,1) &
2683 * ( drdx(i,1,1) * ur(i,1,1) &
2684 + dsdx(i,1,1) * us(i,1,1) &
2685 + dtdx(i,1,1) * ut(i,1,1) )
2686 uy(i,1,1) = w3(i,1,1) &
2687 * ( dsdy(i,1,1) * us(i,1,1) &
2688 + drdy(i,1,1) * ur(i,1,1) &
2689 + dtdy(i,1,1) * ut(i,1,1) )
2690 uz(i,1,1) = w3(i,1,1) &
2691 * ( dtdz(i,1,1) * ut(i,1,1) &
2692 + drdz(i,1,1) * ur(i,1,1) &
2693 + dsdz(i,1,1) * us(i,1,1) )
2694 end do
2695 end subroutine cpu_opgrad_lx11_single
2696
2697 !OCL SERIAL
2698 subroutine cpu_opgrad_lx10_single(ux, uy, uz, u, dx, dy, dz, &
2699 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2700 integer, parameter :: lx = 10
2701 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2702 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2703 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2704 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2705 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2706 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2707 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2708 real(kind=rp) :: ur(lx, lx, lx)
2709 real(kind=rp) :: us(lx, lx, lx)
2710 real(kind=rp) :: ut(lx, lx, lx)
2711 integer :: i, j, k
2712
2713 do j = 1, lx * lx
2714 do i = 1, lx
2715 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2716 + dx(i,2) * u(2,j,1) &
2717 + dx(i,3) * u(3,j,1) &
2718 + dx(i,4) * u(4,j,1) &
2719 + dx(i,5) * u(5,j,1) &
2720 + dx(i,6) * u(6,j,1) &
2721 + dx(i,7) * u(7,j,1) &
2722 + dx(i,8) * u(8,j,1) &
2723 + dx(i,9) * u(9,j,1) &
2724 + dx(i,10) * u(10,j,1)
2725 end do
2726 end do
2727
2728 do k = 1, lx
2729 do j = 1, lx
2730 do i = 1, lx
2731 us(i,j,k) = dy(j,1) * u(i,1,k) &
2732 + dy(j,2) * u(i,2,k) &
2733 + dy(j,3) * u(i,3,k) &
2734 + dy(j,4) * u(i,4,k) &
2735 + dy(j,5) * u(i,5,k) &
2736 + dy(j,6) * u(i,6,k) &
2737 + dy(j,7) * u(i,7,k) &
2738 + dy(j,8) * u(i,8,k) &
2739 + dy(j,9) * u(i,9,k) &
2740 + dy(j,10) * u(i,10,k)
2741 end do
2742 end do
2743 end do
2744
2745 do k = 1, lx
2746 do i = 1, lx*lx
2747 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2748 + dz(k,2) * u(i,1,2) &
2749 + dz(k,3) * u(i,1,3) &
2750 + dz(k,4) * u(i,1,4) &
2751 + dz(k,5) * u(i,1,5) &
2752 + dz(k,6) * u(i,1,6) &
2753 + dz(k,7) * u(i,1,7) &
2754 + dz(k,8) * u(i,1,8) &
2755 + dz(k,9) * u(i,1,9) &
2756 + dz(k,10) * u(i,1,10)
2757 end do
2758 end do
2759
2760 do i = 1, lx * lx * lx
2761 ux(i,1,1) = w3(i,1,1) &
2762 * ( drdx(i,1,1) * ur(i,1,1) &
2763 + dsdx(i,1,1) * us(i,1,1) &
2764 + dtdx(i,1,1) * ut(i,1,1) )
2765 uy(i,1,1) = w3(i,1,1) &
2766 * ( dsdy(i,1,1) * us(i,1,1) &
2767 + drdy(i,1,1) * ur(i,1,1) &
2768 + dtdy(i,1,1) * ut(i,1,1) )
2769 uz(i,1,1) = w3(i,1,1) &
2770 * ( dtdz(i,1,1) * ut(i,1,1) &
2771 + drdz(i,1,1) * ur(i,1,1) &
2772 + dsdz(i,1,1) * us(i,1,1) )
2773 end do
2774 end subroutine cpu_opgrad_lx10_single
2775
2776 !OCL SERIAL
2777 subroutine cpu_opgrad_lx9_single(ux, uy, uz, u, dx, dy, dz, &
2778 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2779 integer, parameter :: lx = 9
2780 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2781 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2782 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2783 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2784 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2785 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2786 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2787 real(kind=rp) :: ur(lx, lx, lx)
2788 real(kind=rp) :: us(lx, lx, lx)
2789 real(kind=rp) :: ut(lx, lx, lx)
2790 integer :: i, j, k
2791
2792 do j = 1, lx * lx
2793 do i = 1, lx
2794 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2795 + dx(i,2) * u(2,j,1) &
2796 + dx(i,3) * u(3,j,1) &
2797 + dx(i,4) * u(4,j,1) &
2798 + dx(i,5) * u(5,j,1) &
2799 + dx(i,6) * u(6,j,1) &
2800 + dx(i,7) * u(7,j,1) &
2801 + dx(i,8) * u(8,j,1) &
2802 + dx(i,9) * u(9,j,1)
2803 end do
2804 end do
2805
2806 do k = 1, lx
2807 do j = 1, lx
2808 do i = 1, lx
2809 us(i,j,k) = dy(j,1) * u(i,1,k) &
2810 + dy(j,2) * u(i,2,k) &
2811 + dy(j,3) * u(i,3,k) &
2812 + dy(j,4) * u(i,4,k) &
2813 + dy(j,5) * u(i,5,k) &
2814 + dy(j,6) * u(i,6,k) &
2815 + dy(j,7) * u(i,7,k) &
2816 + dy(j,8) * u(i,8,k) &
2817 + dy(j,9) * u(i,9,k)
2818 end do
2819 end do
2820 end do
2821
2822 do k = 1, lx
2823 do i = 1, lx*lx
2824 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2825 + dz(k,2) * u(i,1,2) &
2826 + dz(k,3) * u(i,1,3) &
2827 + dz(k,4) * u(i,1,4) &
2828 + dz(k,5) * u(i,1,5) &
2829 + dz(k,6) * u(i,1,6) &
2830 + dz(k,7) * u(i,1,7) &
2831 + dz(k,8) * u(i,1,8) &
2832 + dz(k,9) * u(i,1,9)
2833 end do
2834 end do
2835
2836 do i = 1, lx * lx * lx
2837 ux(i,1,1) = w3(i,1,1) &
2838 * ( drdx(i,1,1) * ur(i,1,1) &
2839 + dsdx(i,1,1) * us(i,1,1) &
2840 + dtdx(i,1,1) * ut(i,1,1) )
2841 uy(i,1,1) = w3(i,1,1) &
2842 * ( dsdy(i,1,1) * us(i,1,1) &
2843 + drdy(i,1,1) * ur(i,1,1) &
2844 + dtdy(i,1,1) * ut(i,1,1) )
2845 uz(i,1,1) = w3(i,1,1) &
2846 * ( dtdz(i,1,1) * ut(i,1,1) &
2847 + drdz(i,1,1) * ur(i,1,1) &
2848 + dsdz(i,1,1) * us(i,1,1) )
2849 end do
2850 end subroutine cpu_opgrad_lx9_single
2851
2852 !OCL SERIAL
2853 subroutine cpu_opgrad_lx8_single(ux, uy, uz, u, dx, dy, dz, &
2854 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2855 integer, parameter :: lx = 8
2856 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2857 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2858 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2859 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2860 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2861 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2862 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2863 real(kind=rp) :: ur(lx, lx, lx)
2864 real(kind=rp) :: us(lx, lx, lx)
2865 real(kind=rp) :: ut(lx, lx, lx)
2866 integer :: i, j, k
2867
2868 do j = 1, lx * lx
2869 do i = 1, lx
2870 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2871 + dx(i,2) * u(2,j,1) &
2872 + dx(i,3) * u(3,j,1) &
2873 + dx(i,4) * u(4,j,1) &
2874 + dx(i,5) * u(5,j,1) &
2875 + dx(i,6) * u(6,j,1) &
2876 + dx(i,7) * u(7,j,1) &
2877 + dx(i,8) * u(8,j,1)
2878 end do
2879 end do
2880
2881 do k = 1, lx
2882 do j = 1, lx
2883 do i = 1, lx
2884 us(i,j,k) = dy(j,1) * u(i,1,k) &
2885 + dy(j,2) * u(i,2,k) &
2886 + dy(j,3) * u(i,3,k) &
2887 + dy(j,4) * u(i,4,k) &
2888 + dy(j,5) * u(i,5,k) &
2889 + dy(j,6) * u(i,6,k) &
2890 + dy(j,7) * u(i,7,k) &
2891 + dy(j,8) * u(i,8,k)
2892 end do
2893 end do
2894 end do
2895
2896 do k = 1, lx
2897 do i = 1, lx*lx
2898 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2899 + dz(k,2) * u(i,1,2) &
2900 + dz(k,3) * u(i,1,3) &
2901 + dz(k,4) * u(i,1,4) &
2902 + dz(k,5) * u(i,1,5) &
2903 + dz(k,6) * u(i,1,6) &
2904 + dz(k,7) * u(i,1,7) &
2905 + dz(k,8) * u(i,1,8)
2906 end do
2907 end do
2908
2909 do i = 1, lx * lx * lx
2910 ux(i,1,1) = w3(i,1,1) &
2911 * ( drdx(i,1,1) * ur(i,1,1) &
2912 + dsdx(i,1,1) * us(i,1,1) &
2913 + dtdx(i,1,1) * ut(i,1,1) )
2914 uy(i,1,1) = w3(i,1,1) &
2915 * ( dsdy(i,1,1) * us(i,1,1) &
2916 + drdy(i,1,1) * ur(i,1,1) &
2917 + dtdy(i,1,1) * ut(i,1,1) )
2918 uz(i,1,1) = w3(i,1,1) &
2919 * ( dtdz(i,1,1) * ut(i,1,1) &
2920 + drdz(i,1,1) * ur(i,1,1) &
2921 + dsdz(i,1,1) * us(i,1,1) )
2922 end do
2923 end subroutine cpu_opgrad_lx8_single
2924
2925 !OCL SERIAL
2926 subroutine cpu_opgrad_lx7_single(ux, uy, uz, u, dx, dy, dz, &
2927 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2928 integer, parameter :: lx = 7
2929 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
2930 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
2931 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
2932 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
2933 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
2934 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
2935 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
2936 real(kind=rp) :: ur(lx, lx, lx)
2937 real(kind=rp) :: us(lx, lx, lx)
2938 real(kind=rp) :: ut(lx, lx, lx)
2939 integer :: i, j, k
2940
2941 do j = 1, lx * lx
2942 do i = 1, lx
2943 ur(i,j,1) = dx(i,1) * u(1,j,1) &
2944 + dx(i,2) * u(2,j,1) &
2945 + dx(i,3) * u(3,j,1) &
2946 + dx(i,4) * u(4,j,1) &
2947 + dx(i,5) * u(5,j,1) &
2948 + dx(i,6) * u(6,j,1) &
2949 + dx(i,7) * u(7,j,1)
2950 end do
2951 end do
2952
2953 do k = 1, lx
2954 do j = 1, lx
2955 do i = 1, lx
2956 us(i,j,k) = dy(j,1) * u(i,1,k) &
2957 + dy(j,2) * u(i,2,k) &
2958 + dy(j,3) * u(i,3,k) &
2959 + dy(j,4) * u(i,4,k) &
2960 + dy(j,5) * u(i,5,k) &
2961 + dy(j,6) * u(i,6,k) &
2962 + dy(j,7) * u(i,7,k)
2963 end do
2964 end do
2965 end do
2966
2967 do k = 1, lx
2968 do i = 1, lx*lx
2969 ut(i,1,k) = dz(k,1) * u(i,1,1) &
2970 + dz(k,2) * u(i,1,2) &
2971 + dz(k,3) * u(i,1,3) &
2972 + dz(k,4) * u(i,1,4) &
2973 + dz(k,5) * u(i,1,5) &
2974 + dz(k,6) * u(i,1,6) &
2975 + dz(k,7) * u(i,1,7)
2976 end do
2977 end do
2978
2979 do i = 1, lx * lx * lx
2980 ux(i,1,1) = w3(i,1,1) &
2981 * ( drdx(i,1,1) * ur(i,1,1) &
2982 + dsdx(i,1,1) * us(i,1,1) &
2983 + dtdx(i,1,1) * ut(i,1,1) )
2984 uy(i,1,1) = w3(i,1,1) &
2985 * ( dsdy(i,1,1) * us(i,1,1) &
2986 + drdy(i,1,1) * ur(i,1,1) &
2987 + dtdy(i,1,1) * ut(i,1,1) )
2988 uz(i,1,1) = w3(i,1,1) &
2989 * ( dtdz(i,1,1) * ut(i,1,1) &
2990 + drdz(i,1,1) * ur(i,1,1) &
2991 + dsdz(i,1,1) * us(i,1,1) )
2992 end do
2993 end subroutine cpu_opgrad_lx7_single
2994
2995 !OCL SERIAL
2996 subroutine cpu_opgrad_lx6_single(ux, uy, uz, u, dx, dy, dz, &
2997 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
2998 integer, parameter :: lx = 6
2999 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
3000 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
3001 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
3002 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
3003 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
3004 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
3005 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
3006 real(kind=rp) :: ur(lx, lx, lx)
3007 real(kind=rp) :: us(lx, lx, lx)
3008 real(kind=rp) :: ut(lx, lx, lx)
3009 integer :: i, j, k
3010
3011 do j = 1, lx * lx
3012 do i = 1, lx
3013 ur(i,j,1) = dx(i,1) * u(1,j,1) &
3014 + dx(i,2) * u(2,j,1) &
3015 + dx(i,3) * u(3,j,1) &
3016 + dx(i,4) * u(4,j,1) &
3017 + dx(i,5) * u(5,j,1) &
3018 + dx(i,6) * u(6,j,1)
3019 end do
3020 end do
3021
3022 do k = 1, lx
3023 do j = 1, lx
3024 do i = 1, lx
3025 us(i,j,k) = dy(j,1) * u(i,1,k) &
3026 + dy(j,2) * u(i,2,k) &
3027 + dy(j,3) * u(i,3,k) &
3028 + dy(j,4) * u(i,4,k) &
3029 + dy(j,5) * u(i,5,k) &
3030 + dy(j,6) * u(i,6,k)
3031 end do
3032 end do
3033 end do
3034
3035 do k = 1, lx
3036 do i = 1, lx*lx
3037 ut(i,1,k) = dz(k,1) * u(i,1,1) &
3038 + dz(k,2) * u(i,1,2) &
3039 + dz(k,3) * u(i,1,3) &
3040 + dz(k,4) * u(i,1,4) &
3041 + dz(k,5) * u(i,1,5) &
3042 + dz(k,6) * u(i,1,6)
3043 end do
3044 end do
3045
3046 do i = 1, lx * lx * lx
3047 ux(i,1,1) = w3(i,1,1) &
3048 * ( drdx(i,1,1) * ur(i,1,1) &
3049 + dsdx(i,1,1) * us(i,1,1) &
3050 + dtdx(i,1,1) * ut(i,1,1) )
3051 uy(i,1,1) = w3(i,1,1) &
3052 * ( dsdy(i,1,1) * us(i,1,1) &
3053 + drdy(i,1,1) * ur(i,1,1) &
3054 + dtdy(i,1,1) * ut(i,1,1) )
3055 uz(i,1,1) = w3(i,1,1) &
3056 * ( dtdz(i,1,1) * ut(i,1,1) &
3057 + drdz(i,1,1) * ur(i,1,1) &
3058 + dsdz(i,1,1) * us(i,1,1) )
3059 end do
3060 end subroutine cpu_opgrad_lx6_single
3061
3062 !OCL SERIAL
3063 subroutine cpu_opgrad_lx5_single(ux, uy, uz, u, dx, dy, dz, &
3064 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
3065 integer, parameter :: lx = 5
3066 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
3067 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
3068 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
3069 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
3070 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
3071 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
3072 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
3073 real(kind=rp) :: ur(lx, lx, lx)
3074 real(kind=rp) :: us(lx, lx, lx)
3075 real(kind=rp) :: ut(lx, lx, lx)
3076 integer :: i, j, k
3077
3078 do j = 1, lx * lx
3079 do i = 1, lx
3080 ur(i,j,1) = dx(i,1) * u(1,j,1) &
3081 + dx(i,2) * u(2,j,1) &
3082 + dx(i,3) * u(3,j,1) &
3083 + dx(i,4) * u(4,j,1) &
3084 + dx(i,5) * u(5,j,1)
3085 end do
3086 end do
3087
3088 do k = 1, lx
3089 do j = 1, lx
3090 do i = 1, lx
3091 us(i,j,k) = dy(j,1) * u(i,1,k) &
3092 + dy(j,2) * u(i,2,k) &
3093 + dy(j,3) * u(i,3,k) &
3094 + dy(j,4) * u(i,4,k) &
3095 + dy(j,5) * u(i,5,k)
3096 end do
3097 end do
3098 end do
3099
3100 do k = 1, lx
3101 do i = 1, lx*lx
3102 ut(i,1,k) = dz(k,1) * u(i,1,1) &
3103 + dz(k,2) * u(i,1,2) &
3104 + dz(k,3) * u(i,1,3) &
3105 + dz(k,4) * u(i,1,4) &
3106 + dz(k,5) * u(i,1,5)
3107 end do
3108 end do
3109
3110 do i = 1, lx * lx * lx
3111 ux(i,1,1) = w3(i,1,1) &
3112 * ( drdx(i,1,1) * ur(i,1,1) &
3113 + dsdx(i,1,1) * us(i,1,1) &
3114 + dtdx(i,1,1) * ut(i,1,1) )
3115 uy(i,1,1) = w3(i,1,1) &
3116 * ( dsdy(i,1,1) * us(i,1,1) &
3117 + drdy(i,1,1) * ur(i,1,1) &
3118 + dtdy(i,1,1) * ut(i,1,1) )
3119 uz(i,1,1) = w3(i,1,1) &
3120 * ( dtdz(i,1,1) * ut(i,1,1) &
3121 + drdz(i,1,1) * ur(i,1,1) &
3122 + dsdz(i,1,1) * us(i,1,1) )
3123 end do
3124 end subroutine cpu_opgrad_lx5_single
3125
3126 !OCL SERIAL
3127 subroutine cpu_opgrad_lx4_single(ux, uy, uz, u, dx, dy, dz, &
3128 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
3129 integer, parameter :: lx = 4
3130 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
3131 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
3132 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
3133 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
3134 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
3135 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
3136 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
3137 real(kind=rp) :: ur(lx, lx, lx)
3138 real(kind=rp) :: us(lx, lx, lx)
3139 real(kind=rp) :: ut(lx, lx, lx)
3140 integer :: i, j, k
3141
3142 do j = 1, lx * lx
3143 do i = 1, lx
3144 ur(i,j,1) = dx(i,1) * u(1,j,1) &
3145 + dx(i,2) * u(2,j,1) &
3146 + dx(i,3) * u(3,j,1) &
3147 + dx(i,4) * u(4,j,1)
3148 end do
3149 end do
3150
3151 do k = 1, lx
3152 do j = 1, lx
3153 do i = 1, lx
3154 us(i,j,k) = dy(j,1) * u(i,1,k) &
3155 + dy(j,2) * u(i,2,k) &
3156 + dy(j,3) * u(i,3,k) &
3157 + dy(j,4) * u(i,4,k)
3158 end do
3159 end do
3160 end do
3161
3162 do k = 1, lx
3163 do i = 1, lx*lx
3164 ut(i,1,k) = dz(k,1) * u(i,1,1) &
3165 + dz(k,2) * u(i,1,2) &
3166 + dz(k,3) * u(i,1,3) &
3167 + dz(k,4) * u(i,1,4)
3168 end do
3169 end do
3170
3171 do i = 1, lx * lx * lx
3172 ux(i,1,1) = w3(i,1,1) &
3173 * ( drdx(i,1,1) * ur(i,1,1) &
3174 + dsdx(i,1,1) * us(i,1,1) &
3175 + dtdx(i,1,1) * ut(i,1,1) )
3176 uy(i,1,1) = w3(i,1,1) &
3177 * ( dsdy(i,1,1) * us(i,1,1) &
3178 + drdy(i,1,1) * ur(i,1,1) &
3179 + dtdy(i,1,1) * ut(i,1,1) )
3180 uz(i,1,1) = w3(i,1,1) &
3181 * ( dtdz(i,1,1) * ut(i,1,1) &
3182 + drdz(i,1,1) * ur(i,1,1) &
3183 + dsdz(i,1,1) * us(i,1,1) )
3184 end do
3185 end subroutine cpu_opgrad_lx4_single
3186
3187 !OCL SERIAL
3188 subroutine cpu_opgrad_lx3_single(ux, uy, uz, u, dx, dy, dz, &
3189 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
3190 integer, parameter :: lx = 3
3191 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
3192 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
3193 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
3194 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
3195 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
3196 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
3197 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
3198 real(kind=rp) :: ur(lx, lx, lx)
3199 real(kind=rp) :: us(lx, lx, lx)
3200 real(kind=rp) :: ut(lx, lx, lx)
3201 integer :: i, j, k
3202
3203 do j = 1, lx * lx
3204 do i = 1, lx
3205 ur(i,j,1) = dx(i,1) * u(1,j,1) &
3206 + dx(i,2) * u(2,j,1) &
3207 + dx(i,3) * u(3,j,1)
3208 end do
3209 end do
3210
3211 do k = 1, lx
3212 do j = 1, lx
3213 do i = 1, lx
3214 us(i,j,k) = dy(j,1) * u(i,1,k) &
3215 + dy(j,2) * u(i,2,k) &
3216 + dy(j,3) * u(i,3,k)
3217 end do
3218 end do
3219 end do
3220
3221 do k = 1, lx
3222 do i = 1, lx*lx
3223 ut(i,1,k) = dz(k,1) * u(i,1,1) &
3224 + dz(k,2) * u(i,1,2) &
3225 + dz(k,3) * u(i,1,3)
3226 end do
3227 end do
3228
3229 do i = 1, lx * lx * lx
3230 ux(i,1,1) = w3(i,1,1) &
3231 * ( drdx(i,1,1) * ur(i,1,1) &
3232 + dsdx(i,1,1) * us(i,1,1) &
3233 + dtdx(i,1,1) * ut(i,1,1) )
3234 uy(i,1,1) = w3(i,1,1) &
3235 * ( dsdy(i,1,1) * us(i,1,1) &
3236 + drdy(i,1,1) * ur(i,1,1) &
3237 + dtdy(i,1,1) * ut(i,1,1) )
3238 uz(i,1,1) = w3(i,1,1) &
3239 * ( dtdz(i,1,1) * ut(i,1,1) &
3240 + drdz(i,1,1) * ur(i,1,1) &
3241 + dsdz(i,1,1) * us(i,1,1) )
3242 end do
3243 end subroutine cpu_opgrad_lx3_single
3244
3245 !OCL SERIAL
3246 subroutine cpu_opgrad_lx2_single(ux, uy, uz, u, dx, dy, dz, &
3247 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3)
3248 integer, parameter :: lx = 2
3249 real(kind=rp), dimension(lx, lx, lx), intent(inout) :: ux, uy, uz
3250 real(kind=rp), dimension(lx, lx, lx), intent(in) :: u
3251 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
3252 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdx, dsdx, dtdx
3253 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdy, dsdy, dtdy
3254 real(kind=rp), dimension(lx, lx, lx), intent(in) :: drdz, dsdz, dtdz
3255 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
3256 real(kind=rp) :: ur(lx, lx, lx)
3257 real(kind=rp) :: us(lx, lx, lx)
3258 real(kind=rp) :: ut(lx, lx, lx)
3259 integer :: i, j, k
3260
3261 do j = 1, lx * lx
3262 do i = 1, lx
3263 ur(i,j,1) = dx(i,1) * u(1,j,1) &
3264 + dx(i,2) * u(2,j,1)
3265 end do
3266 end do
3267
3268 do k = 1, lx
3269 do j = 1, lx
3270 do i = 1, lx
3271 us(i,j,k) = dy(j,1) * u(i,1,k) &
3272 + dy(j,2) * u(i,2,k)
3273 end do
3274 end do
3275 end do
3276
3277 do k = 1, lx
3278 do i = 1, lx*lx
3279 ut(i,1,k) = dz(k,1) * u(i,1,1) &
3280 + dz(k,2) * u(i,1,2)
3281 end do
3282 end do
3283
3284 do i = 1, lx * lx * lx
3285 ux(i,1,1) = w3(i,1,1) &
3286 * ( drdx(i,1,1) * ur(i,1,1) &
3287 + dsdx(i,1,1) * us(i,1,1) &
3288 + dtdx(i,1,1) * ut(i,1,1) )
3289 uy(i,1,1) = w3(i,1,1) &
3290 * ( dsdy(i,1,1) * us(i,1,1) &
3291 + drdy(i,1,1) * ur(i,1,1) &
3292 + dtdy(i,1,1) * ut(i,1,1) )
3293 uz(i,1,1) = w3(i,1,1) &
3294 * ( dtdz(i,1,1) * ut(i,1,1) &
3295 + drdz(i,1,1) * ur(i,1,1) &
3296 + dsdz(i,1,1) * us(i,1,1) )
3297 end do
3298
3299 end subroutine cpu_opgrad_lx2_single
3300
3301end submodule cpu_opgrad
Operators CPU backend.
Definition opr_cpu.f90:34