Neko  0.8.1
A portable framework for high-order spectral element flow simulations
sx_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 sx_cdtp
35  use num_types, only : rp
36  use math
37  implicit none
38  private
39 
44 
45 contains
46 
47  subroutine sx_cdtp_lx(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd, lx)
48  integer, intent(in) :: nel, nd, lx
49  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
50  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
51  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
52  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
53  real(kind=rp) :: tmp
54  integer :: e, i, j, k, kk, jj
55 
56  call col3(wx, b, x, nd)
57  call invcol2(wx, jac, nd)
58  call col3(ta1, wx, dr, nd)
59 
60  do i = 1, lx
61  do jj = 1, lx * lx * nel
62  tmp = 0d0
63  do kk = 1, lx
64  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
65  end do
66  dtx(i,jj,1,1) = tmp
67  end do
68  end do
69 
70  call col3(ta1, wx, ds, nd)
71 
72  do k = 1, lx
73  do i = 1, lx
74  do j = 1, lx
75  do e = 1, nel
76  tmp = 0d0
77  !NEC$ unroll_completely
78  do kk = 1, lx
79  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
80  end do
81  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
82  end do
83  end do
84  end do
85  end do
86 
87  call col3(ta1, wx, dt, nd)
88 
89  do j = 1, lx
90  do i = 1, lx
91  do k = 1, lx
92  do e = 1, nel
93  tmp = 0d0
94  !NEC$ unroll_completely
95  do kk=1, lx
96  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
97  end do
98  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
99  end do
100  end do
101  end do
102  end do
103 
104  end subroutine sx_cdtp_lx
105 
106  subroutine sx_cdtp_lx14(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
107  integer, parameter :: lx = 14
108  integer, intent(in) :: nel, nd
109  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
110  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
111  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
112  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
113  real(kind=rp) :: tmp
114  integer :: e, i, j, k, kk, jj
115 
116  call col3(wx, b, x, nd)
117  call invcol2(wx, jac, nd)
118  call col3(ta1, wx, dr, nd)
119 
120  do i = 1, lx
121  do jj = 1, lx * lx * nel
122  tmp = 0d0
123  do kk = 1, lx
124  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
125  end do
126  dtx(i,jj,1,1) = tmp
127  end do
128  end do
129 
130  call col3(ta1, wx, ds, nd)
131 
132  do k = 1, lx
133  do i = 1, lx
134  do j = 1, lx
135  do e = 1, nel
136  tmp = 0d0
137  !NEC$ unroll_completely
138  do kk = 1, lx
139  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
140  end do
141  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
142  end do
143  end do
144  end do
145  end do
146 
147  call col3(ta1, wx, dt, nd)
148 
149  do j = 1, lx
150  do i = 1, lx
151  do k = 1, lx
152  do e = 1, nel
153  tmp = 0d0
154  !NEC$ unroll_completely
155  do kk=1, lx
156  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
157  end do
158  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
159  end do
160  end do
161  end do
162  end do
163 
164  end subroutine sx_cdtp_lx14
165 
166  subroutine sx_cdtp_lx13(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
167  integer, parameter :: lx = 13
168  integer, intent(in) :: nel, nd
169  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
170  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
171  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
172  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
173  real(kind=rp) :: tmp
174  integer :: e, i, j, k, kk, jj
175 
176  call col3(wx, b, x, nd)
177  call invcol2(wx, jac, nd)
178  call col3(ta1, wx, dr, nd)
179 
180  do i = 1, lx
181  do jj = 1, lx * lx * nel
182  tmp = 0d0
183  do kk = 1, lx
184  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
185  end do
186  dtx(i,jj,1,1) = tmp
187  end do
188  end do
189 
190  call col3(ta1, wx, ds, nd)
191 
192  do k = 1, lx
193  do i = 1, lx
194  do j = 1, lx
195  do e = 1, nel
196  tmp = 0d0
197  !NEC$ unroll_completely
198  do kk = 1, lx
199  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
200  end do
201  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
202  end do
203  end do
204  end do
205  end do
206 
207  call col3(ta1, wx, dt, nd)
208 
209  do j = 1, lx
210  do i = 1, lx
211  do k = 1, lx
212  do e = 1, nel
213  tmp = 0d0
214  !NEC$ unroll_completely
215  do kk=1, lx
216  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
217  end do
218  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
219  end do
220  end do
221  end do
222  end do
223 
224  end subroutine sx_cdtp_lx13
225 
226  subroutine sx_cdtp_lx12(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
227  integer, parameter :: lx = 12
228  integer, intent(in) :: nel, nd
229  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
230  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
231  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
232  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
233  real(kind=rp) :: tmp
234  integer :: e, i, j, k, kk, jj
235 
236  call col3(wx, b, x, nd)
237  call invcol2(wx, jac, nd)
238  call col3(ta1, wx, dr, nd)
239 
240  do i = 1, lx
241  do jj = 1, lx * lx * nel
242  tmp = 0d0
243  do kk = 1, lx
244  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
245  end do
246  dtx(i,jj,1,1) = tmp
247  end do
248  end do
249 
250  call col3(ta1, wx, ds, nd)
251 
252  do k = 1, lx
253  do i = 1, lx
254  do j = 1, lx
255  do e = 1, nel
256  tmp = 0d0
257  !NEC$ unroll_completely
258  do kk = 1, lx
259  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
260  end do
261  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
262  end do
263  end do
264  end do
265  end do
266 
267  call col3(ta1, wx, dt, nd)
268 
269  do j = 1, lx
270  do i = 1, lx
271  do k = 1, lx
272  do e = 1, nel
273  tmp = 0d0
274  !NEC$ unroll_completely
275  do kk=1, lx
276  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
277  end do
278  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
279  end do
280  end do
281  end do
282  end do
283 
284  end subroutine sx_cdtp_lx12
285 
286  subroutine sx_cdtp_lx11(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
287  integer, parameter :: lx = 11
288  integer, intent(in) :: nel, nd
289  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
290  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
291  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
292  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
293  real(kind=rp) :: tmp
294  integer :: e, i, j, k, kk, jj
295 
296  call col3(wx, b, x, nd)
297  call invcol2(wx, jac, nd)
298  call col3(ta1, wx, dr, nd)
299 
300  do i = 1, lx
301  do jj = 1, lx * lx * nel
302  tmp = 0d0
303  do kk = 1, lx
304  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
305  end do
306  dtx(i,jj,1,1) = tmp
307  end do
308  end do
309 
310  call col3(ta1, wx, ds, nd)
311 
312  do k = 1, lx
313  do i = 1, lx
314  do j = 1, lx
315  do e = 1, nel
316  tmp = 0d0
317  !NEC$ unroll_completely
318  do kk = 1, lx
319  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
320  end do
321  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
322  end do
323  end do
324  end do
325  end do
326 
327  call col3(ta1, wx, dt, nd)
328 
329  do j = 1, lx
330  do i = 1, lx
331  do k = 1, lx
332  do e = 1, nel
333  tmp = 0d0
334  !NEC$ unroll_completely
335  do kk = 1, lx
336  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
337  end do
338  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
339  end do
340  end do
341  end do
342  end do
343 
344  end subroutine sx_cdtp_lx11
345 
346  subroutine sx_cdtp_lx10(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
347  integer, parameter :: lx = 10
348  integer, intent(in) :: nel, nd
349  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
350  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
351  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
352  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
353  real(kind=rp) :: tmp
354  integer :: e, i, j, k, kk, jj
355 
356  call col3(wx, b, x, nd)
357  call invcol2(wx, jac, nd)
358  call col3(ta1, wx, dr, nd)
359 
360  do i = 1, lx
361  do jj = 1, lx * lx * nel
362  tmp = 0d0
363  do kk = 1, lx
364  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
365  end do
366  dtx(i,jj,1,1) = tmp
367  end do
368  end do
369 
370  call col3(ta1, wx, ds, nd)
371 
372  do k = 1, lx
373  do i = 1, lx
374  do j = 1, lx
375  do e = 1, nel
376  tmp = 0d0
377  !NEC$ unroll_completely
378  do kk = 1, lx
379  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
380  end do
381  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
382  end do
383  end do
384  end do
385  end do
386 
387  call col3(ta1, wx, dt, nd)
388 
389  do j = 1, lx
390  do i = 1, lx
391  do k = 1, lx
392  do e = 1, nel
393  tmp = 0d0
394  !NEC$ unroll_completely
395  do kk=1, lx
396  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
397  end do
398  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
399  end do
400  end do
401  end do
402  end do
403 
404  end subroutine sx_cdtp_lx10
405 
406  subroutine sx_cdtp_lx9(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
407  integer, parameter :: lx = 9
408  integer, intent(in) :: nel, nd
409  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
410  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
411  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
412  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
413  real(kind=rp) :: tmp
414  integer :: e, i, j, k, kk, jj
415 
416  call col3(wx, b, x, nd)
417  call invcol2(wx, jac, nd)
418  call col3(ta1, wx, dr, nd)
419 
420  do i = 1, lx
421  do jj = 1, lx * lx * nel
422  tmp = 0d0
423  do kk = 1, lx
424  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
425  end do
426  dtx(i,jj,1,1) = tmp
427  end do
428  end do
429 
430  call col3(ta1, wx, ds, nd)
431 
432  do k = 1, lx
433  do i = 1, lx
434  do j = 1, lx
435  do e = 1, nel
436  tmp = 0d0
437  !NEC$ unroll_completely
438  do kk=1, lx
439  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
440  end do
441  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
442  end do
443  end do
444  end do
445  end do
446 
447  call col3(ta1, wx, dt, nd)
448 
449  do j = 1, lx
450  do i = 1, lx
451  do k = 1, lx
452  do e = 1, nel
453  tmp = 0d0
454  !NEC$ unroll_completely
455  do kk=1, lx
456  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
457  end do
458  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
459  end do
460  end do
461  end do
462  end do
463 
464  end subroutine sx_cdtp_lx9
465 
466  subroutine sx_cdtp_lx8(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
467  integer, parameter :: lx = 8
468  integer, intent(in) :: nel, nd
469  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
470  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
471  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
472  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
473  real(kind=rp) :: tmp
474  integer :: e, i, j, k, kk, jj
475 
476  call col3(wx, b, x, nd)
477  call invcol2(wx, jac, nd)
478  call col3(ta1, wx, dr, nd)
479 
480  do i = 1, lx
481  do jj = 1, lx * lx * nel
482  tmp = 0d0
483  do kk = 1, lx
484  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
485  end do
486  dtx(i,jj,1,1) = tmp
487  end do
488  end do
489 
490  call col3(ta1, wx, ds, nd)
491 
492  do k = 1, lx
493  do i = 1, lx
494  do j = 1, lx
495  do e = 1, nel
496  tmp = 0d0
497  !NEC$ unroll_completely
498  do kk = 1, lx
499  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
500  end do
501  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
502  end do
503  end do
504  end do
505  end do
506 
507  call col3(ta1, wx, dt, nd)
508 
509  do j = 1, lx
510  do i = 1, lx
511  do k = 1, lx
512  do e = 1, nel
513  tmp = 0d0
514  !NEC$ unroll_completely
515  do kk=1, lx
516  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
517  end do
518  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
519  end do
520  end do
521  end do
522  end do
523 
524  end subroutine sx_cdtp_lx8
525 
526  subroutine sx_cdtp_lx7(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
527  integer, parameter :: lx = 7
528  integer, intent(in) :: nel, nd
529  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
530  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
531  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
532  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
533  real(kind=rp) :: tmp
534  integer :: e, i, j, k, kk, jj
535 
536  call col3(wx, b, x, nd)
537  call invcol2(wx, jac, nd)
538  call col3(ta1, wx, dr, nd)
539 
540  do i = 1, lx
541  do jj = 1, lx * lx * nel
542  tmp = 0d0
543  do kk = 1, lx
544  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
545  end do
546  dtx(i,jj,1,1) = tmp
547  end do
548  end do
549 
550  call col3(ta1, wx, ds, nd)
551 
552  do k = 1, lx
553  do i = 1, lx
554  do j = 1, lx
555  do e = 1, nel
556  tmp = 0d0
557  !NEC$ unroll_completely
558  do kk=1, lx
559  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
560  end do
561  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
562  end do
563  end do
564  end do
565  end do
566 
567  call col3(ta1, wx, dt, nd)
568 
569  do j = 1, lx
570  do i = 1, lx
571  do k = 1, lx
572  do e = 1, nel
573  tmp = 0d0
574  !NEC$ unroll_completely
575  do kk=1, lx
576  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
577  end do
578  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
579  end do
580  end do
581  end do
582  end do
583 
584  end subroutine sx_cdtp_lx7
585 
586  subroutine sx_cdtp_lx6(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
587  integer, parameter :: lx = 6
588  integer, intent(in) :: nel, nd
589  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
590  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
591  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
592  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
593  real(kind=rp) :: tmp
594  integer :: e, i, j, k, kk, jj
595 
596  call col3(wx, b, x, nd)
597  call invcol2(wx, jac, nd)
598  call col3(ta1, wx, dr, nd)
599 
600  do i = 1, lx
601  do jj = 1, lx * lx * nel
602  tmp = 0d0
603  do kk = 1, lx
604  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
605  end do
606  dtx(i,jj,1,1) = tmp
607  end do
608  end do
609 
610  call col3(ta1, wx, ds, nd)
611 
612  do k = 1, lx
613  do i = 1, lx
614  do j = 1, lx
615  do e = 1, nel
616  tmp = 0d0
617  !NEC$ unroll_completely
618  do kk=1, lx
619  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
620  end do
621  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
622  end do
623  end do
624  end do
625  end do
626 
627  call col3(ta1, wx, dt, nd)
628 
629  do j = 1, lx
630  do i = 1, lx
631  do k = 1, lx
632  do e = 1, nel
633  tmp = 0d0
634  !NEC$ unroll_completely
635  do kk = 1, lx
636  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
637  end do
638  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
639  end do
640  end do
641  end do
642  end do
643 
644  end subroutine sx_cdtp_lx6
645 
646  subroutine sx_cdtp_lx5(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
647  integer, parameter :: lx = 5
648  integer, intent(in) :: nel, nd
649  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
650  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
651  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
652  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
653  real(kind=rp) :: tmp
654  integer :: e, i, j, k, kk, jj
655 
656  call col3(wx, b, x, nd)
657  call invcol2(wx, jac, nd)
658  call col3(ta1, wx, dr, nd)
659 
660  do i = 1, lx
661  do jj = 1, lx * lx * nel
662  tmp = 0d0
663  do kk = 1, lx
664  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
665  end do
666  dtx(i,jj,1,1) = tmp
667  end do
668  end do
669 
670  call col3(ta1, wx, ds, nd)
671 
672  do k = 1, lx
673  do i = 1, lx
674  do j = 1, lx
675  do e = 1, nel
676  tmp = 0d0
677  !NEC$ unroll_completely
678  do kk = 1, lx
679  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
680  end do
681  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
682  end do
683  end do
684  end do
685  end do
686 
687  call col3(ta1, wx, dt, nd)
688 
689  do j = 1, lx
690  do i = 1, lx
691  do k = 1, lx
692  do e = 1, nel
693  tmp = 0d0
694  !NEC$ unroll_completely
695  do kk=1, lx
696  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
697  end do
698  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
699  end do
700  end do
701  end do
702  end do
703 
704  end subroutine sx_cdtp_lx5
705 
706  subroutine sx_cdtp_lx4(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
707  integer, parameter :: lx = 4
708  integer, intent(in) :: nel, nd
709  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
710  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
711  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
712  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
713  real(kind=rp) :: tmp
714  integer :: e, i, j, k, kk, jj
715 
716  call col3(wx, b, x, nd)
717  call invcol2(wx, jac, nd)
718  call col3(ta1, wx, dr, nd)
719 
720  do i = 1, lx
721  do jj = 1, lx * lx * nel
722  tmp = 0d0
723  do kk=1, lx
724  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
725  end do
726  dtx(i,jj,1,1) = tmp
727  end do
728  end do
729 
730  call col3(ta1, wx, ds, nd)
731 
732  do k = 1, lx
733  do i = 1, lx
734  do j = 1, lx
735  do e = 1, nel
736  tmp = 0d0
737  !NEC$ unroll_completely
738  do kk=1, lx
739  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
740  end do
741  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
742  end do
743  end do
744  end do
745  end do
746 
747  call col3(ta1, wx, dt, nd)
748 
749  do j = 1, lx
750  do i = 1, lx
751  do k = 1, lx
752  do e = 1, nel
753  tmp = 0d0
754  !NEC$ unroll_completely
755  do kk = 1, lx
756  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
757  end do
758  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
759  end do
760  end do
761  end do
762  end do
763 
764  end subroutine sx_cdtp_lx4
765 
766  subroutine sx_cdtp_lx3(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
767  integer, parameter :: lx = 3
768  integer, intent(in) :: nel, nd
769  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
770  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
771  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
772  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
773  real(kind=rp) :: tmp
774  integer :: e, i, j, k, kk, jj
775 
776  call col3(wx, b, x, nd)
777  call invcol2(wx, jac, nd)
778  call col3(ta1, wx, dr, nd)
779 
780  do i = 1, lx
781  do jj = 1, lx * lx * nel
782  tmp = 0d0
783  do kk = 1, lx
784  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
785  end do
786  dtx(i,jj,1,1) = tmp
787  end do
788  end do
789 
790  call col3(ta1, wx, ds, nd)
791 
792  do k = 1, lx
793  do i = 1, lx
794  do j = 1, lx
795  do e = 1, nel
796  tmp = 0d0
797  !NEC$ unroll_completely
798  do kk = 1, lx
799  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
800  end do
801  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
802  end do
803  end do
804  end do
805  end do
806 
807  call col3(ta1, wx, dt, nd)
808 
809  do j = 1, lx
810  do i = 1, lx
811  do k = 1, lx
812  do e = 1, nel
813  tmp = 0d0
814  !NEC$ unroll_completely
815  do kk = 1, lx
816  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
817  end do
818  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
819  end do
820  end do
821  end do
822  end do
823 
824  end subroutine sx_cdtp_lx3
825 
826  subroutine sx_cdtp_lx2(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
827  integer, parameter :: lx = 2
828  integer, intent(in) :: nel, nd
829  real(kind=rp), dimension(lx,lx,lx,nel), intent(inout) :: dtx
830  real(kind=rp), dimension(lx,lx,lx,nel), intent(in) :: x, dr, ds, dt, jac, b
831  real(kind=rp), intent(in) :: dxt(lx,lx), dyt(lx,lx), dzt(lx,lx)
832  real(kind=rp), dimension(lx,lx,lx,nel) :: wx, ta1
833  real(kind=rp) :: tmp
834  integer :: e, i, j, k, kk, jj
835 
836  call col3(wx, b, x, nd)
837  call invcol2(wx, jac, nd)
838  call col3(ta1, wx, dr, nd)
839 
840  do i = 1, lx
841  do jj = 1, lx * lx * nel
842  tmp = 0d0
843  do kk = 1, lx
844  tmp = tmp + dxt(i,kk) * ta1(kk,jj,1,1)
845  end do
846  dtx(i,jj,1,1) = tmp
847  end do
848  end do
849 
850  call col3(ta1, wx, ds, nd)
851 
852  do k = 1, lx
853  do i = 1, lx
854  do j = 1, lx
855  do e = 1, nel
856  tmp = 0d0
857  !NEC$ unroll_completely
858  do kk = 1, lx
859  tmp = tmp + dyt(j, kk) * ta1(i,kk,k,e)
860  end do
861  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
862  end do
863  end do
864  end do
865  end do
866 
867  call col3(ta1, wx, dt, nd)
868 
869  do j = 1, lx
870  do i = 1, lx
871  do k = 1, lx
872  do e = 1, nel
873  tmp = 0d0
874  !NEC$ unroll_completely
875  do kk = 1, lx
876  tmp = tmp + dzt(k, kk)*ta1(i,j,kk,e)
877  end do
878  dtx(i,j,k,e) = dtx(i,j,k,e) + tmp
879  end do
880  end do
881  end do
882  end do
883 
884  end subroutine sx_cdtp_lx2
885 
886 
887 end module sx_cdtp
Definition: math.f90:60
subroutine, public invcol2(a, b, n)
Vector division .
Definition: math.f90:631
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition: math.f90:658
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
DT*X kernels for SX-Aurora.
Definition: sx_cdtp.f90:34
subroutine, public sx_cdtp_lx2(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:827
subroutine, public sx_cdtp_lx8(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:467
subroutine, public sx_cdtp_lx9(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:407
subroutine, public sx_cdtp_lx6(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:587
subroutine, public sx_cdtp_lx12(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:227
subroutine, public sx_cdtp_lx7(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:527
subroutine, public sx_cdtp_lx10(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:347
subroutine, public sx_cdtp_lx11(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:287
subroutine, public sx_cdtp_lx4(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:707
subroutine, public sx_cdtp_lx3(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:767
subroutine, public sx_cdtp_lx(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd, lx)
Definition: sx_cdtp.f90:48
subroutine, public sx_cdtp_lx14(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:107
subroutine, public sx_cdtp_lx5(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:647
subroutine, public sx_cdtp_lx13(dtx, x, dr, ds, dt, dxt, dyt, dzt, B, jac, nel, nd)
Definition: sx_cdtp.f90:167