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