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