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