Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
cpu_dudxyz.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_dudxyz
35 implicit none
36
37contains
38
39 module subroutine opr_cpu_dudxyz(du, u, dr, ds, dt, coef)
40 type(coef_t), intent(in), target :: coef
41 real(kind=rp), dimension(coef%Xh%lx, coef%Xh%ly, & coef%Xh%lz, coef%msh%nelv), intent(inout) :: du
42 real(kind=rp), dimension(coef%Xh%lx, coef%Xh%ly, & coef%Xh%lz, coef%msh%nelv), intent(in) :: u, dr, ds, dt
43
44 associate(xh => coef%Xh, msh => coef%msh, dof => coef%dof)
45 !$omp parallel
46 select case (coef%Xh%lx)
47 case (14)
48 call cpu_dudxyz_lx14(du, u, dr, ds, dt, &
49 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
50 case (13)
51 call cpu_dudxyz_lx13(du, u, dr, ds, dt, &
52 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
53 case (12)
54 call cpu_dudxyz_lx12(du, u, dr, ds, dt, &
55 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
56 case (11)
57 call cpu_dudxyz_lx11(du, u, dr, ds, dt, &
58 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
59 case (10)
60 call cpu_dudxyz_lx10(du, u, dr, ds, dt, &
61 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
62 case (9)
63 call cpu_dudxyz_lx9(du, u, dr, ds, dt, &
64 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
65 case (8)
66 call cpu_dudxyz_lx8(du, u, dr, ds, dt, &
67 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
68 case (7)
69 call cpu_dudxyz_lx7(du, u, dr, ds, dt, &
70 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
71 case (6)
72 call cpu_dudxyz_lx6(du, u, dr, ds, dt, &
73 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
74 case (5)
75 call cpu_dudxyz_lx5(du, u, dr, ds, dt, &
76 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
77 case (4)
78 call cpu_dudxyz_lx4(du, u, dr, ds, dt, &
79 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
80 case (3)
81 call cpu_dudxyz_lx3(du, u, dr, ds, dt, &
82 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
83 case (2)
84 call cpu_dudxyz_lx2(du, u, dr, ds, dt, &
85 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
86 case default
87 call cpu_dudxyz_lx(du, u, dr, ds, dt, &
88 xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv, xh%lx)
89 end select
90 !$omp end parallel
91 end associate
92
93 end subroutine opr_cpu_dudxyz
94
95 subroutine cpu_dudxyz_lx(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, lx)
96 integer, intent(in) :: nel, lx
97 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
98 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
99 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
100 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
101 real(kind=rp), dimension(lx, lx, lx) :: drst
102 real(kind=rp) :: tmp
103 integer :: e, i, j, k, l
104
105 !$omp do
106 do e = 1, nel
107 do j = 1, lx * lx
108 do i = 1, lx
109 tmp = 0.0_rp
110 do k = 1, lx
111 tmp = tmp + dx(i,k) * u(k,j,1,e)
112 end do
113 du(i,j,1,e) = tmp
114 end do
115 end do
116
117 do i = 1, lx * lx * lx
118 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
119 end do
120
121 do k = 1, lx
122 do j = 1, lx
123 do i = 1, lx
124 tmp = 0.0_rp
125 do l = 1, lx
126 tmp = tmp + dy(j,l) * u(i,l,k,e)
127 end do
128 drst(i,j,k) = tmp
129 end do
130 end do
131 end do
132
133 do i = 1, lx * lx * lx
134 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
135 end do
136
137 do k = 1, lx
138 do i = 1, lx*lx
139 tmp = 0.0_rp
140 do l = 1, lx
141 tmp = tmp + dz(k,l) * u(i,1,l,e)
142 end do
143 drst(i,1,k) = tmp
144 end do
145 end do
146
147 do i = 1, lx * lx * lx
148 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
149 end do
150
151 do i = 1, lx * lx * lx
152 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
153 end do
154
155 end do
156 !$omp end do
157 end subroutine cpu_dudxyz_lx
158
159 subroutine cpu_dudxyz_lx14(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
160 integer, parameter :: lx = 14
161 integer, intent(in) :: nel
162 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
163 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
164 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
165 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
166 real(kind=rp), dimension(lx, lx, lx) :: drst
167 integer :: e, i, j, k
168
169 !$omp do
170 do e = 1, nel
171 do j = 1, lx * lx
172 do i = 1, lx
173 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
174 + dx(i,2) * u(2,j,1,e) &
175 + dx(i,3) * u(3,j,1,e) &
176 + dx(i,4) * u(4,j,1,e) &
177 + dx(i,5) * u(5,j,1,e) &
178 + dx(i,6) * u(6,j,1,e) &
179 + dx(i,7) * u(7,j,1,e) &
180 + dx(i,8) * u(8,j,1,e) &
181 + dx(i,9) * u(9,j,1,e) &
182 + dx(i,10) * u(10,j,1,e) &
183 + dx(i,11) * u(11,j,1,e) &
184 + dx(i,12) * u(12,j,1,e) &
185 + dx(i,13) * u(13,j,1,e) &
186 + dx(i,14) * u(14,j,1,e)
187 end do
188 end do
189
190 do i = 1, lx * lx * lx
191 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
192 end do
193
194 do k = 1, lx
195 do j = 1, lx
196 do i = 1, lx
197 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
198 + dy(j,2) * u(i,2,k,e) &
199 + dy(j,3) * u(i,3,k,e) &
200 + dy(j,4) * u(i,4,k,e) &
201 + dy(j,5) * u(i,5,k,e) &
202 + dy(j,6) * u(i,6,k,e) &
203 + dy(j,7) * u(i,7,k,e) &
204 + dy(j,8) * u(i,8,k,e) &
205 + dy(j,9) * u(i,9,k,e) &
206 + dy(j,10) * u(i,10,k,e) &
207 + dy(j,11) * u(i,11,k,e) &
208 + dy(j,12) * u(i,12,k,e) &
209 + dy(j,13) * u(i,13,k,e) &
210 + dy(j,14) * u(i,14,k,e)
211 end do
212 end do
213 end do
214
215 do i = 1, lx * lx * lx
216 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
217 end do
218
219 do k = 1, lx
220 do i = 1, lx*lx
221 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
222 + dz(k,2) * u(i,1,2,e) &
223 + dz(k,3) * u(i,1,3,e) &
224 + dz(k,4) * u(i,1,4,e) &
225 + dz(k,5) * u(i,1,5,e) &
226 + dz(k,6) * u(i,1,6,e) &
227 + dz(k,7) * u(i,1,7,e) &
228 + dz(k,8) * u(i,1,8,e) &
229 + dz(k,9) * u(i,1,9,e) &
230 + dz(k,10) * u(i,1,10,e) &
231 + dz(k,11) * u(i,1,11,e) &
232 + dz(k,12) * u(i,1,12,e) &
233 + dz(k,13) * u(i,1,13,e) &
234 + dz(k,14) * u(i,1,14,e)
235 end do
236 end do
237
238 do i = 1, lx * lx * lx
239 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
240 end do
241
242 do i = 1, lx * lx * lx
243 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
244 end do
245
246 end do
247 !$omp end do
248 end subroutine cpu_dudxyz_lx14
249
250 subroutine cpu_dudxyz_lx13(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
251 integer, parameter :: lx = 13
252 integer, intent(in) :: nel
253 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
254 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
255 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
256 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
257 real(kind=rp), dimension(lx, lx, lx) :: drst
258 integer :: e, i, j, k
259
260 !$omp do
261 do e = 1, nel
262 do j = 1, lx * lx
263 do i = 1, lx
264 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
265 + dx(i,2) * u(2,j,1,e) &
266 + dx(i,3) * u(3,j,1,e) &
267 + dx(i,4) * u(4,j,1,e) &
268 + dx(i,5) * u(5,j,1,e) &
269 + dx(i,6) * u(6,j,1,e) &
270 + dx(i,7) * u(7,j,1,e) &
271 + dx(i,8) * u(8,j,1,e) &
272 + dx(i,9) * u(9,j,1,e) &
273 + dx(i,10) * u(10,j,1,e) &
274 + dx(i,11) * u(11,j,1,e) &
275 + dx(i,12) * u(12,j,1,e) &
276 + dx(i,13) * u(13,j,1,e)
277 end do
278 end do
279
280 do i = 1, lx * lx * lx
281 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
282 end do
283
284 do k = 1, lx
285 do j = 1, lx
286 do i = 1, lx
287 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
288 + dy(j,2) * u(i,2,k,e) &
289 + dy(j,3) * u(i,3,k,e) &
290 + dy(j,4) * u(i,4,k,e) &
291 + dy(j,5) * u(i,5,k,e) &
292 + dy(j,6) * u(i,6,k,e) &
293 + dy(j,7) * u(i,7,k,e) &
294 + dy(j,8) * u(i,8,k,e) &
295 + dy(j,9) * u(i,9,k,e) &
296 + dy(j,10) * u(i,10,k,e) &
297 + dy(j,11) * u(i,11,k,e) &
298 + dy(j,12) * u(i,12,k,e) &
299 + dy(j,13) * u(i,13,k,e)
300 end do
301 end do
302 end do
303
304 do i = 1, lx * lx * lx
305 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
306 end do
307
308 do k = 1, lx
309 do i = 1, lx*lx
310 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
311 + dz(k,2) * u(i,1,2,e) &
312 + dz(k,3) * u(i,1,3,e) &
313 + dz(k,4) * u(i,1,4,e) &
314 + dz(k,5) * u(i,1,5,e) &
315 + dz(k,6) * u(i,1,6,e) &
316 + dz(k,7) * u(i,1,7,e) &
317 + dz(k,8) * u(i,1,8,e) &
318 + dz(k,9) * u(i,1,9,e) &
319 + dz(k,10) * u(i,1,10,e) &
320 + dz(k,11) * u(i,1,11,e) &
321 + dz(k,12) * u(i,1,12,e) &
322 + dz(k,13) * u(i,1,13,e)
323 end do
324 end do
325
326 do i = 1, lx * lx * lx
327 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
328 end do
329
330 do i = 1, lx * lx * lx
331 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
332 end do
333
334 end do
335 !$omp end do
336 end subroutine cpu_dudxyz_lx13
337
338 subroutine cpu_dudxyz_lx12(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
339 integer, parameter :: lx = 12
340 integer, intent(in) :: nel
341 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
342 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
343 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
344 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
345 real(kind=rp), dimension(lx, lx, lx) :: drst
346 integer :: e, i, j, k
347 !$omp do
348 do e = 1, nel
349 do j = 1, lx * lx
350 do i = 1, lx
351 du(i,j,1,e) = 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 end do
364 end do
365
366 do i = 1, lx * lx * lx
367 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
368 end do
369
370 do k = 1, lx
371 do j = 1, lx
372 do i = 1, lx
373 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
374 + dy(j,2) * u(i,2,k,e) &
375 + dy(j,3) * u(i,3,k,e) &
376 + dy(j,4) * u(i,4,k,e) &
377 + dy(j,5) * u(i,5,k,e) &
378 + dy(j,6) * u(i,6,k,e) &
379 + dy(j,7) * u(i,7,k,e) &
380 + dy(j,8) * u(i,8,k,e) &
381 + dy(j,9) * u(i,9,k,e) &
382 + dy(j,10) * u(i,10,k,e) &
383 + dy(j,11) * u(i,11,k,e) &
384 + dy(j,12) * u(i,12,k,e)
385 end do
386 end do
387 end do
388
389 do i = 1, lx * lx * lx
390 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
391 end do
392
393 do k = 1, lx
394 do i = 1, lx*lx
395 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
396 + dz(k,2) * u(i,1,2,e) &
397 + dz(k,3) * u(i,1,3,e) &
398 + dz(k,4) * u(i,1,4,e) &
399 + dz(k,5) * u(i,1,5,e) &
400 + dz(k,6) * u(i,1,6,e) &
401 + dz(k,7) * u(i,1,7,e) &
402 + dz(k,8) * u(i,1,8,e) &
403 + dz(k,9) * u(i,1,9,e) &
404 + dz(k,10) * u(i,1,10,e) &
405 + dz(k,11) * u(i,1,11,e) &
406 + dz(k,12) * u(i,1,12,e)
407 end do
408 end do
409
410 do i = 1, lx * lx * lx
411 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
412 end do
413
414 do i = 1, lx * lx * lx
415 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
416 end do
417
418 end do
419 !$omp end do
420 end subroutine cpu_dudxyz_lx12
421
422 subroutine cpu_dudxyz_lx11(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
423 integer, parameter :: lx = 11
424 integer, intent(in) :: nel
425 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
426 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
427 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
428 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
429 real(kind=rp), dimension(lx, lx, lx) :: drst
430 integer :: e, i, j, k
431 !$omp do
432 do e = 1, nel
433 do j = 1, lx * lx
434 do i = 1, lx
435 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
436 + dx(i,2) * u(2,j,1,e) &
437 + dx(i,3) * u(3,j,1,e) &
438 + dx(i,4) * u(4,j,1,e) &
439 + dx(i,5) * u(5,j,1,e) &
440 + dx(i,6) * u(6,j,1,e) &
441 + dx(i,7) * u(7,j,1,e) &
442 + dx(i,8) * u(8,j,1,e) &
443 + dx(i,9) * u(9,j,1,e) &
444 + dx(i,10) * u(10,j,1,e) &
445 + dx(i,11) * u(11,j,1,e)
446 end do
447 end do
448
449 do i = 1, lx * lx * lx
450 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
451 end do
452
453 do k = 1, lx
454 do j = 1, lx
455 do i = 1, lx
456 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
457 + dy(j,2) * u(i,2,k,e) &
458 + dy(j,3) * u(i,3,k,e) &
459 + dy(j,4) * u(i,4,k,e) &
460 + dy(j,5) * u(i,5,k,e) &
461 + dy(j,6) * u(i,6,k,e) &
462 + dy(j,7) * u(i,7,k,e) &
463 + dy(j,8) * u(i,8,k,e) &
464 + dy(j,9) * u(i,9,k,e) &
465 + dy(j,10) * u(i,10,k,e) &
466 + dy(j,11) * u(i,11,k,e)
467 end do
468 end do
469 end do
470
471 do i = 1, lx * lx * lx
472 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
473 end do
474
475 do k = 1, lx
476 do i = 1, lx*lx
477 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
478 + dz(k,2) * u(i,1,2,e) &
479 + dz(k,3) * u(i,1,3,e) &
480 + dz(k,4) * u(i,1,4,e) &
481 + dz(k,5) * u(i,1,5,e) &
482 + dz(k,6) * u(i,1,6,e) &
483 + dz(k,7) * u(i,1,7,e) &
484 + dz(k,8) * u(i,1,8,e) &
485 + dz(k,9) * u(i,1,9,e) &
486 + dz(k,10) * u(i,1,10,e) &
487 + dz(k,11) * u(i,1,11,e)
488 end do
489 end do
490
491 do i = 1, lx * lx * lx
492 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
493 end do
494
495 do i = 1, lx * lx * lx
496 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
497 end do
498
499 end do
500 !$omp end do
501 end subroutine cpu_dudxyz_lx11
502
503 subroutine cpu_dudxyz_lx10(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
504 integer, parameter :: lx = 10
505 integer, intent(in) :: nel
506 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
507 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
508 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
509 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
510 real(kind=rp), dimension(lx, lx, lx) :: drst
511 integer :: e, i, j, k
512 !$omp do
513 do e = 1, nel
514 do j = 1, lx * lx
515 do i = 1, lx
516 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
517 + dx(i,2) * u(2,j,1,e) &
518 + dx(i,3) * u(3,j,1,e) &
519 + dx(i,4) * u(4,j,1,e) &
520 + dx(i,5) * u(5,j,1,e) &
521 + dx(i,6) * u(6,j,1,e) &
522 + dx(i,7) * u(7,j,1,e) &
523 + dx(i,8) * u(8,j,1,e) &
524 + dx(i,9) * u(9,j,1,e) &
525 + dx(i,10) * u(10,j,1,e)
526 end do
527 end do
528
529 do i = 1, lx * lx * lx
530 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
531 end do
532
533 do k = 1, lx
534 do j = 1, lx
535 do i = 1, lx
536 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
537 + dy(j,2) * u(i,2,k,e) &
538 + dy(j,3) * u(i,3,k,e) &
539 + dy(j,4) * u(i,4,k,e) &
540 + dy(j,5) * u(i,5,k,e) &
541 + dy(j,6) * u(i,6,k,e) &
542 + dy(j,7) * u(i,7,k,e) &
543 + dy(j,8) * u(i,8,k,e) &
544 + dy(j,9) * u(i,9,k,e) &
545 + dy(j,10) * u(i,10,k,e)
546 end do
547 end do
548 end do
549
550 do i = 1, lx * lx * lx
551 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
552 end do
553
554 do k = 1, lx
555 do i = 1, lx*lx
556 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
557 + dz(k,2) * u(i,1,2,e) &
558 + dz(k,3) * u(i,1,3,e) &
559 + dz(k,4) * u(i,1,4,e) &
560 + dz(k,5) * u(i,1,5,e) &
561 + dz(k,6) * u(i,1,6,e) &
562 + dz(k,7) * u(i,1,7,e) &
563 + dz(k,8) * u(i,1,8,e) &
564 + dz(k,9) * u(i,1,9,e) &
565 + dz(k,10) * u(i,1,10,e)
566 end do
567 end do
568
569 do i = 1, lx * lx * lx
570 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
571 end do
572
573 do i = 1, lx * lx * lx
574 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
575 end do
576
577 end do
578 !$omp end do
579 end subroutine cpu_dudxyz_lx10
580
581 subroutine cpu_dudxyz_lx9(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
582 integer, parameter :: lx = 9
583 integer, intent(in) :: nel
584 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
585 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
586 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
587 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
588 real(kind=rp), dimension(lx, lx, lx) :: drst
589 integer :: e, i, j, k
590 !$omp do
591 do e = 1, nel
592 do j = 1, lx * lx
593 do i = 1, lx
594 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
595 + dx(i,2) * u(2,j,1,e) &
596 + dx(i,3) * u(3,j,1,e) &
597 + dx(i,4) * u(4,j,1,e) &
598 + dx(i,5) * u(5,j,1,e) &
599 + dx(i,6) * u(6,j,1,e) &
600 + dx(i,7) * u(7,j,1,e) &
601 + dx(i,8) * u(8,j,1,e) &
602 + dx(i,9) * u(9,j,1,e)
603 end do
604 end do
605
606 do i = 1, lx * lx * lx
607 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
608 end do
609
610 do k = 1, lx
611 do j = 1, lx
612 do i = 1, lx
613 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
614 + dy(j,2) * u(i,2,k,e) &
615 + dy(j,3) * u(i,3,k,e) &
616 + dy(j,4) * u(i,4,k,e) &
617 + dy(j,5) * u(i,5,k,e) &
618 + dy(j,6) * u(i,6,k,e) &
619 + dy(j,7) * u(i,7,k,e) &
620 + dy(j,8) * u(i,8,k,e) &
621 + dy(j,9) * u(i,9,k,e)
622 end do
623 end do
624 end do
625
626 do i = 1, lx * lx * lx
627 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
628 end do
629
630 do k = 1, lx
631 do i = 1, lx*lx
632 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
633 + dz(k,2) * u(i,1,2,e) &
634 + dz(k,3) * u(i,1,3,e) &
635 + dz(k,4) * u(i,1,4,e) &
636 + dz(k,5) * u(i,1,5,e) &
637 + dz(k,6) * u(i,1,6,e) &
638 + dz(k,7) * u(i,1,7,e) &
639 + dz(k,8) * u(i,1,8,e) &
640 + dz(k,9) * u(i,1,9,e)
641 end do
642 end do
643
644 do i = 1, lx * lx * lx
645 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
646 end do
647
648 do i = 1, lx * lx * lx
649 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
650 end do
651
652 end do
653 !$omp end do
654 end subroutine cpu_dudxyz_lx9
655
656 subroutine cpu_dudxyz_lx8(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
657 integer, parameter :: lx = 8
658 integer, intent(in) :: nel
659 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
660 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
661 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
662 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
663 real(kind=rp), dimension(lx, lx, lx) :: drst
664 integer :: e, i, j, k
665 !$omp do
666 do e = 1, nel
667 do j = 1, lx * lx
668 do i = 1, lx
669 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
670 + dx(i,2) * u(2,j,1,e) &
671 + dx(i,3) * u(3,j,1,e) &
672 + dx(i,4) * u(4,j,1,e) &
673 + dx(i,5) * u(5,j,1,e) &
674 + dx(i,6) * u(6,j,1,e) &
675 + dx(i,7) * u(7,j,1,e) &
676 + dx(i,8) * u(8,j,1,e)
677 end do
678 end do
679
680 do i = 1, lx * lx * lx
681 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
682 end do
683
684 do k = 1, lx
685 do j = 1, lx
686 do i = 1, lx
687 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
688 + dy(j,2) * u(i,2,k,e) &
689 + dy(j,3) * u(i,3,k,e) &
690 + dy(j,4) * u(i,4,k,e) &
691 + dy(j,5) * u(i,5,k,e) &
692 + dy(j,6) * u(i,6,k,e) &
693 + dy(j,7) * u(i,7,k,e) &
694 + dy(j,8) * u(i,8,k,e)
695 end do
696 end do
697 end do
698
699 do i = 1, lx * lx * lx
700 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
701 end do
702
703 do k = 1, lx
704 do i = 1, lx*lx
705 drst(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 end do
714 end do
715
716 do i = 1, lx * lx * lx
717 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
718 end do
719
720 do i = 1, lx * lx * lx
721 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
722 end do
723
724 end do
725 !$omp end do
726 end subroutine cpu_dudxyz_lx8
727
728 subroutine cpu_dudxyz_lx7(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
729 integer, parameter :: lx = 7
730 integer, intent(in) :: nel
731 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
732 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
733 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
734 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
735 real(kind=rp), dimension(lx, lx, lx) :: drst
736 integer :: e, i, j, k
737 !$omp do
738 do e = 1, nel
739 do j = 1, lx * lx
740 do i = 1, lx
741 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
742 + dx(i,2) * u(2,j,1,e) &
743 + dx(i,3) * u(3,j,1,e) &
744 + dx(i,4) * u(4,j,1,e) &
745 + dx(i,5) * u(5,j,1,e) &
746 + dx(i,6) * u(6,j,1,e) &
747 + dx(i,7) * u(7,j,1,e)
748 end do
749 end do
750
751 do i = 1, lx * lx * lx
752 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
753 end do
754
755 do k = 1, lx
756 do j = 1, lx
757 do i = 1, lx
758 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
759 + dy(j,2) * u(i,2,k,e) &
760 + dy(j,3) * u(i,3,k,e) &
761 + dy(j,4) * u(i,4,k,e) &
762 + dy(j,5) * u(i,5,k,e) &
763 + dy(j,6) * u(i,6,k,e) &
764 + dy(j,7) * u(i,7,k,e)
765 end do
766 end do
767 end do
768
769 do i = 1, lx * lx * lx
770 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
771 end do
772
773 do k = 1, lx
774 do i = 1, lx*lx
775 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
776 + dz(k,2) * u(i,1,2,e) &
777 + dz(k,3) * u(i,1,3,e) &
778 + dz(k,4) * u(i,1,4,e) &
779 + dz(k,5) * u(i,1,5,e) &
780 + dz(k,6) * u(i,1,6,e) &
781 + dz(k,7) * u(i,1,7,e)
782 end do
783 end do
784
785 do i = 1, lx * lx * lx
786 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
787 end do
788
789 do i = 1, lx * lx * lx
790 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
791 end do
792
793 end do
794 !$omp end do
795 end subroutine cpu_dudxyz_lx7
796
797 subroutine cpu_dudxyz_lx6(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
798 integer, parameter :: lx = 6
799 integer, intent(in) :: nel
800 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
801 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
802 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
803 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
804 real(kind=rp), dimension(lx, lx, lx) :: drst
805 integer :: e, i, j, k
806 !$omp do
807 do e = 1, nel
808 do j = 1, lx * lx
809 do i = 1, lx
810 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
811 + dx(i,2) * u(2,j,1,e) &
812 + dx(i,3) * u(3,j,1,e) &
813 + dx(i,4) * u(4,j,1,e) &
814 + dx(i,5) * u(5,j,1,e) &
815 + dx(i,6) * u(6,j,1,e)
816 end do
817 end do
818
819 do i = 1, lx * lx * lx
820 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
821 end do
822
823 do k = 1, lx
824 do j = 1, lx
825 do i = 1, lx
826 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
827 + dy(j,2) * u(i,2,k,e) &
828 + dy(j,3) * u(i,3,k,e) &
829 + dy(j,4) * u(i,4,k,e) &
830 + dy(j,5) * u(i,5,k,e) &
831 + dy(j,6) * u(i,6,k,e)
832 end do
833 end do
834 end do
835
836 do i = 1, lx * lx * lx
837 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
838 end do
839
840 do k = 1, lx
841 do i = 1, lx*lx
842 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
843 + dz(k,2) * u(i,1,2,e) &
844 + dz(k,3) * u(i,1,3,e) &
845 + dz(k,4) * u(i,1,4,e) &
846 + dz(k,5) * u(i,1,5,e) &
847 + dz(k,6) * u(i,1,6,e)
848 end do
849 end do
850
851 do i = 1, lx * lx * lx
852 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
853 end do
854
855 do i = 1, lx * lx * lx
856 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
857 end do
858
859 end do
860 !$omp end do
861 end subroutine cpu_dudxyz_lx6
862
863 subroutine cpu_dudxyz_lx5(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
864 integer, parameter :: lx = 5
865 integer, intent(in) :: nel
866 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
867 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
868 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
869 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
870 real(kind=rp), dimension(lx, lx, lx) :: drst
871 integer :: e, i, j, k
872 !$omp do
873 do e = 1, nel
874 do j = 1, lx * lx
875 do i = 1, lx
876 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
877 + dx(i,2) * u(2,j,1,e) &
878 + dx(i,3) * u(3,j,1,e) &
879 + dx(i,4) * u(4,j,1,e) &
880 + dx(i,5) * u(5,j,1,e)
881 end do
882 end do
883
884 do i = 1, lx * lx * lx
885 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
886 end do
887
888 do k = 1, lx
889 do j = 1, lx
890 do i = 1, lx
891 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
892 + dy(j,2) * u(i,2,k,e) &
893 + dy(j,3) * u(i,3,k,e) &
894 + dy(j,4) * u(i,4,k,e) &
895 + dy(j,5) * u(i,5,k,e)
896 end do
897 end do
898 end do
899
900 do i = 1, lx * lx * lx
901 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
902 end do
903
904 do k = 1, lx
905 do i = 1, lx*lx
906 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
907 + dz(k,2) * u(i,1,2,e) &
908 + dz(k,3) * u(i,1,3,e) &
909 + dz(k,4) * u(i,1,4,e) &
910 + dz(k,5) * u(i,1,5,e)
911 end do
912 end do
913
914 do i = 1, lx * lx * lx
915 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
916 end do
917
918 do i = 1, lx * lx * lx
919 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
920 end do
921
922 end do
923 !$omp end do
924 end subroutine cpu_dudxyz_lx5
925
926 subroutine cpu_dudxyz_lx4(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
927 integer, parameter :: lx = 4
928 integer, intent(in) :: nel
929 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
930 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
931 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
932 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
933 real(kind=rp), dimension(lx, lx, lx) :: drst
934 integer :: e, i, j, k
935 !$omp do
936 do e = 1, nel
937 do j = 1, lx * lx
938 do i = 1, lx
939 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
940 + dx(i,2) * u(2,j,1,e) &
941 + dx(i,3) * u(3,j,1,e) &
942 + dx(i,4) * u(4,j,1,e)
943 end do
944 end do
945
946 do i = 1, lx * lx * lx
947 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
948 end do
949
950 do k = 1, lx
951 do j = 1, lx
952 do i = 1, lx
953 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
954 + dy(j,2) * u(i,2,k,e) &
955 + dy(j,3) * u(i,3,k,e) &
956 + dy(j,4) * u(i,4,k,e)
957 end do
958 end do
959 end do
960
961 do i = 1, lx * lx * lx
962 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
963 end do
964
965 do k = 1, lx
966 do i = 1, lx*lx
967 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
968 + dz(k,2) * u(i,1,2,e) &
969 + dz(k,3) * u(i,1,3,e) &
970 + dz(k,4) * u(i,1,4,e)
971 end do
972 end do
973
974 do i = 1, lx * lx * lx
975 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
976 end do
977
978 do i = 1, lx * lx * lx
979 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
980 end do
981
982 end do
983 !$omp end do
984 end subroutine cpu_dudxyz_lx4
985
986 subroutine cpu_dudxyz_lx3(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
987 integer, parameter :: lx = 3
988 integer, intent(in) :: nel
989 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
990 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
991 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
992 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
993 real(kind=rp), dimension(lx, lx, lx) :: drst
994 integer :: e, i, j, k
995 !$omp do
996 do e = 1, nel
997 do j = 1, lx * lx
998 do i = 1, lx
999 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
1000 + dx(i,2) * u(2,j,1,e) &
1001 + dx(i,3) * u(3,j,1,e)
1002 end do
1003 end do
1004
1005 do i = 1, lx * lx * lx
1006 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
1007 end do
1008
1009 do k = 1, lx
1010 do j = 1, lx
1011 do i = 1, lx
1012 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
1013 + dy(j,2) * u(i,2,k,e) &
1014 + dy(j,3) * u(i,3,k,e)
1015 end do
1016 end do
1017 end do
1018
1019 do i = 1, lx * lx * lx
1020 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
1021 end do
1022
1023 do k = 1, lx
1024 do i = 1, lx*lx
1025 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
1026 + dz(k,2) * u(i,1,2,e) &
1027 + dz(k,3) * u(i,1,3,e)
1028 end do
1029 end do
1030
1031 do i = 1, lx * lx * lx
1032 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
1033 end do
1034
1035 do i = 1, lx * lx * lx
1036 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
1037 end do
1038
1039 end do
1040 !$omp end do
1041 end subroutine cpu_dudxyz_lx3
1042
1043 subroutine cpu_dudxyz_lx2(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
1044 integer, parameter :: lx = 2
1045 integer, intent(in) :: nel
1046 real(kind=rp), dimension(lx, lx, lx, nel), intent(inout) :: du
1047 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: u, dr, ds, dt
1048 real(kind=rp), dimension(lx, lx, lx, nel), intent(in) :: jacinv
1049 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1050 real(kind=rp), dimension(lx, lx, lx) :: drst
1051 integer :: e, i, j, k
1052 !$omp do
1053 do e = 1, nel
1054 do j = 1, lx * lx
1055 do i = 1, lx
1056 du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
1057 + dx(i,2) * u(2,j,1,e)
1058 end do
1059 end do
1060
1061 do i = 1, lx * lx * lx
1062 du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
1063 end do
1064
1065 do k = 1, lx
1066 do j = 1, lx
1067 do i = 1, lx
1068 drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
1069 + dy(j,2) * u(i,2,k,e)
1070 end do
1071 end do
1072 end do
1073
1074 do i = 1, lx * lx * lx
1075 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
1076 end do
1077
1078 do k = 1, lx
1079 do i = 1, lx*lx
1080 drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
1081 + dz(k,2) * u(i,1,2,e)
1082 end do
1083 end do
1084
1085 do i = 1, lx * lx * lx
1086 du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
1087 end do
1088
1089 do i = 1, lx * lx * lx
1090 du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
1091 end do
1092
1093 end do
1094 !$omp end do
1095 end subroutine cpu_dudxyz_lx2
1096
1097end submodule cpu_dudxyz
1098
Operators CPU backend.
Definition opr_cpu.f90:34