Neko  0.8.1
A portable framework for high-order spectral element flow simulations
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 !
34 module cpu_cdtp
35  use num_types, only : rp
36  implicit none
37 
38 contains
39 
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
46  real(kind=rp) :: tmp
47  integer :: e, i, j, k, l
48 
49  do e = 1, nel
50 
51  do i = 1, lx*lx*lx
52  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
53  end do
54 
55  do i = 1, lx*lx*lx
56  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
57  end do
58 
59  do j = 1, lx * lx
60  do i = 1, lx
61  tmp = 0.0_rp
62  !DIR$ LOOP_INFO MIN_TRIPS(15)
63  do k = 1, lx
64  tmp = tmp + dxt(i,k) * ta1(k,j,1)
65  end do
66  dtx(i,j,1,e) = tmp
67  end do
68  end do
69 
70  do i = 1, lx*lx*lx
71  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
72  end do
73 
74  do k = 1, lx
75  do j = 1, lx
76  do i = 1, lx
77  tmp = 0.0_rp
78  !DIR$ LOOP_INFO MIN_TRIPS(15)
79  do l = 1, lx
80  tmp = tmp + dyt(j,l) * ta1(i,l,k)
81  end do
82  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
83  end do
84  end do
85  end do
86 
87  do i = 1, lx*lx*lx
88  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
89  end do
90 
91  do k = 1, lx
92  do i = 1, lx*lx
93  tmp = 0.0_rp
94  !DIR$ LOOP_INFO MIN_TRIPS(15)
95  do l = 1, lx
96  tmp = tmp + dzt(k,l) * ta1(i,1,l)
97  end do
98  dtx(i,1,k,e) = dtx(i,1,k,e) + tmp
99  end do
100  end do
101 
102  end do
103  end subroutine cpu_cdtp_lx
104 
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
113 
114  do e = 1, nel
115 
116  do i = 1, lx*lx*lx
117  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
118  end do
119 
120  do i = 1, lx*lx*lx
121  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
122  end do
123 
124  do j = 1, lx * lx
125  do i = 1, lx
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)
140  end do
141  end do
142 
143  do i = 1, lx*lx*lx
144  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
145  end do
146 
147  do k = 1, lx
148  do j = 1, lx
149  do i = 1, lx
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)
165  end do
166  end do
167  end do
168 
169  do i = 1, lx*lx*lx
170  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
171  end do
172 
173  do k = 1, lx
174  do i = 1, lx*lx
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)
190  end do
191  end do
192 
193  end do
194  end subroutine cpu_cdtp_lx14
195 
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
204 
205  do e = 1, nel
206 
207  do i = 1, lx*lx*lx
208  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
209  end do
210 
211  do i = 1, lx*lx*lx
212  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
213  end do
214 
215  do j = 1, lx * lx
216  do i = 1, lx
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)
230  end do
231  end do
232 
233  do i = 1, lx*lx*lx
234  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
235  end do
236 
237  do k = 1, lx
238  do j = 1, lx
239  do i = 1, lx
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)
254  end do
255  end do
256  end do
257 
258  do i = 1, lx*lx*lx
259  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
260  end do
261 
262  do k = 1, lx
263  do i = 1, lx*lx
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)
278  end do
279  end do
280 
281  end do
282  end subroutine cpu_cdtp_lx13
283 
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
292 
293  do e = 1, nel
294 
295  do i = 1, lx*lx*lx
296  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
297  end do
298 
299  do i = 1, lx*lx*lx
300  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
301  end do
302 
303  do j = 1, lx * lx
304  do i = 1, lx
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)
317  end do
318  end do
319 
320  do i = 1, lx*lx*lx
321  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
322  end do
323 
324  do k = 1, lx
325  do j = 1, lx
326  do i = 1, lx
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)
340  end do
341  end do
342  end do
343 
344  do i = 1, lx*lx*lx
345  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
346  end do
347 
348  do k = 1, lx
349  do i = 1, lx*lx
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)
363  end do
364  end do
365 
366  end do
367  end subroutine cpu_cdtp_lx12
368 
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
377 
378  do e = 1, nel
379 
380  do i = 1, lx*lx*lx
381  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
382  end do
383 
384  do i = 1, lx*lx*lx
385  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
386  end do
387 
388  do j = 1, lx * lx
389  do i = 1, lx
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)
401  end do
402  end do
403 
404  do i = 1, lx*lx*lx
405  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
406  end do
407 
408  do k = 1, lx
409  do j = 1, lx
410  do i = 1, lx
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)
423  end do
424  end do
425  end do
426 
427  do i = 1, lx*lx*lx
428  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
429  end do
430 
431  do k = 1, lx
432  do i = 1, lx*lx
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)
445  end do
446  end do
447 
448  end do
449  end subroutine cpu_cdtp_lx11
450 
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
459 
460  do e = 1, nel
461 
462  do i = 1, lx*lx*lx
463  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
464  end do
465 
466  do i = 1, lx*lx*lx
467  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
468  end do
469 
470  do j = 1, lx * lx
471  do i = 1, lx
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)
482  end do
483  end do
484 
485  do i = 1, lx*lx*lx
486  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
487  end do
488 
489  do k = 1, lx
490  do j = 1, lx
491  do i = 1, lx
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)
503 
504  end do
505  end do
506  end do
507 
508  do i = 1, lx*lx*lx
509  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
510  end do
511 
512  do k = 1, lx
513  do i = 1, lx*lx
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)
525  end do
526  end do
527 
528  end do
529  end subroutine cpu_cdtp_lx10
530 
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
539 
540  do e = 1, nel
541 
542  do i = 1, lx*lx*lx
543  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
544  end do
545 
546  do i = 1, lx*lx*lx
547  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
548  end do
549 
550  do j = 1, lx * lx
551  do i = 1, lx
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)
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  end do
582  end do
583  end do
584 
585  do i = 1, lx*lx*lx
586  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
587  end do
588 
589  do k = 1, lx
590  do i = 1, lx*lx
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)
601  end do
602  end do
603 
604  end do
605  end subroutine cpu_cdtp_lx9
606 
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
615 
616  do e = 1, nel
617 
618  do i = 1, lx*lx*lx
619  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
620  end do
621 
622  do i = 1, lx*lx*lx
623  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
624  end do
625 
626  do j = 1, lx * lx
627  do i = 1, lx
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)
636  end do
637  end do
638 
639  do i = 1, lx*lx*lx
640  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
641  end do
642 
643  do k = 1, lx
644  do j = 1, lx
645  do i = 1, lx
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)
655  end do
656  end do
657  end do
658 
659  do i = 1, lx*lx*lx
660  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
661  end do
662 
663  do k = 1, lx
664  do i = 1, lx*lx
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)
674  end do
675  end do
676 
677  end do
678  end subroutine cpu_cdtp_lx8
679 
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
688 
689  do e = 1, nel
690 
691  do i = 1, lx*lx*lx
692  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
693  end do
694 
695  do i = 1, lx*lx*lx
696  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
697  end do
698 
699  do j = 1, lx * lx
700  do i = 1, lx
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)
708  end do
709  end do
710 
711  do i = 1, lx*lx*lx
712  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
713  end do
714 
715  do k = 1, lx
716  do j = 1, lx
717  do i = 1, lx
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)
726  end do
727  end do
728  end do
729 
730  do i = 1, lx*lx*lx
731  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
732  end do
733 
734  do k = 1, lx
735  do i = 1, lx*lx
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)
744  end do
745  end do
746 
747  end do
748  end subroutine cpu_cdtp_lx7
749 
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
758 
759  do e = 1, nel
760 
761  do i = 1, lx*lx*lx
762  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
763  end do
764 
765  do i = 1, lx*lx*lx
766  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
767  end do
768 
769  do j = 1, lx * lx
770  do i = 1, lx
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)
777  end do
778  end do
779 
780  do i = 1, lx*lx*lx
781  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
782  end do
783 
784  do k = 1, lx
785  do j = 1, lx
786  do i = 1, lx
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)
794  end do
795  end do
796  end do
797 
798  do i = 1, lx*lx*lx
799  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
800  end do
801 
802  do k = 1, lx
803  do i = 1, lx*lx
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)
811  end do
812  end do
813 
814  end do
815  end subroutine cpu_cdtp_lx6
816 
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
825 
826  do e = 1, nel
827 
828  do i = 1, lx*lx*lx
829  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
830  end do
831 
832  do i = 1, lx*lx*lx
833  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
834  end do
835 
836  do j = 1, lx * lx
837  do i = 1, lx
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)
843  end do
844  end do
845 
846  do i = 1, lx*lx*lx
847  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
848  end do
849 
850  do k = 1, lx
851  do j = 1, lx
852  do i = 1, lx
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)
859  end do
860  end do
861  end do
862 
863  do i = 1, lx*lx*lx
864  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
865  end do
866 
867  do k = 1, lx
868  do i = 1, lx*lx
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)
875  end do
876  end do
877 
878  end do
879  end subroutine cpu_cdtp_lx5
880 
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
889 
890  do e = 1, nel
891 
892  do i = 1, lx*lx*lx
893  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
894  end do
895 
896  do i = 1, lx*lx*lx
897  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
898  end do
899 
900  do j = 1, lx * lx
901  do i = 1, lx
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)
906  end do
907  end do
908 
909  do i = 1, lx*lx*lx
910  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
911  end do
912 
913  do k = 1, lx
914  do j = 1, lx
915  do i = 1, lx
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)
921  end do
922  end do
923  end do
924 
925  do i = 1, lx*lx*lx
926  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
927  end do
928 
929  do k = 1, lx
930  do i = 1, lx*lx
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)
936  end do
937  end do
938 
939  end do
940  end subroutine cpu_cdtp_lx4
941 
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
950 
951  do e = 1, nel
952 
953  do i = 1, lx*lx*lx
954  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
955  end do
956 
957  do i = 1, lx*lx*lx
958  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
959  end do
960 
961  do j = 1, lx * lx
962  do i = 1, lx
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)
966  end do
967  end do
968 
969  do i = 1, lx*lx*lx
970  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
971  end do
972 
973  do k = 1, lx
974  do j = 1, lx
975  do i = 1, lx
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)
980  end do
981  end do
982  end do
983 
984  do i = 1, lx*lx*lx
985  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
986  end do
987 
988  do k = 1, lx
989  do i = 1, lx*lx
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)
994  end do
995  end do
996 
997  end do
998  end subroutine cpu_cdtp_lx3
999 
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
1008 
1009  do e = 1, nel
1010 
1011  do i = 1, lx*lx*lx
1012  wx(i,1,1) = ( b(i,1,1,e) * x(i,1,1,e) ) / jac(i,1,1,e)
1013  end do
1014 
1015  do i = 1, lx*lx*lx
1016  ta1(i,1,1) = wx(i,1,1) * dr(i,1,1,e)
1017  end do
1018 
1019  do j = 1, lx * lx
1020  do i = 1, lx
1021  dtx(i,j,1,e) = dxt(i,1) * ta1(1,j,1) &
1022  + dxt(i,2) * ta1(2,j,1)
1023  end do
1024  end do
1025 
1026  do i = 1, lx*lx*lx
1027  ta1(i,1,1) = wx(i,1,1) * ds(i,1,1,e)
1028  end do
1029 
1030  do k = 1, lx
1031  do j = 1, lx
1032  do i = 1, lx
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)
1036  end do
1037  end do
1038  end do
1039 
1040  do i = 1, lx*lx*lx
1041  ta1(i,1,1) = wx(i,1,1) * dt(i,1,1,e)
1042  end do
1043 
1044  do k = 1, lx
1045  do i = 1, lx*lx
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)
1049  end do
1050  end do
1051 
1052  end do
1053  end subroutine cpu_cdtp_lx2
1054 
1055 end module cpu_cdtp
DT*X kernels.
Definition: cdtp.f90:34
subroutine cpu_cdtp_lx8(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:608
subroutine cpu_cdtp_lx5(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:818
subroutine cpu_cdtp_lx13(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:197
subroutine cpu_cdtp_lx2(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:1001
subroutine cpu_cdtp_lx9(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:532
subroutine cpu_cdtp_lx(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, lx)
Definition: cdtp.f90:41
subroutine cpu_cdtp_lx14(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:106
subroutine cpu_cdtp_lx10(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:452
subroutine cpu_cdtp_lx11(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:370
subroutine cpu_cdtp_lx12(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:285
subroutine cpu_cdtp_lx7(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:681
subroutine cpu_cdtp_lx3(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:943
subroutine cpu_cdtp_lx6(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:751
subroutine cpu_cdtp_lx4(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel)
Definition: cdtp.f90:882
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12