Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
sx_opgrad.f90
Go to the documentation of this file.
1! Copyright (c) 2021, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
34submodule(opr_sx) sx_opgrad
35 implicit none
36
37contains
38
39 module subroutine opr_sx_opgrad(ux, uy, uz, u, coef)
40 type(coef_t), intent(in) :: coef
41 real(kind=rp), intent(inout) :: ux(coef%Xh%lxyz, coef%msh%nelv)
42 real(kind=rp), intent(inout) :: uy(coef%Xh%lxyz, coef%msh%nelv)
43 real(kind=rp), intent(inout) :: uz(coef%Xh%lxyz, coef%msh%nelv)
44 real(kind=rp), intent(in) :: u(coef%Xh%lxyz, coef%msh%nelv)
45
46 associate(xh => coef%Xh, msh => coef%msh)
47 select case (xh%lx)
48 case (18)
49 call sx_opgrad_lx18(ux, uy, uz, u, &
50 xh%dx, xh%dy, xh%dz, &
51 coef%drdx, coef%dsdx, coef%dtdx, &
52 coef%drdy, coef%dsdy, coef%dtdy, &
53 coef%drdz, coef%dsdz, coef%dtdz, &
54 xh%w3, msh%nelv)
55 case (17)
56 call sx_opgrad_lx17(ux, uy, uz, u, &
57 xh%dx, xh%dy, xh%dz, &
58 coef%drdx, coef%dsdx, coef%dtdx, &
59 coef%drdy, coef%dsdy, coef%dtdy, &
60 coef%drdz, coef%dsdz, coef%dtdz, &
61 xh%w3, msh%nelv)
62 case (16)
63 call sx_opgrad_lx16(ux, uy, uz, u, &
64 xh%dx, xh%dy, xh%dz, &
65 coef%drdx, coef%dsdx, coef%dtdx, &
66 coef%drdy, coef%dsdy, coef%dtdy, &
67 coef%drdz, coef%dsdz, coef%dtdz, &
68 xh%w3, msh%nelv)
69 case (15)
70 call sx_opgrad_lx15(ux, uy, uz, u, &
71 xh%dx, xh%dy, xh%dz, &
72 coef%drdx, coef%dsdx, coef%dtdx, &
73 coef%drdy, coef%dsdy, coef%dtdy, &
74 coef%drdz, coef%dsdz, coef%dtdz, &
75 xh%w3, msh%nelv)
76 case (14)
77 call sx_opgrad_lx14(ux, uy, uz, u, &
78 xh%dx, xh%dy, xh%dz, &
79 coef%drdx, coef%dsdx, coef%dtdx, &
80 coef%drdy, coef%dsdy, coef%dtdy, &
81 coef%drdz, coef%dsdz, coef%dtdz, &
82 xh%w3, msh%nelv)
83 case (13)
84 call sx_opgrad_lx13(ux, uy, uz, u, &
85 xh%dx, xh%dy, xh%dz, &
86 coef%drdx, coef%dsdx, coef%dtdx, &
87 coef%drdy, coef%dsdy, coef%dtdy, &
88 coef%drdz, coef%dsdz, coef%dtdz, &
89 xh%w3, msh%nelv)
90 case (12)
91 call sx_opgrad_lx12(ux, uy, uz, u, &
92 xh%dx, xh%dy, xh%dz, &
93 coef%drdx, coef%dsdx, coef%dtdx, &
94 coef%drdy, coef%dsdy, coef%dtdy, &
95 coef%drdz, coef%dsdz, coef%dtdz, &
96 xh%w3, msh%nelv)
97 case (11)
98 call sx_opgrad_lx11(ux, uy, uz, u, &
99 xh%dx, xh%dy, xh%dz, &
100 coef%drdx, coef%dsdx, coef%dtdx, &
101 coef%drdy, coef%dsdy, coef%dtdy, &
102 coef%drdz, coef%dsdz, coef%dtdz, &
103 xh%w3, msh%nelv)
104 case (10)
105 call sx_opgrad_lx10(ux, uy, uz, u, &
106 xh%dx, xh%dy, xh%dz, &
107 coef%drdx, coef%dsdx, coef%dtdx, &
108 coef%drdy, coef%dsdy, coef%dtdy, &
109 coef%drdz, coef%dsdz, coef%dtdz, &
110 xh%w3, msh%nelv)
111 case (9)
112 call sx_opgrad_lx9(ux, uy, uz, u, &
113 xh%dx, xh%dy, xh%dz, &
114 coef%drdx, coef%dsdx, coef%dtdx, &
115 coef%drdy, coef%dsdy, coef%dtdy, &
116 coef%drdz, coef%dsdz, coef%dtdz, &
117 xh%w3, msh%nelv)
118 case (8)
119 call sx_opgrad_lx8(ux, uy, uz, u, &
120 xh%dx, xh%dy, xh%dz, &
121 coef%drdx, coef%dsdx, coef%dtdx, &
122 coef%drdy, coef%dsdy, coef%dtdy, &
123 coef%drdz, coef%dsdz, coef%dtdz, &
124 xh%w3, msh%nelv)
125 case (7)
126 call sx_opgrad_lx7(ux, uy, uz, u, &
127 xh%dx, xh%dy, xh%dz, &
128 coef%drdx, coef%dsdx, coef%dtdx, &
129 coef%drdy, coef%dsdy, coef%dtdy, &
130 coef%drdz, coef%dsdz, coef%dtdz, &
131 xh%w3, msh%nelv)
132 case (6)
133 call sx_opgrad_lx6(ux, uy, uz, u, &
134 xh%dx, xh%dy, xh%dz, &
135 coef%drdx, coef%dsdx, coef%dtdx, &
136 coef%drdy, coef%dsdy, coef%dtdy, &
137 coef%drdz, coef%dsdz, coef%dtdz, &
138 xh%w3, msh%nelv)
139 case (5)
140 call sx_opgrad_lx5(ux, uy, uz, u, &
141 xh%dx, xh%dy, xh%dz, &
142 coef%drdx, coef%dsdx, coef%dtdx, &
143 coef%drdy, coef%dsdy, coef%dtdy, &
144 coef%drdz, coef%dsdz, coef%dtdz, &
145 xh%w3, msh%nelv)
146 case (4)
147 call sx_opgrad_lx4(ux, uy, uz, u, &
148 xh%dx, xh%dy, xh%dz, &
149 coef%drdx, coef%dsdx, coef%dtdx, &
150 coef%drdy, coef%dsdy, coef%dtdy, &
151 coef%drdz, coef%dsdz, coef%dtdz, &
152 xh%w3, msh%nelv)
153 case (3)
154 call sx_opgrad_lx3(ux, uy, uz, u, &
155 xh%dx, xh%dy, xh%dz, &
156 coef%drdx, coef%dsdx, coef%dtdx, &
157 coef%drdy, coef%dsdy, coef%dtdy, &
158 coef%drdz, coef%dsdz, coef%dtdz, &
159 xh%w3, msh%nelv)
160 case (2)
161 call sx_opgrad_lx2(ux, uy, uz, u, &
162 xh%dx, xh%dy, xh%dz, &
163 coef%drdx, coef%dsdx, coef%dtdx, &
164 coef%drdy, coef%dsdy, coef%dtdy, &
165 coef%drdz, coef%dsdz, coef%dtdz, &
166 xh%w3, msh%nelv)
167 case default
168 call sx_opgrad_lx(ux, uy, uz, u, &
169 xh%dx, xh%dy, xh%dz, &
170 coef%drdx, coef%dsdx, coef%dtdx, &
171 coef%drdy, coef%dsdy, coef%dtdy, &
172 coef%drdz, coef%dsdz, coef%dtdz, &
173 xh%w3, msh%nelv, xh%lx)
174 end select
175 end associate
176
177 end subroutine opr_sx_opgrad
178
179 subroutine sx_opgrad_lx(ux, uy, uz, u, dx, dy, dz, &
180 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n, lx)
181 integer, intent(in) :: n, lx
182 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
183 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
184 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
185 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
186 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
187 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
188 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
189 real(kind=rp) :: ur(lx, lx, lx, n)
190 real(kind=rp) :: us(lx, lx, lx, n)
191 real(kind=rp) :: ut(lx, lx, lx, n)
192 real(kind=rp) :: wr, ws, wt
193 integer :: e, i, j, k, jj, kk
194
195 do i = 1, lx
196 do jj = 1, lx * lx * n
197 wr = 0d0
198 do kk = 1, lx
199 wr = wr + dx(i, kk) * u(kk, jj,1,1)
200 end do
201 ur(i, jj,1,1) = wr
202 end do
203 end do
204
205 do k = 1, lx
206 do i = 1, lx
207 do j = 1, lx
208 do e = 1, n
209 ws = 0d0
210 !NEC$ unroll_completely
211 do kk = 1, lx
212 ws = ws + dy(j, kk) * u(i, kk, k,e)
213 end do
214 us(i, j, k,e) = ws
215 end do
216 end do
217 end do
218 end do
219
220 do j = 1, lx
221 do i = 1, lx
222 do k = 1, lx
223 do e = 1, n
224 wt = 0d0
225 !NEC$ unroll_completely
226 do kk = 1, lx
227 wt = wt + dz(k, kk) * u(i, j, kk, e)
228 end do
229 ut(i, j, k, e) = wt
230 end do
231 end do
232 end do
233 end do
234
235 do i = 1, lx * lx * lx
236 do e = 1, n
237 ux(i,1,1,e) = w3(i,1,1) * &
238 ( drdx(i,1,1,e) * ur(i,1,1,e) &
239 + dsdx(i,1,1,e) * us(i,1,1,e) &
240 + dtdx(i,1,1,e) * ut(i,1,1,e))
241
242 uy(i,1,1,e) = w3(i,1,1) * &
243 ( dsdy(i,1,1,e) * us(i,1,1,e) &
244 + drdy(i,1,1,e) * ur(i,1,1,e) &
245 + dtdy(i,1,1,e) * ut(i,1,1,e) )
246
247 uz(i,1,1,e) = w3(i,1,1) * &
248 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
249 + drdz(i,1,1,e) * ur(i,1,1,e) &
250 + dsdz(i,1,1,e) * us(i,1,1,e))
251 end do
252 end do
253 end subroutine sx_opgrad_lx
254
255 subroutine sx_opgrad_lx18(ux, uy, uz, u, dx, dy, dz, &
256 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
257 integer, parameter :: lx = 18
258 integer, intent(in) :: n
259 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
260 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
261 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
262 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
263 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
264 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
265 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
266 real(kind=rp) :: ur(lx, lx, lx, n)
267 real(kind=rp) :: us(lx, lx, lx, n)
268 real(kind=rp) :: ut(lx, lx, lx, n)
269 real(kind=rp) :: wr, ws, wt
270 integer :: e, i, j, k, jj, kk
271
272 do i = 1, lx
273 do jj = 1, lx * lx * n
274 wr = 0d0
275 do kk = 1, lx
276 wr = wr + dx(i, kk) * u(kk, jj,1,1)
277 end do
278 ur(i, jj,1,1) = wr
279 end do
280 end do
281
282 do k = 1, lx
283 do i = 1, lx
284 do j = 1, lx
285 do e = 1, n
286 ws = 0d0
287 !NEC$ unroll_completely
288 do kk = 1, lx
289 ws = ws + dy(j, kk) * u(i, kk, k,e)
290 end do
291 us(i, j, k,e) = ws
292 end do
293 end do
294 end do
295 end do
296
297 do j = 1, lx
298 do i = 1, lx
299 do k = 1, lx
300 do e = 1, n
301 wt = 0d0
302 !NEC$ unroll_completely
303 do kk = 1, lx
304 wt = wt + dz(k, kk) * u(i, j, kk,e)
305 end do
306 ut(i, j, k,e) = wt
307 end do
308 end do
309 end do
310 end do
311
312 do i = 1, lx * lx * lx
313 do e = 1, n
314 ux(i,1,1,e) = w3(i,1,1) * &
315 ( drdx(i,1,1,e) * ur(i,1,1,e) &
316 + dsdx(i,1,1,e) * us(i,1,1,e) &
317 + dtdx(i,1,1,e) * ut(i,1,1,e))
318
319 uy(i,1,1,e) = w3(i,1,1) * &
320 ( dsdy(i,1,1,e) * us(i,1,1,e) &
321 + drdy(i,1,1,e) * ur(i,1,1,e) &
322 + dtdy(i,1,1,e) * ut(i,1,1,e) )
323
324 uz(i,1,1,e) = w3(i,1,1) * &
325 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
326 + drdz(i,1,1,e) * ur(i,1,1,e) &
327 + dsdz(i,1,1,e) * us(i,1,1,e))
328 end do
329 end do
330 end subroutine sx_opgrad_lx18
331
332 subroutine sx_opgrad_lx17(ux, uy, uz, u, dx, dy, dz, &
333 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
334 integer, parameter :: lx = 17
335 integer, intent(in) :: n
336 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
337 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
338 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
339 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
340 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
341 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
342 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
343 real(kind=rp) :: ur(lx, lx, lx, n)
344 real(kind=rp) :: us(lx, lx, lx, n)
345 real(kind=rp) :: ut(lx, lx, lx, n)
346 real(kind=rp) :: wr, ws, wt
347 integer :: e, i, j, k, jj, kk
348
349 do i = 1, lx
350 do jj = 1, lx * lx * n
351 wr = 0d0
352 do kk = 1, lx
353 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
354 end do
355 ur(i, jj, 1, 1) = wr
356 end do
357 end do
358
359 do k = 1, lx
360 do i = 1, lx
361 do j = 1, lx
362 do e = 1, n
363 ws = 0d0
364 !NEC$ unroll_completely
365 do kk = 1, lx
366 ws = ws + dy(j, kk) * u(i, kk, k,e)
367 end do
368 us(i, j, k,e) = ws
369 end do
370 end do
371 end do
372 end do
373
374 do j = 1, lx
375 do i = 1, lx
376 do k = 1, lx
377 do e = 1, n
378 wt = 0d0
379 !NEC$ unroll_completely
380 do kk = 1, lx
381 wt = wt + dz(k, kk) * u(i, j, kk, e)
382 end do
383 ut(i, j, k, e) = wt
384 end do
385 end do
386 end do
387 end do
388
389 do i = 1, lx * lx * lx
390 do e = 1, n
391 ux(i,1,1,e) = w3(i,1,1) * &
392 ( drdx(i,1,1,e) * ur(i,1,1,e) &
393 + dsdx(i,1,1,e) * us(i,1,1,e) &
394 + dtdx(i,1,1,e) * ut(i,1,1,e))
395
396 uy(i,1,1,e) = w3(i,1,1) * &
397 ( dsdy(i,1,1,e) * us(i,1,1,e) &
398 + drdy(i,1,1,e) * ur(i,1,1,e) &
399 + dtdy(i,1,1,e) * ut(i,1,1,e) )
400
401 uz(i,1,1,e) = w3(i,1,1) * &
402 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
403 + drdz(i,1,1,e) * ur(i,1,1,e) &
404 + dsdz(i,1,1,e) * us(i,1,1,e))
405 end do
406 end do
407 end subroutine sx_opgrad_lx17
408
409 subroutine sx_opgrad_lx16(ux, uy, uz, u, dx, dy, dz, &
410 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
411 integer, parameter :: lx = 16
412 integer, intent(in) :: n
413 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
414 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
415 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
416 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
417 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
418 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
419 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
420 real(kind=rp) :: ur(lx, lx, lx, n)
421 real(kind=rp) :: us(lx, lx, lx, n)
422 real(kind=rp) :: ut(lx, lx, lx, n)
423 real(kind=rp) :: wr, ws, wt
424 integer :: e, i, j, k, jj, kk
425
426 do i = 1, lx
427 do jj = 1, lx * lx * n
428 wr = 0d0
429 do kk = 1, lx
430 wr = wr + dx(i, kk) * u(kk, jj,1,1)
431 end do
432 ur(i, jj,1,1) = wr
433 end do
434 end do
435
436 do k = 1, lx
437 do i = 1, lx
438 do j = 1, lx
439 do e = 1, n
440 ws = 0d0
441 !NEC$ unroll_completely
442 do kk = 1, lx
443 ws = ws + dy(j, kk) * u(i, kk, k,e)
444 end do
445 us(i, j, k,e) = ws
446 end do
447 end do
448 end do
449 end do
450
451 do j = 1, lx
452 do i = 1, lx
453 do k = 1, lx
454 do e = 1, n
455 wt = 0d0
456 !NEC$ unroll_completely
457 do kk = 1, lx
458 wt = wt + dz(k, kk) * u(i, j, kk,e)
459 end do
460 ut(i, j, k,e) = wt
461 end do
462 end do
463 end do
464 end do
465
466 do i = 1, lx * lx * lx
467 do e = 1, n
468 ux(i,1,1,e) = w3(i,1,1) * &
469 ( drdx(i,1,1,e) * ur(i,1,1,e) &
470 + dsdx(i,1,1,e) * us(i,1,1,e) &
471 + dtdx(i,1,1,e) * ut(i,1,1,e))
472
473 uy(i,1,1,e) = w3(i,1,1) * &
474 ( dsdy(i,1,1,e) * us(i,1,1,e) &
475 + drdy(i,1,1,e) * ur(i,1,1,e) &
476 + dtdy(i,1,1,e) * ut(i,1,1,e) )
477
478 uz(i,1,1,e) = w3(i,1,1) * &
479 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
480 + drdz(i,1,1,e) * ur(i,1,1,e) &
481 + dsdz(i,1,1,e) * us(i,1,1,e))
482 end do
483 end do
484 end subroutine sx_opgrad_lx16
485
486 subroutine sx_opgrad_lx15(ux, uy, uz, u, dx, dy, dz, &
487 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
488 integer, parameter :: lx = 15
489 integer, intent(in) :: n
490 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
491 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
492 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
493 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
494 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
495 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
496 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
497 real(kind=rp) :: ur(lx, lx, lx, n)
498 real(kind=rp) :: us(lx, lx, lx, n)
499 real(kind=rp) :: ut(lx, lx, lx, n)
500 real(kind=rp) :: wr, ws, wt
501 integer :: e, i, j, k, jj, kk
502
503 do i = 1, lx
504 do jj = 1, lx * lx * n
505 wr = 0d0
506 do kk = 1, lx
507 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
508 end do
509 ur(i, jj,1,1) = wr
510 end do
511 end do
512
513 do k = 1, lx
514 do i = 1, lx
515 do j = 1, lx
516 do e = 1, n
517 ws = 0d0
518 !NEC$ unroll_completely
519 do kk = 1, lx
520 ws = ws + dy(j, kk) * u(i, kk, k, e)
521 end do
522 us(i, j, k, e) = ws
523 end do
524 end do
525 end do
526 end do
527
528 do j = 1, lx
529 do i = 1, lx
530 do k = 1, lx
531 do e = 1, n
532 wt = 0d0
533 !NEC$ unroll_completely
534 do kk = 1, lx
535 wt = wt + dz(k, kk) * u(i, j, kk, e)
536 end do
537 ut(i, j, k, e) = wt
538 end do
539 end do
540 end do
541 end do
542
543 do i = 1, lx * lx * lx
544 do e = 1, n
545 ux(i,1,1,e) = w3(i,1,1) * &
546 ( drdx(i,1,1,e) * ur(i,1,1,e) &
547 + dsdx(i,1,1,e) * us(i,1,1,e) &
548 + dtdx(i,1,1,e) * ut(i,1,1,e))
549
550 uy(i,1,1,e) = w3(i,1,1) * &
551 ( dsdy(i,1,1,e) * us(i,1,1,e) &
552 + drdy(i,1,1,e) * ur(i,1,1,e) &
553 + dtdy(i,1,1,e) * ut(i,1,1,e) )
554
555 uz(i,1,1,e) = w3(i,1,1) * &
556 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
557 + drdz(i,1,1,e) * ur(i,1,1,e) &
558 + dsdz(i,1,1,e) * us(i,1,1,e))
559 end do
560 end do
561 end subroutine sx_opgrad_lx15
562
563 subroutine sx_opgrad_lx14(ux, uy, uz, u, dx, dy, dz, &
564 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
565 integer, parameter :: lx = 14
566 integer, intent(in) :: n
567 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
568 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
569 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
570 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
571 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
572 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
573 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
574 real(kind=rp) :: ur(lx, lx, lx, n)
575 real(kind=rp) :: us(lx, lx, lx, n)
576 real(kind=rp) :: ut(lx, lx, lx, n)
577 real(kind=rp) :: wr, ws, wt
578 integer :: e, i, j, k, jj, kk
579
580 do i = 1, lx
581 do jj = 1, lx * lx * n
582 wr = 0d0
583 do kk = 1, lx
584 wr = wr + dx(i, kk) * u(kk, jj,1,1)
585 end do
586 ur(i, jj,1,1) = wr
587 end do
588 end do
589
590 do k = 1, lx
591 do i = 1, lx
592 do j = 1, lx
593 do e = 1, n
594 ws = 0d0
595 !NEC$ unroll_completely
596 do kk = 1, lx
597 ws = ws + dy(j, kk) * u(i, kk, k, e)
598 end do
599 us(i, j, k,e) = ws
600 end do
601 end do
602 end do
603 end do
604
605 do j = 1, lx
606 do i = 1, lx
607 do k = 1, lx
608 do e = 1, n
609 wt = 0d0
610 !NEC$ unroll_completely
611 do kk = 1, lx
612 wt = wt + dz(k, kk) * u(i, j, kk,e)
613 end do
614 ut(i, j, k,e) = wt
615 end do
616 end do
617 end do
618 end do
619
620 do i = 1, lx * lx * lx
621 do e = 1, n
622 ux(i,1,1,e) = w3(i,1,1) * &
623 ( drdx(i,1,1,e) * ur(i,1,1,e) &
624 + dsdx(i,1,1,e) * us(i,1,1,e) &
625 + dtdx(i,1,1,e) * ut(i,1,1,e))
626
627 uy(i,1,1,e) = w3(i,1,1) * &
628 ( dsdy(i,1,1,e) * us(i,1,1,e) &
629 + drdy(i,1,1,e) * ur(i,1,1,e) &
630 + dtdy(i,1,1,e) * ut(i,1,1,e) )
631
632 uz(i,1,1,e) = w3(i,1,1) * &
633 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
634 + drdz(i,1,1,e) * ur(i,1,1,e) &
635 + dsdz(i,1,1,e) * us(i,1,1,e))
636 end do
637 end do
638 end subroutine sx_opgrad_lx14
639
640 subroutine sx_opgrad_lx13(ux, uy, uz, u, dx, dy, dz, &
641 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
642 integer, parameter :: lx = 13
643 integer, intent(in) :: n
644 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
645 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
646 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
647 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
648 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
649 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
650 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
651 real(kind=rp) :: ur(lx, lx, lx, n)
652 real(kind=rp) :: us(lx, lx, lx, n)
653 real(kind=rp) :: ut(lx, lx, lx, n)
654 real(kind=rp) :: wr, ws, wt
655 integer :: e, i, j, k, jj, kk
656
657 do i = 1, lx
658 do jj = 1, lx * lx * n
659 wr = 0d0
660 do kk = 1, lx
661 wr = wr + dx(i, kk) * u(kk, jj,1,1)
662 end do
663 ur(i, jj,1,1) = wr
664 end do
665 end do
666
667 do k = 1, lx
668 do i = 1, lx
669 do j = 1, lx
670 do e = 1, n
671 ws = 0d0
672 !NEC$ unroll_completely
673 do kk = 1, lx
674 ws = ws + dy(j, kk) * u(i, kk, k,e)
675 end do
676 us(i, j, k,e) = ws
677 end do
678 end do
679 end do
680 end do
681
682 do j = 1, lx
683 do i = 1, lx
684 do k = 1, lx
685 do e = 1, n
686 wt = 0d0
687 !NEC$ unroll_completely
688 do kk = 1, lx
689 wt = wt + dz(k, kk) * u(i, j, kk,e)
690 end do
691 ut(i, j, k,e) = wt
692 end do
693 end do
694 end do
695 end do
696
697 do i = 1, lx * lx * lx
698 do e = 1, n
699 ux(i,1,1,e) = w3(i,1,1) * &
700 ( drdx(i,1,1,e) * ur(i,1,1,e) &
701 + dsdx(i,1,1,e) * us(i,1,1,e) &
702 + dtdx(i,1,1,e) * ut(i,1,1,e))
703
704 uy(i,1,1,e) = w3(i,1,1) * &
705 ( dsdy(i,1,1,e) * us(i,1,1,e) &
706 + drdy(i,1,1,e) * ur(i,1,1,e) &
707 + dtdy(i,1,1,e) * ut(i,1,1,e) )
708
709 uz(i,1,1,e) = w3(i,1,1) * &
710 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
711 + drdz(i,1,1,e) * ur(i,1,1,e) &
712 + dsdz(i,1,1,e) * us(i,1,1,e))
713 end do
714 end do
715 end subroutine sx_opgrad_lx13
716
717 subroutine sx_opgrad_lx12(ux, uy, uz, u, dx, dy, dz, &
718 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
719 integer, parameter :: lx = 12
720 integer, intent(in) :: n
721 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
722 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
723 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
724 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
725 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
726 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
727 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
728 real(kind=rp) :: ur(lx, lx, lx, n)
729 real(kind=rp) :: us(lx, lx, lx, n)
730 real(kind=rp) :: ut(lx, lx, lx, n)
731 real(kind=rp) :: wr, ws, wt
732 integer :: e, i, j, k, jj, kk
733
734 do i = 1, lx
735 do jj = 1, lx * lx * n
736 wr = 0d0
737 do kk = 1, lx
738 wr = wr + dx(i, kk) * u(kk, jj,1,1)
739 end do
740 ur(i, jj,1,1) = wr
741 end do
742 end do
743
744 do k = 1, lx
745 do i = 1, lx
746 do j = 1, lx
747 do e = 1, n
748 ws = 0d0
749 !NEC$ unroll_completely
750 do kk = 1, lx
751 ws = ws + dy(j, kk) * u(i, kk, k,e)
752 end do
753 us(i, j, k,e) = ws
754 end do
755 end do
756 end do
757 end do
758
759 do j = 1, lx
760 do i = 1, lx
761 do k = 1, lx
762 do e = 1, n
763 wt = 0d0
764 !NEC$ unroll_completely
765 do kk = 1, lx
766 wt = wt + dz(k, kk) * u(i, j, kk,e)
767 end do
768 ut(i, j, k,e) = wt
769 end do
770 end do
771 end do
772 end do
773
774 do i = 1, lx * lx * lx
775 do e = 1, n
776 ux(i,1,1,e) = w3(i,1,1) * &
777 ( drdx(i,1,1,e) * ur(i,1,1,e) &
778 + dsdx(i,1,1,e) * us(i,1,1,e) &
779 + dtdx(i,1,1,e) * ut(i,1,1,e))
780
781 uy(i,1,1,e) = w3(i,1,1) * &
782 ( dsdy(i,1,1,e) * us(i,1,1,e) &
783 + drdy(i,1,1,e) * ur(i,1,1,e) &
784 + dtdy(i,1,1,e) * ut(i,1,1,e) )
785
786 uz(i,1,1,e) = w3(i,1,1) * &
787 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
788 + drdz(i,1,1,e) * ur(i,1,1,e) &
789 + dsdz(i,1,1,e) * us(i,1,1,e))
790 end do
791 end do
792 end subroutine sx_opgrad_lx12
793
794 subroutine sx_opgrad_lx11(ux, uy, uz, u, dx, dy, dz, &
795 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
796 integer, parameter :: lx = 11
797 integer, intent(in) :: n
798 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
799 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
800 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
801 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
802 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
803 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
804 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
805 real(kind=rp) :: ur(lx, lx, lx, n)
806 real(kind=rp) :: us(lx, lx, lx, n)
807 real(kind=rp) :: ut(lx, lx, lx, n)
808 real(kind=rp) :: wr, ws, wt
809 integer :: e, i, j, k, jj, kk
810
811 do i = 1, lx
812 do jj = 1, lx * lx * n
813 wr = 0d0
814 do kk = 1, lx
815 wr = wr + dx(i, kk) * u(kk, jj,1,1)
816 end do
817 ur(i, jj,1,1) = wr
818 end do
819 end do
820
821 do k = 1, lx
822 do i = 1, lx
823 do j = 1, lx
824 do e = 1, n
825 ws = 0d0
826 !NEC$ unroll_completely
827 do kk = 1, lx
828 ws = ws + dy(j, kk) * u(i, kk, k,e)
829 end do
830 us(i, j, k,e) = ws
831 end do
832 end do
833 end do
834 end do
835
836 do j = 1, lx
837 do i = 1, lx
838 do k = 1, lx
839 do e = 1, n
840 wt = 0d0
841 !NEC$ unroll_completely
842 do kk = 1, lx
843 wt = wt + dz(k, kk) * u(i, j, kk,e)
844 end do
845 ut(i, j, k,e) = wt
846 end do
847 end do
848 end do
849 end do
850
851 do i = 1, lx * lx * lx
852 do e = 1, n
853 ux(i,1,1,e) = w3(i,1,1) * &
854 ( drdx(i,1,1,e) * ur(i,1,1,e) &
855 + dsdx(i,1,1,e) * us(i,1,1,e) &
856 + dtdx(i,1,1,e) * ut(i,1,1,e))
857
858 uy(i,1,1,e) = w3(i,1,1) * &
859 ( dsdy(i,1,1,e) * us(i,1,1,e) &
860 + drdy(i,1,1,e) * ur(i,1,1,e) &
861 + dtdy(i,1,1,e) * ut(i,1,1,e) )
862
863 uz(i,1,1,e) = w3(i,1,1) * &
864 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
865 + drdz(i,1,1,e) * ur(i,1,1,e) &
866 + dsdz(i,1,1,e) * us(i,1,1,e))
867 end do
868 end do
869 end subroutine sx_opgrad_lx11
870
871 subroutine sx_opgrad_lx10(ux, uy, uz, u, dx, dy, dz, &
872 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
873 integer, parameter :: lx = 10
874 integer, intent(in) :: n
875 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
876 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
877 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
878 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
879 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
880 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
881 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
882 real(kind=rp) :: ur(lx, lx, lx, n)
883 real(kind=rp) :: us(lx, lx, lx, n)
884 real(kind=rp) :: ut(lx, lx, lx, n)
885 real(kind=rp) :: wr, ws, wt
886 integer :: e, i, j, k, jj, kk
887
888 do i = 1, lx
889 do jj = 1, lx * lx * n
890 wr = 0d0
891 do kk = 1, lx
892 wr = wr + dx(i, kk) * u(kk, jj,1,1)
893 end do
894 ur(i, jj,1,1) = wr
895 end do
896 end do
897
898 do k = 1, lx
899 do i = 1, lx
900 do j = 1, lx
901 do e = 1, n
902 ws = 0d0
903 !NEC$ unroll_completely
904 do kk = 1, lx
905 ws = ws + dy(j, kk) * u(i, kk, k,e)
906 end do
907 us(i, j, k,e) = ws
908 end do
909 end do
910 end do
911 end do
912
913 do j = 1, lx
914 do i = 1, lx
915 do k = 1, lx
916 do e = 1, n
917 wt = 0d0
918 !NEC$ unroll_completely
919 do kk = 1, lx
920 wt = wt + dz(k, kk) * u(i, j, kk,e)
921 end do
922 ut(i, j, k,e) = wt
923 end do
924 end do
925 end do
926 end do
927
928 do i = 1, lx * lx * lx
929 do e = 1, n
930 ux(i,1,1,e) = w3(i,1,1) * &
931 ( drdx(i,1,1,e) * ur(i,1,1,e) &
932 + dsdx(i,1,1,e) * us(i,1,1,e) &
933 + dtdx(i,1,1,e) * ut(i,1,1,e))
934
935 uy(i,1,1,e) = w3(i,1,1) * &
936 ( dsdy(i,1,1,e) * us(i,1,1,e) &
937 + drdy(i,1,1,e) * ur(i,1,1,e) &
938 + dtdy(i,1,1,e) * ut(i,1,1,e) )
939
940 uz(i,1,1,e) = w3(i,1,1) * &
941 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
942 + drdz(i,1,1,e) * ur(i,1,1,e) &
943 + dsdz(i,1,1,e) * us(i,1,1,e))
944 end do
945 end do
946 end subroutine sx_opgrad_lx10
947
948 subroutine sx_opgrad_lx9(ux, uy, uz, u, dx, dy, dz, &
949 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
950 integer, parameter :: lx = 9
951 integer, intent(in) :: n
952 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
953 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
954 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
955 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
956 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
957 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
958 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
959 real(kind=rp) :: ur(lx, lx, lx, n)
960 real(kind=rp) :: us(lx, lx, lx, n)
961 real(kind=rp) :: ut(lx, lx, lx, n)
962 real(kind=rp) :: wr, ws, wt
963 integer :: e, i, j, k, jj, kk
964
965 do i = 1, lx
966 do jj = 1, lx * lx * n
967 wr = 0d0
968 do kk = 1, lx
969 wr = wr + dx(i, kk) * u(kk, jj,1,1)
970 end do
971 ur(i, jj,1,1) = wr
972 end do
973 end do
974
975 do k = 1, lx
976 do i = 1, lx
977 do j = 1, lx
978 do e = 1, n
979 ws = 0d0
980 !NEC$ unroll_completely
981 do kk = 1, lx
982 ws = ws + dy(j, kk) * u(i, kk, k,e)
983 end do
984 us(i, j, k,e) = ws
985 end do
986 end do
987 end do
988 end do
989
990 do j = 1, lx
991 do i = 1, lx
992 do k = 1, lx
993 do e = 1, n
994 wt = 0d0
995 !NEC$ unroll_completely
996 do kk = 1, lx
997 wt = wt + dz(k, kk) * u(i, j, kk,e)
998 end do
999 ut(i, j, k,e) = wt
1000 end do
1001 end do
1002 end do
1003 end do
1004
1005 do i = 1, lx * lx * lx
1006 do e = 1, n
1007 ux(i,1,1,e) = w3(i,1,1) * &
1008 ( drdx(i,1,1,e) * ur(i,1,1,e) &
1009 + dsdx(i,1,1,e) * us(i,1,1,e) &
1010 + dtdx(i,1,1,e) * ut(i,1,1,e))
1011
1012 uy(i,1,1,e) = w3(i,1,1) * &
1013 ( dsdy(i,1,1,e) * us(i,1,1,e) &
1014 + drdy(i,1,1,e) * ur(i,1,1,e) &
1015 + dtdy(i,1,1,e) * ut(i,1,1,e) )
1016
1017 uz(i,1,1,e) = w3(i,1,1) * &
1018 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
1019 + drdz(i,1,1,e) * ur(i,1,1,e) &
1020 + dsdz(i,1,1,e) * us(i,1,1,e))
1021 end do
1022 end do
1023 end subroutine sx_opgrad_lx9
1024
1025 subroutine sx_opgrad_lx8(ux, uy, uz, u, dx, dy, dz, &
1026 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1027 integer, parameter :: lx = 8
1028 integer, intent(in) :: n
1029 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1030 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1031 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1032 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1033 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1034 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1035 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1036 real(kind=rp) :: ur(lx, lx, lx, n)
1037 real(kind=rp) :: us(lx, lx, lx, n)
1038 real(kind=rp) :: ut(lx, lx, lx, n)
1039 real(kind=rp) :: wr, ws, wt
1040 integer :: e, i, j, k, jj, kk
1041
1042 do i = 1, lx
1043 do jj = 1, lx * lx * n
1044 wr = 0d0
1045 do kk = 1, lx
1046 wr = wr + dx(i, kk) * u(kk, jj,1,1)
1047 end do
1048 ur(i, jj,1,1) = wr
1049 end do
1050 end do
1051
1052 do k = 1, lx
1053 do i = 1, lx
1054 do j = 1, lx
1055 do e = 1, n
1056 ws = 0d0
1057 !NEC$ unroll_completely
1058 do kk = 1, lx
1059 ws = ws + dy(j, kk) * u(i, kk, k,e)
1060 end do
1061 us(i, j, k,e) = ws
1062 end do
1063 end do
1064 end do
1065 end do
1066
1067 do j = 1, lx
1068 do i = 1, lx
1069 do k = 1, lx
1070 do e = 1, n
1071 wt = 0d0
1072 !NEC$ unroll_completely
1073 do kk = 1, lx
1074 wt = wt + dz(k, kk) * u(i, j, kk,e)
1075 end do
1076 ut(i, j, k,e) = wt
1077 end do
1078 end do
1079 end do
1080 end do
1081
1082 do i = 1, lx * lx * lx
1083 do e = 1, n
1084 ux(i,1,1,e) = w3(i,1,1) * &
1085 ( drdx(i,1,1,e) * ur(i,1,1,e) &
1086 + dsdx(i,1,1,e) * us(i,1,1,e) &
1087 + dtdx(i,1,1,e) * ut(i,1,1,e))
1088
1089 uy(i,1,1,e) = w3(i,1,1) * &
1090 ( dsdy(i,1,1,e) * us(i,1,1,e) &
1091 + drdy(i,1,1,e) * ur(i,1,1,e) &
1092 + dtdy(i,1,1,e) * ut(i,1,1,e) )
1093
1094 uz(i,1,1,e) = w3(i,1,1) * &
1095 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
1096 + drdz(i,1,1,e) * ur(i,1,1,e) &
1097 + dsdz(i,1,1,e) * us(i,1,1,e))
1098 end do
1099 end do
1100 end subroutine sx_opgrad_lx8
1101
1102 subroutine sx_opgrad_lx7(ux, uy, uz, u, dx, dy, dz, &
1103 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1104 integer, parameter :: lx = 7
1105 integer, intent(in) :: n
1106 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1107 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1108 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1109 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1110 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1111 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1112 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1113 real(kind=rp) :: ur(lx, lx, lx, n)
1114 real(kind=rp) :: us(lx, lx, lx, n)
1115 real(kind=rp) :: ut(lx, lx, lx, n)
1116 real(kind=rp) :: wr, ws, wt
1117 integer :: e, i, j, k, jj, kk
1118
1119 do i = 1, lx
1120 do jj = 1, lx * lx * n
1121 wr = 0d0
1122 do kk = 1, lx
1123 wr = wr + dx(i, kk) * u(kk, jj,1,1)
1124 end do
1125 ur(i, jj,1,1) = wr
1126 end do
1127 end do
1128
1129 do k = 1, lx
1130 do i = 1, lx
1131 do j = 1, lx
1132 do e = 1, n
1133 ws = 0d0
1134 !NEC$ unroll_completely
1135 do kk = 1, lx
1136 ws = ws + dy(j, kk) * u(i, kk, k,e)
1137 end do
1138 us(i, j, k,e) = ws
1139 end do
1140 end do
1141 end do
1142 end do
1143
1144 do j = 1, lx
1145 do i = 1, lx
1146 do k = 1, lx
1147 do e = 1, n
1148 wt = 0d0
1149 !NEC$ unroll_completely
1150 do kk = 1, lx
1151 wt = wt + dz(k, kk) * u(i, j, kk,e)
1152 end do
1153 ut(i, j, k,e) = wt
1154 end do
1155 end do
1156 end do
1157 end do
1158
1159 do i = 1, lx * lx * lx
1160 do e = 1, n
1161 ux(i,1,1,e) = w3(i,1,1) * &
1162 ( drdx(i,1,1,e) * ur(i,1,1,e) &
1163 + dsdx(i,1,1,e) * us(i,1,1,e) &
1164 + dtdx(i,1,1,e) * ut(i,1,1,e))
1165
1166 uy(i,1,1,e) = w3(i,1,1) * &
1167 ( dsdy(i,1,1,e) * us(i,1,1,e) &
1168 + drdy(i,1,1,e) * ur(i,1,1,e) &
1169 + dtdy(i,1,1,e) * ut(i,1,1,e) )
1170
1171 uz(i,1,1,e) = w3(i,1,1) * &
1172 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
1173 + drdz(i,1,1,e) * ur(i,1,1,e) &
1174 + dsdz(i,1,1,e) * us(i,1,1,e))
1175 end do
1176 end do
1177 end subroutine sx_opgrad_lx7
1178
1179 subroutine sx_opgrad_lx6(ux, uy, uz, u, dx, dy, dz, &
1180 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1181 integer, parameter :: lx = 6
1182 integer, intent(in) :: n
1183 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1184 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1185 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1186 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1187 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1188 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1189 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1190 real(kind=rp) :: ur(lx, lx, lx, n)
1191 real(kind=rp) :: us(lx, lx, lx, n)
1192 real(kind=rp) :: ut(lx, lx, lx, n)
1193 real(kind=rp) :: wr, ws, wt
1194 integer :: e, i, j, k, jj, kk
1195
1196 do i = 1, lx
1197 do jj = 1, lx * lx * n
1198 wr = 0d0
1199 do kk = 1, lx
1200 wr = wr + dx(i, kk) * u(kk, jj,1,1)
1201 end do
1202 ur(i, jj,1,1) = wr
1203 end do
1204 end do
1205
1206 do k = 1, lx
1207 do i = 1, lx
1208 do j = 1, lx
1209 do e = 1, n
1210 ws = 0d0
1211 !NEC$ unroll_completely
1212 do kk = 1, lx
1213 ws = ws + dy(j, kk) * u(i, kk, k,e)
1214 end do
1215 us(i, j, k,e) = ws
1216 end do
1217 end do
1218 end do
1219 end do
1220
1221 do j = 1, lx
1222 do i = 1, lx
1223 do k = 1, lx
1224 do e = 1, n
1225 wt = 0d0
1226 !NEC$ unroll_completely
1227 do kk = 1, lx
1228 wt = wt + dz(k, kk) * u(i, j, kk,e)
1229 end do
1230 ut(i, j, k,e) = wt
1231 end do
1232 end do
1233 end do
1234 end do
1235
1236 do i = 1, lx * lx * lx
1237 do e = 1, n
1238 ux(i,1,1,e) = w3(i,1,1) * &
1239 ( drdx(i,1,1,e) * ur(i,1,1,e) &
1240 + dsdx(i,1,1,e) * us(i,1,1,e) &
1241 + dtdx(i,1,1,e) * ut(i,1,1,e))
1242
1243 uy(i,1,1,e) = w3(i,1,1) * &
1244 ( dsdy(i,1,1,e) * us(i,1,1,e) &
1245 + drdy(i,1,1,e) * ur(i,1,1,e) &
1246 + dtdy(i,1,1,e) * ut(i,1,1,e) )
1247
1248 uz(i,1,1,e) = w3(i,1,1) * &
1249 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
1250 + drdz(i,1,1,e) * ur(i,1,1,e) &
1251 + dsdz(i,1,1,e) * us(i,1,1,e))
1252 end do
1253 end do
1254 end subroutine sx_opgrad_lx6
1255
1256 subroutine sx_opgrad_lx5(ux, uy, uz, u, dx, dy, dz, &
1257 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1258 integer, parameter :: lx = 5
1259 integer, intent(in) :: n
1260 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1261 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1262 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1263 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1264 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1265 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1266 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1267 real(kind=rp) :: ur(lx, lx, lx, n)
1268 real(kind=rp) :: us(lx, lx, lx, n)
1269 real(kind=rp) :: ut(lx, lx, lx, n)
1270 real(kind=rp) :: wr, ws, wt
1271 integer :: e, i, j, k, jj, kk
1272
1273 do i = 1, lx
1274 do jj = 1, lx * lx * n
1275 wr = 0d0
1276 do kk = 1, lx
1277 wr = wr + dx(i, kk) * u(kk, jj,1,1)
1278 end do
1279 ur(i, jj,1,1) = wr
1280 end do
1281 end do
1282
1283 do k = 1, lx
1284 do i = 1, lx
1285 do j = 1, lx
1286 do e = 1, n
1287 ws = 0d0
1288 !NEC$ unroll_completely
1289 do kk = 1, lx
1290 ws = ws + dy(j, kk) * u(i, kk, k,e)
1291 end do
1292 us(i, j, k,e) = ws
1293 end do
1294 end do
1295 end do
1296 end do
1297
1298 do j = 1, lx
1299 do i = 1, lx
1300 do k = 1, lx
1301 do e = 1, n
1302 wt = 0d0
1303 !NEC$ unroll_completely
1304 do kk = 1, lx
1305 wt = wt + dz(k, kk) * u(i, j, kk,e)
1306 end do
1307 ut(i, j, k,e) = wt
1308 end do
1309 end do
1310 end do
1311 end do
1312
1313 do i = 1, lx * lx * lx
1314 do e = 1, n
1315 ux(i,1,1,e) = w3(i,1,1) * &
1316 ( drdx(i,1,1,e) * ur(i,1,1,e) &
1317 + dsdx(i,1,1,e) * us(i,1,1,e) &
1318 + dtdx(i,1,1,e) * ut(i,1,1,e))
1319
1320 uy(i,1,1,e) = w3(i,1,1) * &
1321 ( dsdy(i,1,1,e) * us(i,1,1,e) &
1322 + drdy(i,1,1,e) * ur(i,1,1,e) &
1323 + dtdy(i,1,1,e) * ut(i,1,1,e) )
1324
1325 uz(i,1,1,e) = w3(i,1,1) * &
1326 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
1327 + drdz(i,1,1,e) * ur(i,1,1,e) &
1328 + dsdz(i,1,1,e) * us(i,1,1,e))
1329 end do
1330 end do
1331 end subroutine sx_opgrad_lx5
1332
1333 subroutine sx_opgrad_lx4(ux, uy, uz, u, dx, dy, dz, &
1334 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1335 integer, parameter :: lx = 4
1336 integer, intent(in) :: n
1337 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1338 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1339 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1340 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1341 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1342 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1343 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1344 real(kind=rp) :: ur(lx, lx, lx, n)
1345 real(kind=rp) :: us(lx, lx, lx, n)
1346 real(kind=rp) :: ut(lx, lx, lx, n)
1347 real(kind=rp) :: wr, ws, wt
1348 integer :: e, i, j, k, jj, kk
1349
1350 do i = 1, lx
1351 do jj = 1, lx * lx * n
1352 wr = 0d0
1353 do kk = 1, lx
1354 wr = wr + dx(i, kk) * u(kk, jj,1,1)
1355 end do
1356 ur(i, jj,1,1) = wr
1357 end do
1358 end do
1359
1360 do k = 1, lx
1361 do i = 1, lx
1362 do j = 1, lx
1363 do e = 1, n
1364 ws = 0d0
1365 !NEC$ unroll_completely
1366 do kk = 1, lx
1367 ws = ws + dy(j, kk) * u(i, kk, k,e)
1368 end do
1369 us(i, j, k,e) = ws
1370 end do
1371 end do
1372 end do
1373 end do
1374
1375 do j = 1, lx
1376 do i = 1, lx
1377 do k = 1, lx
1378 do e = 1, n
1379 wt = 0d0
1380 !NEC$ unroll_completely
1381 do kk = 1, lx
1382 wt = wt + dz(k, kk) * u(i, j, kk,e)
1383 end do
1384 ut(i, j, k,e) = wt
1385 end do
1386 end do
1387 end do
1388 end do
1389
1390 do i = 1, lx * lx * lx
1391 do e = 1, n
1392 ux(i,1,1,e) = w3(i,1,1) * &
1393 ( drdx(i,1,1,e) * ur(i,1,1,e) &
1394 + dsdx(i,1,1,e) * us(i,1,1,e) &
1395 + dtdx(i,1,1,e) * ut(i,1,1,e))
1396
1397 uy(i,1,1,e) = w3(i,1,1) * &
1398 ( dsdy(i,1,1,e) * us(i,1,1,e) &
1399 + drdy(i,1,1,e) * ur(i,1,1,e) &
1400 + dtdy(i,1,1,e) * ut(i,1,1,e) )
1401
1402 uz(i,1,1,e) = w3(i,1,1) * &
1403 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
1404 + drdz(i,1,1,e) * ur(i,1,1,e) &
1405 + dsdz(i,1,1,e) * us(i,1,1,e))
1406 end do
1407 end do
1408 end subroutine sx_opgrad_lx4
1409
1410 subroutine sx_opgrad_lx3(ux, uy, uz, u, dx, dy, dz, &
1411 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1412 integer, parameter :: lx = 3
1413 integer, intent(in) :: n
1414 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1415 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1416 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1417 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1418 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1419 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1420 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1421 real(kind=rp) :: ur(lx, lx, lx, n)
1422 real(kind=rp) :: us(lx, lx, lx, n)
1423 real(kind=rp) :: ut(lx, lx, lx, n)
1424 real(kind=rp) :: wr, ws, wt
1425 integer :: e, i, j, k, jj, kk
1426
1427 do i = 1, lx
1428 do jj = 1, lx * lx * n
1429 wr = 0d0
1430 do kk = 1, lx
1431 wr = wr + dx(i, kk) * u(kk, jj,1,1)
1432 end do
1433 ur(i, jj,1,1) = wr
1434 end do
1435 end do
1436
1437 do k = 1, lx
1438 do i = 1, lx
1439 do j = 1, lx
1440 do e = 1, n
1441 ws = 0d0
1442 !NEC$ unroll_completely
1443 do kk = 1, lx
1444 ws = ws + dy(j, kk) * u(i, kk, k,e)
1445 end do
1446 us(i, j, k,e) = ws
1447 end do
1448 end do
1449 end do
1450 end do
1451
1452 do j = 1, lx
1453 do i = 1, lx
1454 do k = 1, lx
1455 do e = 1, n
1456 wt = 0d0
1457 !NEC$ unroll_completely
1458 do kk = 1, lx
1459 wt = wt + dz(k, kk) * u(i, j, kk,e)
1460 end do
1461 ut(i, j, k,e) = wt
1462 end do
1463 end do
1464 end do
1465 end do
1466
1467 do i = 1, lx * lx * lx
1468 do e = 1, n
1469 ux(i,1,1,e) = w3(i,1,1) * &
1470 ( drdx(i,1,1,e) * ur(i,1,1,e) &
1471 + dsdx(i,1,1,e) * us(i,1,1,e) &
1472 + dtdx(i,1,1,e) * ut(i,1,1,e))
1473
1474 uy(i,1,1,e) = w3(i,1,1) * &
1475 ( dsdy(i,1,1,e) * us(i,1,1,e) &
1476 + drdy(i,1,1,e) * ur(i,1,1,e) &
1477 + dtdy(i,1,1,e) * ut(i,1,1,e) )
1478
1479 uz(i,1,1,e) = w3(i,1,1) * &
1480 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
1481 + drdz(i,1,1,e) * ur(i,1,1,e) &
1482 + dsdz(i,1,1,e) * us(i,1,1,e))
1483 end do
1484 end do
1485 end subroutine sx_opgrad_lx3
1486
1487 subroutine sx_opgrad_lx2(ux, uy, uz, u, dx, dy, dz, &
1488 drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1489 integer, parameter :: lx = 2
1490 integer, intent(in) :: n
1491 real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: ux, uy, uz
1492 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: u
1493 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1494 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
1495 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
1496 real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
1497 real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
1498 real(kind=rp) :: ur(lx, lx, lx, n)
1499 real(kind=rp) :: us(lx, lx, lx, n)
1500 real(kind=rp) :: ut(lx, lx, lx, n)
1501 real(kind=rp) :: wr, ws, wt
1502 integer :: e, i, j, k, jj, kk
1503
1504 do i = 1, lx
1505 do jj = 1, lx * lx * n
1506 wr = 0d0
1507 do kk = 1, lx
1508 wr = wr + dx(i, kk) * u(kk, jj,1,1)
1509 end do
1510 ur(i, jj,1,1) = wr
1511 end do
1512 end do
1513
1514 do k = 1, lx
1515 do i = 1, lx
1516 do j = 1, lx
1517 do e = 1, n
1518 ws = 0d0
1519 !NEC$ unroll_completely
1520 do kk = 1, lx
1521 ws = ws + dy(j, kk) * u(i, kk, k,e)
1522 end do
1523 us(i, j, k,e) = ws
1524 end do
1525 end do
1526 end do
1527 end do
1528
1529 do j = 1, lx
1530 do i = 1, lx
1531 do k = 1, lx
1532 do e = 1, n
1533 wt = 0d0
1534 !NEC$ unroll_completely
1535 do kk = 1, lx
1536 wt = wt + dz(k, kk) * u(i, j, kk,e)
1537 end do
1538 ut(i, j, k,e) = wt
1539 end do
1540 end do
1541 end do
1542 end do
1543
1544 do i = 1, lx * lx * lx
1545 do e = 1, n
1546 ux(i,1,1,e) = w3(i,1,1) * &
1547 ( drdx(i,1,1,e) * ur(i,1,1,e) &
1548 + dsdx(i,1,1,e) * us(i,1,1,e) &
1549 + dtdx(i,1,1,e) * ut(i,1,1,e))
1550
1551 uy(i,1,1,e) = w3(i,1,1) * &
1552 ( dsdy(i,1,1,e) * us(i,1,1,e) &
1553 + drdy(i,1,1,e) * ur(i,1,1,e) &
1554 + dtdy(i,1,1,e) * ut(i,1,1,e) )
1555
1556 uz(i,1,1,e) = w3(i,1,1) * &
1557 ( dtdz(i,1,1,e) * ut(i,1,1,e) &
1558 + drdz(i,1,1,e) * ur(i,1,1,e) &
1559 + dsdz(i,1,1,e) * us(i,1,1,e))
1560 end do
1561 end do
1562 end subroutine sx_opgrad_lx2
1563
1564end submodule sx_opgrad
Operators SX-Aurora backend.
Definition opr_sx.f90:2