40 subroutine cpu_cdtp_lx(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, lx)
41 integer,
intent(in) :: nel, lx
42 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
43 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
44 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
45 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
47 integer :: e, i, j, k, l
52 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
56 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
64 tmp = tmp + dxt(i,k) * ta1(k,j,1)
71 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
80 tmp = tmp + dyt(j,l) * ta1(i,l,k)
82 dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
88 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
96 tmp = tmp + dzt(k,l) * ta1(i,1,l)
98 dtx(i,1,k,e) = dtx(i,1,k,e) + tmp
105 subroutine cpu_cdtp_lx14(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
106 integer,
parameter :: lx = 14
107 integer,
intent(in) :: nel
108 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
109 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
110 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
111 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
112 integer :: e, i, j, k
117 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
121 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
126 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
127 + dxt(i,2) * ta1(2,j,1) &
128 + dxt(i,3) * ta1(3,j,1) &
129 + dxt(i,4) * ta1(4,j,1) &
130 + dxt(i,5) * ta1(5,j,1) &
131 + dxt(i,6) * ta1(6,j,1) &
132 + dxt(i,7) * ta1(7,j,1) &
133 + dxt(i,8) * ta1(8,j,1) &
134 + dxt(i,9) * ta1(9,j,1) &
135 + dxt(i,10) * ta1(10,j,1) &
136 + dxt(i,11) * ta1(11,j,1) &
137 + dxt(i,12) * ta1(12,j,1) &
138 + dxt(i,13) * ta1(13,j,1) &
139 + dxt(i,14) * ta1(14,j,1)
144 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
150 dtx(i,j,k,e) = dtx(i,j,k,e) &
151 + dyt(j,1) * ta1(i,1,k) &
152 + dyt(j,2) * ta1(i,2,k) &
153 + dyt(j,3) * ta1(i,3,k) &
154 + dyt(j,4) * ta1(i,4,k) &
155 + dyt(j,5) * ta1(i,5,k) &
156 + dyt(j,6) * ta1(i,6,k) &
157 + dyt(j,7) * ta1(i,7,k) &
158 + dyt(j,8) * ta1(i,8,k) &
159 + dyt(j,9) * ta1(i,9,k) &
160 + dyt(j,10) * ta1(i,10,k) &
161 + dyt(j,11) * ta1(i,11,k) &
162 + dyt(j,12) * ta1(i,12,k) &
163 + dyt(j,13) * ta1(i,13,k) &
164 + dyt(j,14) * ta1(i,14,k)
170 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
175 dtx(i,1,k,e) = dtx(i,1,k,e) &
176 + dzt(k,1) * ta1(i,1,1) &
177 + dzt(k,2) * ta1(i,1,2) &
178 + dzt(k,3) * ta1(i,1,3) &
179 + dzt(k,4) * ta1(i,1,4) &
180 + dzt(k,5) * ta1(i,1,5) &
181 + dzt(k,6) * ta1(i,1,6) &
182 + dzt(k,7) * ta1(i,1,7) &
183 + dzt(k,8) * ta1(i,1,8) &
184 + dzt(k,9) * ta1(i,1,9) &
185 + dzt(k,10) * ta1(i,1,10) &
186 + dzt(k,11) * ta1(i,1,11) &
187 + dzt(k,12) * ta1(i,1,12) &
188 + dzt(k,13) * ta1(i,1,13) &
189 + dzt(k,14) * ta1(i,1,14)
196 subroutine cpu_cdtp_lx13(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
197 integer,
parameter :: lx = 13
198 integer,
intent(in) :: nel
199 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
200 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
201 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
202 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
203 integer :: e, i, j, k
208 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
212 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
217 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
218 + dxt(i,2) * ta1(2,j,1) &
219 + dxt(i,3) * ta1(3,j,1) &
220 + dxt(i,4) * ta1(4,j,1) &
221 + dxt(i,5) * ta1(5,j,1) &
222 + dxt(i,6) * ta1(6,j,1) &
223 + dxt(i,7) * ta1(7,j,1) &
224 + dxt(i,8) * ta1(8,j,1) &
225 + dxt(i,9) * ta1(9,j,1) &
226 + dxt(i,10) * ta1(10,j,1) &
227 + dxt(i,11) * ta1(11,j,1) &
228 + dxt(i,12) * ta1(12,j,1) &
229 + dxt(i,13) * ta1(13,j,1)
234 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
240 dtx(i,j,k,e) = dtx(i,j,k,e) &
241 + dyt(j,1) * ta1(i,1,k) &
242 + dyt(j,2) * ta1(i,2,k) &
243 + dyt(j,3) * ta1(i,3,k) &
244 + dyt(j,4) * ta1(i,4,k) &
245 + dyt(j,5) * ta1(i,5,k) &
246 + dyt(j,6) * ta1(i,6,k) &
247 + dyt(j,7) * ta1(i,7,k) &
248 + dyt(j,8) * ta1(i,8,k) &
249 + dyt(j,9) * ta1(i,9,k) &
250 + dyt(j,10) * ta1(i,10,k) &
251 + dyt(j,11) * ta1(i,11,k) &
252 + dyt(j,12) * ta1(i,12,k) &
253 + dyt(j,13) * ta1(i,13,k)
259 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
264 dtx(i,1,k,e) = dtx(i,1,k,e) &
265 + dzt(k,1) * ta1(i,1,1) &
266 + dzt(k,2) * ta1(i,1,2) &
267 + dzt(k,3) * ta1(i,1,3) &
268 + dzt(k,4) * ta1(i,1,4) &
269 + dzt(k,5) * ta1(i,1,5) &
270 + dzt(k,6) * ta1(i,1,6) &
271 + dzt(k,7) * ta1(i,1,7) &
272 + dzt(k,8) * ta1(i,1,8) &
273 + dzt(k,9) * ta1(i,1,9) &
274 + dzt(k,10) * ta1(i,1,10) &
275 + dzt(k,11) * ta1(i,1,11) &
276 + dzt(k,12) * ta1(i,1,12) &
277 + dzt(k,13) * ta1(i,1,13)
284 subroutine cpu_cdtp_lx12(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
285 integer,
parameter :: lx = 12
286 integer,
intent(in) :: nel
287 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
288 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
289 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
290 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
291 integer :: e, i, j, k
296 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
300 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
305 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
306 + dxt(i,2) * ta1(2,j,1) &
307 + dxt(i,3) * ta1(3,j,1) &
308 + dxt(i,4) * ta1(4,j,1) &
309 + dxt(i,5) * ta1(5,j,1) &
310 + dxt(i,6) * ta1(6,j,1) &
311 + dxt(i,7) * ta1(7,j,1) &
312 + dxt(i,8) * ta1(8,j,1) &
313 + dxt(i,9) * ta1(9,j,1) &
314 + dxt(i,10) * ta1(10,j,1) &
315 + dxt(i,11) * ta1(11,j,1) &
316 + dxt(i,12) * ta1(12,j,1)
321 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
327 dtx(i,j,k,e) = dtx(i,j,k,e) &
328 + dyt(j,1) * ta1(i,1,k) &
329 + dyt(j,2) * ta1(i,2,k) &
330 + dyt(j,3) * ta1(i,3,k) &
331 + dyt(j,4) * ta1(i,4,k) &
332 + dyt(j,5) * ta1(i,5,k) &
333 + dyt(j,6) * ta1(i,6,k) &
334 + dyt(j,7) * ta1(i,7,k) &
335 + dyt(j,8) * ta1(i,8,k) &
336 + dyt(j,9) * ta1(i,9,k) &
337 + dyt(j,10) * ta1(i,10,k) &
338 + dyt(j,11) * ta1(i,11,k) &
339 + dyt(j,12) * ta1(i,12,k)
345 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
350 dtx(i,1,k,e) = dtx(i,1,k,e) &
351 + dzt(k,1) * ta1(i,1,1) &
352 + dzt(k,2) * ta1(i,1,2) &
353 + dzt(k,3) * ta1(i,1,3) &
354 + dzt(k,4) * ta1(i,1,4) &
355 + dzt(k,5) * ta1(i,1,5) &
356 + dzt(k,6) * ta1(i,1,6) &
357 + dzt(k,7) * ta1(i,1,7) &
358 + dzt(k,8) * ta1(i,1,8) &
359 + dzt(k,9) * ta1(i,1,9) &
360 + dzt(k,10) * ta1(i,1,10) &
361 + dzt(k,11) * ta1(i,1,11) &
362 + dzt(k,12) * ta1(i,1,12)
369 subroutine cpu_cdtp_lx11(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
370 integer,
parameter :: lx = 11
371 integer,
intent(in) :: nel
372 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
373 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
374 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
375 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
376 integer :: e, i, j, k
381 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
385 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
390 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
391 + dxt(i,2) * ta1(2,j,1) &
392 + dxt(i,3) * ta1(3,j,1) &
393 + dxt(i,4) * ta1(4,j,1) &
394 + dxt(i,5) * ta1(5,j,1) &
395 + dxt(i,6) * ta1(6,j,1) &
396 + dxt(i,7) * ta1(7,j,1) &
397 + dxt(i,8) * ta1(8,j,1) &
398 + dxt(i,9) * ta1(9,j,1) &
399 + dxt(i,10) * ta1(10,j,1) &
400 + dxt(i,11) * ta1(11,j,1)
405 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
411 dtx(i,j,k,e) = dtx(i,j,k,e) &
412 + dyt(j,1) * ta1(i,1,k) &
413 + dyt(j,2) * ta1(i,2,k) &
414 + dyt(j,3) * ta1(i,3,k) &
415 + dyt(j,4) * ta1(i,4,k) &
416 + dyt(j,5) * ta1(i,5,k) &
417 + dyt(j,6) * ta1(i,6,k) &
418 + dyt(j,7) * ta1(i,7,k) &
419 + dyt(j,8) * ta1(i,8,k) &
420 + dyt(j,9) * ta1(i,9,k) &
421 + dyt(j,10) * ta1(i,10,k) &
422 + dyt(j,11) * ta1(i,11,k)
428 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
433 dtx(i,1,k,e) = dtx(i,1,k,e) &
434 + dzt(k,1) * ta1(i,1,1) &
435 + dzt(k,2) * ta1(i,1,2) &
436 + dzt(k,3) * ta1(i,1,3) &
437 + dzt(k,4) * ta1(i,1,4) &
438 + dzt(k,5) * ta1(i,1,5) &
439 + dzt(k,6) * ta1(i,1,6) &
440 + dzt(k,7) * ta1(i,1,7) &
441 + dzt(k,8) * ta1(i,1,8) &
442 + dzt(k,9) * ta1(i,1,9) &
443 + dzt(k,10) * ta1(i,1,10) &
444 + dzt(k,11) * ta1(i,1,11)
451 subroutine cpu_cdtp_lx10(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
452 integer,
parameter :: lx = 10
453 integer,
intent(in) :: nel
454 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
455 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
456 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
457 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
458 integer :: e, i, j, k
463 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
467 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
472 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
473 + dxt(i,2) * ta1(2,j,1) &
474 + dxt(i,3) * ta1(3,j,1) &
475 + dxt(i,4) * ta1(4,j,1) &
476 + dxt(i,5) * ta1(5,j,1) &
477 + dxt(i,6) * ta1(6,j,1) &
478 + dxt(i,7) * ta1(7,j,1) &
479 + dxt(i,8) * ta1(8,j,1) &
480 + dxt(i,9) * ta1(9,j,1) &
481 + dxt(i,10) * ta1(10,j,1)
486 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
492 dtx(i,j,k,e) = dtx(i,j,k,e) &
493 + dyt(j,1) * ta1(i,1,k) &
494 + dyt(j,2) * ta1(i,2,k) &
495 + dyt(j,3) * ta1(i,3,k) &
496 + dyt(j,4) * ta1(i,4,k) &
497 + dyt(j,5) * ta1(i,5,k) &
498 + dyt(j,6) * ta1(i,6,k) &
499 + dyt(j,7) * ta1(i,7,k) &
500 + dyt(j,8) * ta1(i,8,k) &
501 + dyt(j,9) * ta1(i,9,k) &
502 + dyt(j,10) * ta1(i,10,k)
509 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
514 dtx(i,1,k,e) = dtx(i,1,k,e) &
515 + dzt(k,1) * ta1(i,1,1) &
516 + dzt(k,2) * ta1(i,1,2) &
517 + dzt(k,3) * ta1(i,1,3) &
518 + dzt(k,4) * ta1(i,1,4) &
519 + dzt(k,5) * ta1(i,1,5) &
520 + dzt(k,6) * ta1(i,1,6) &
521 + dzt(k,7) * ta1(i,1,7) &
522 + dzt(k,8) * ta1(i,1,8) &
523 + dzt(k,9) * ta1(i,1,9) &
524 + dzt(k,10) * ta1(i,1,10)
531 subroutine cpu_cdtp_lx9(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
532 integer,
parameter :: lx = 9
533 integer,
intent(in) :: nel
534 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
535 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
536 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
537 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
538 integer :: e, i, j, k
543 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
547 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
552 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
553 + dxt(i,2) * ta1(2,j,1) &
554 + dxt(i,3) * ta1(3,j,1) &
555 + dxt(i,4) * ta1(4,j,1) &
556 + dxt(i,5) * ta1(5,j,1) &
557 + dxt(i,6) * ta1(6,j,1) &
558 + dxt(i,7) * ta1(7,j,1) &
559 + dxt(i,8) * ta1(8,j,1) &
560 + dxt(i,9) * ta1(9,j,1)
565 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
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)
586 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
591 dtx(i,1,k,e) = dtx(i,1,k,e) &
592 + dzt(k,1) * ta1(i,1,1) &
593 + dzt(k,2) * ta1(i,1,2) &
594 + dzt(k,3) * ta1(i,1,3) &
595 + dzt(k,4) * ta1(i,1,4) &
596 + dzt(k,5) * ta1(i,1,5) &
597 + dzt(k,6) * ta1(i,1,6) &
598 + dzt(k,7) * ta1(i,1,7) &
599 + dzt(k,8) * ta1(i,1,8) &
600 + dzt(k,9) * ta1(i,1,9)
607 subroutine cpu_cdtp_lx8(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
608 integer,
parameter :: lx = 8
609 integer,
intent(in) :: nel
610 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
611 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
612 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
613 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
614 integer :: e, i, j, k
619 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
623 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
628 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
629 + dxt(i,2) * ta1(2,j,1) &
630 + dxt(i,3) * ta1(3,j,1) &
631 + dxt(i,4) * ta1(4,j,1) &
632 + dxt(i,5) * ta1(5,j,1) &
633 + dxt(i,6) * ta1(6,j,1) &
634 + dxt(i,7) * ta1(7,j,1) &
635 + dxt(i,8) * ta1(8,j,1)
640 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
646 dtx(i,j,k,e) = dtx(i,j,k,e) &
647 + dyt(j,1) * ta1(i,1,k) &
648 + dyt(j,2) * ta1(i,2,k) &
649 + dyt(j,3) * ta1(i,3,k) &
650 + dyt(j,4) * ta1(i,4,k) &
651 + dyt(j,5) * ta1(i,5,k) &
652 + dyt(j,6) * ta1(i,6,k) &
653 + dyt(j,7) * ta1(i,7,k) &
654 + dyt(j,8) * ta1(i,8,k)
660 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
665 dtx(i,1,k,e) = dtx(i,1,k,e) &
666 + dzt(k,1) * ta1(i,1,1) &
667 + dzt(k,2) * ta1(i,1,2) &
668 + dzt(k,3) * ta1(i,1,3) &
669 + dzt(k,4) * ta1(i,1,4) &
670 + dzt(k,5) * ta1(i,1,5) &
671 + dzt(k,6) * ta1(i,1,6) &
672 + dzt(k,7) * ta1(i,1,7) &
673 + dzt(k,8) * ta1(i,1,8)
680 subroutine cpu_cdtp_lx7(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
681 integer,
parameter :: lx = 7
682 integer,
intent(in) :: nel
683 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
684 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
685 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
686 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
687 integer :: e, i, j, k
692 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
696 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
701 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
702 + dxt(i,2) * ta1(2,j,1) &
703 + dxt(i,3) * ta1(3,j,1) &
704 + dxt(i,4) * ta1(4,j,1) &
705 + dxt(i,5) * ta1(5,j,1) &
706 + dxt(i,6) * ta1(6,j,1) &
707 + dxt(i,7) * ta1(7,j,1)
712 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
718 dtx(i,j,k,e) = dtx(i,j,k,e) &
719 + dyt(j,1) * ta1(i,1,k) &
720 + dyt(j,2) * ta1(i,2,k) &
721 + dyt(j,3) * ta1(i,3,k) &
722 + dyt(j,4) * ta1(i,4,k) &
723 + dyt(j,5) * ta1(i,5,k) &
724 + dyt(j,6) * ta1(i,6,k) &
725 + dyt(j,7) * ta1(i,7,k)
731 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
736 dtx(i,1,k,e) = dtx(i,1,k,e) &
737 + dzt(k,1) * ta1(i,1,1) &
738 + dzt(k,2) * ta1(i,1,2) &
739 + dzt(k,3) * ta1(i,1,3) &
740 + dzt(k,4) * ta1(i,1,4) &
741 + dzt(k,5) * ta1(i,1,5) &
742 + dzt(k,6) * ta1(i,1,6) &
743 + dzt(k,7) * ta1(i,1,7)
750 subroutine cpu_cdtp_lx6(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
751 integer,
parameter :: lx = 6
752 integer,
intent(in) :: nel
753 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
754 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
755 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
756 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
757 integer :: e, i, j, k
762 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
766 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
771 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
772 + dxt(i,2) * ta1(2,j,1) &
773 + dxt(i,3) * ta1(3,j,1) &
774 + dxt(i,4) * ta1(4,j,1) &
775 + dxt(i,5) * ta1(5,j,1) &
776 + dxt(i,6) * ta1(6,j,1)
781 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
787 dtx(i,j,k,e) = dtx(i,j,k,e) &
788 + dyt(j,1) * ta1(i,1,k) &
789 + dyt(j,2) * ta1(i,2,k) &
790 + dyt(j,3) * ta1(i,3,k) &
791 + dyt(j,4) * ta1(i,4,k) &
792 + dyt(j,5) * ta1(i,5,k) &
793 + dyt(j,6) * ta1(i,6,k)
799 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
804 dtx(i,1,k,e) = dtx(i,1,k,e) &
805 + dzt(k,1) * ta1(i,1,1) &
806 + dzt(k,2) * ta1(i,1,2) &
807 + dzt(k,3) * ta1(i,1,3) &
808 + dzt(k,4) * ta1(i,1,4) &
809 + dzt(k,5) * ta1(i,1,5) &
810 + dzt(k,6) * ta1(i,1,6)
817 subroutine cpu_cdtp_lx5(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
818 integer,
parameter :: lx = 5
819 integer,
intent(in) :: nel
820 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
821 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
822 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
823 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
824 integer :: e, i, j, k
829 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
833 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
838 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
839 + dxt(i,2) * ta1(2,j,1) &
840 + dxt(i,3) * ta1(3,j,1) &
841 + dxt(i,4) * ta1(4,j,1) &
842 + dxt(i,5) * ta1(5,j,1)
847 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
853 dtx(i,j,k,e) = dtx(i,j,k,e) &
854 + dyt(j,1) * ta1(i,1,k) &
855 + dyt(j,2) * ta1(i,2,k) &
856 + dyt(j,3) * ta1(i,3,k) &
857 + dyt(j,4) * ta1(i,4,k) &
858 + dyt(j,5) * ta1(i,5,k)
864 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
869 dtx(i,1,k,e) = dtx(i,1,k,e) &
870 + dzt(k,1) * ta1(i,1,1) &
871 + dzt(k,2) * ta1(i,1,2) &
872 + dzt(k,3) * ta1(i,1,3) &
873 + dzt(k,4) * ta1(i,1,4) &
874 + dzt(k,5) * ta1(i,1,5)
881 subroutine cpu_cdtp_lx4(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
882 integer,
parameter :: lx = 4
883 integer,
intent(in) :: nel
884 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
885 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
886 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
887 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
888 integer :: e, i, j, k
893 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
897 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
902 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
903 + dxt(i,2) * ta1(2,j,1) &
904 + dxt(i,3) * ta1(3,j,1) &
905 + dxt(i,4) * ta1(4,j,1)
910 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
916 dtx(i,j,k,e) = dtx(i,j,k,e) &
917 + dyt(j,1) * ta1(i,1,k) &
918 + dyt(j,2) * ta1(i,2,k) &
919 + dyt(j,3) * ta1(i,3,k) &
920 + dyt(j,4) * ta1(i,4,k)
926 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
931 dtx(i,1,k,e) = dtx(i,1,k,e) &
932 + dzt(k,1) * ta1(i,1,1) &
933 + dzt(k,2) * ta1(i,1,2) &
934 + dzt(k,3) * ta1(i,1,3) &
935 + dzt(k,4) * ta1(i,1,4)
942 subroutine cpu_cdtp_lx3(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
943 integer,
parameter :: lx = 3
944 integer,
intent(in) :: nel
945 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
946 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
947 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
948 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
949 integer :: e, i, j, k
954 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
958 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
963 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
964 + dxt(i,2) * ta1(2,j,1) &
965 + dxt(i,3) * ta1(3,j,1)
970 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
976 dtx(i,j,k,e) = dtx(i,j,k,e) &
977 + dyt(j,1) * ta1(i,1,k) &
978 + dyt(j,2) * ta1(i,2,k) &
979 + dyt(j,3) * ta1(i,3,k)
985 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
990 dtx(i,1,k,e) = dtx(i,1,k,e) &
991 + dzt(k,1) * ta1(i,1,1) &
992 + dzt(k,2) * ta1(i,1,2) &
993 + dzt(k,3) * ta1(i,1,3)
1000 subroutine cpu_cdtp_lx2(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
1001 integer,
parameter :: lx = 2
1002 integer,
intent(in) :: nel
1003 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(inout) :: dtx
1004 real(kind=
rp),
dimension(lx,lx,lx,nel),
intent(in) :: x, dr, ds, dt, jac, b
1005 real(kind=
rp),
intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
1006 real(kind=
rp),
dimension(lx,lx,lx) :: wx, ta1
1007 integer :: e, i, j, k
1012 wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
1016 ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
1021 dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
1022 + dxt(i,2) * ta1(2,j,1)
1027 ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
1033 dtx(i,j,k,e) = dtx(i,j,k,e) &
1034 + dyt(j,1) * ta1(i,1,k) &
1035 + dyt(j,2) * ta1(i,2,k)
1041 ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
1046 dtx(i,1,k,e) = dtx(i,1,k,e) &
1047 + dzt(k,1) * ta1(i,1,1) &
1048 + dzt(k,2) * ta1(i,1,2)
subroutine cpu_cdtp_lx8(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx5(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx13(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx2(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx9(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, lx)
subroutine cpu_cdtp_lx14(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx10(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx11(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx12(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx7(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx3(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx6(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
subroutine cpu_cdtp_lx4(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
integer, parameter, public rp
Global precision used in computations.