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