Neko  0.8.1
A portable framework for high-order spectral element flow simulations
opgrad.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 module cpu_opgrad
35  use num_types, only : rp
36  implicit none
37 
38 contains
39 
40  subroutine cpu_opgrad_lx(ux, uy, uz, u, dx, dy, dz, &
41  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n, lx)
42  integer, intent(in) :: n, lx
43  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
44  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
45  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
46  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
47  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
48  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
49  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
50  real(kind=rp) :: ur(lx,lx,lx)
51  real(kind=rp) :: us(lx,lx,lx)
52  real(kind=rp) :: ut(lx,lx,lx)
53  real(kind=rp) :: tmp
54  integer :: e, i, j, k, l
55 
56  do e = 1, n
57  do j = 1, lx * lx
58  do i = 1, lx
59  tmp = 0.0_rp
60  do k = 1, lx
61  tmp = tmp + dx(i,k) * u(k,j,1,e)
62  end do
63  ur(i,j,1) = tmp
64  end do
65  end do
66 
67  do k = 1, lx
68  do j = 1, lx
69  do i = 1, lx
70  tmp = 0.0_rp
71  do l = 1, lx
72  tmp = tmp + dy(j,l) * u(i,l,k,e)
73  end do
74  us(i,j,k) = tmp
75  end do
76  end do
77  end do
78 
79  do k = 1, lx
80  do i = 1, lx*lx
81  tmp = 0.0_rp
82  do l = 1, lx
83  tmp = tmp + dz(k,l) * u(i,1,l,e)
84  end do
85  ut(i,1,k) = tmp
86  end do
87  end do
88 
89  do i = 1, lx * lx * lx
90  ux(i,1,1,e) = w3(i,1,1) &
91  * ( drdx(i,1,1,e) * ur(i,1,1) &
92  + dsdx(i,1,1,e) * us(i,1,1) &
93  + dtdx(i,1,1,e) * ut(i,1,1) )
94  uy(i,1,1,e) = w3(i,1,1) &
95  * ( dsdy(i,1,1,e) * us(i,1,1) &
96  + drdy(i,1,1,e) * ur(i,1,1) &
97  + dtdy(i,1,1,e) * ut(i,1,1) )
98  uz(i,1,1,e) = w3(i,1,1) &
99  * ( dtdz(i,1,1,e) * ut(i,1,1) &
100  + drdz(i,1,1,e) * ur(i,1,1) &
101  + dsdz(i,1,1,e) * us(i,1,1) )
102  end do
103  end do
104  end subroutine cpu_opgrad_lx
105 
106  subroutine cpu_opgrad_lx18(ux, uy, uz, u, dx, dy, dz, &
107  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
108  integer, parameter :: lx = 18
109  integer, intent(in) :: n
110  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
111  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
112  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
113  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
114  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
115  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
116  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
117  real(kind=rp) :: ur(lx,lx,lx)
118  real(kind=rp) :: us(lx,lx,lx)
119  real(kind=rp) :: ut(lx,lx,lx)
120  integer :: e, i, j, k
121 
122  do e = 1, n
123  do j = 1, lx * lx
124  do i = 1, lx
125  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
126  + dx(i,2) * u(2,j,1,e) &
127  + dx(i,3) * u(3,j,1,e) &
128  + dx(i,4) * u(4,j,1,e) &
129  + dx(i,5) * u(5,j,1,e) &
130  + dx(i,6) * u(6,j,1,e) &
131  + dx(i,7) * u(7,j,1,e) &
132  + dx(i,8) * u(8,j,1,e) &
133  + dx(i,9) * u(9,j,1,e) &
134  + dx(i,10) * u(10,j,1,e) &
135  + dx(i,11) * u(11,j,1,e) &
136  + dx(i,12) * u(12,j,1,e) &
137  + dx(i,13) * u(13,j,1,e) &
138  + dx(i,14) * u(14,j,1,e) &
139  + dx(i,15) * u(15,j,1,e) &
140  + dx(i,16) * u(16,j,1,e) &
141  + dx(i,17) * u(17,j,1,e) &
142  + dx(i,18) * u(18,j,1,e)
143  end do
144  end do
145 
146  do k = 1, lx
147  do j = 1, lx
148  do i = 1, lx
149  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
150  + dy(j,2) * u(i,2,k,e) &
151  + dy(j,3) * u(i,3,k,e) &
152  + dy(j,4) * u(i,4,k,e) &
153  + dy(j,5) * u(i,5,k,e) &
154  + dy(j,6) * u(i,6,k,e) &
155  + dy(j,7) * u(i,7,k,e) &
156  + dy(j,8) * u(i,8,k,e) &
157  + dy(j,9) * u(i,9,k,e) &
158  + dy(j,10) * u(i,10,k,e) &
159  + dy(j,11) * u(i,11,k,e) &
160  + dy(j,12) * u(i,12,k,e) &
161  + dy(j,13) * u(i,13,k,e) &
162  + dy(j,14) * u(i,14,k,e) &
163  + dy(j,15) * u(i,15,k,e) &
164  + dy(j,16) * u(i,16,k,e) &
165  + dy(j,17) * u(i,17,k,e) &
166  + dy(j,18) * u(i,18,k,e)
167  end do
168  end do
169  end do
170 
171  do k = 1, lx
172  do i = 1, lx*lx
173  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
174  + dz(k,2) * u(i,1,2,e) &
175  + dz(k,3) * u(i,1,3,e) &
176  + dz(k,4) * u(i,1,4,e) &
177  + dz(k,5) * u(i,1,5,e) &
178  + dz(k,6) * u(i,1,6,e) &
179  + dz(k,7) * u(i,1,7,e) &
180  + dz(k,8) * u(i,1,8,e) &
181  + dz(k,9) * u(i,1,9,e) &
182  + dz(k,10) * u(i,1,10,e) &
183  + dz(k,11) * u(i,1,11,e) &
184  + dz(k,12) * u(i,1,12,e) &
185  + dz(k,13) * u(i,1,13,e) &
186  + dz(k,14) * u(i,1,14,e) &
187  + dz(k,15) * u(i,1,15,e) &
188  + dz(k,16) * u(i,1,16,e) &
189  + dz(k,17) * u(i,1,17,e) &
190  + dz(k,18) * u(i,1,18,e)
191  end do
192  end do
193 
194  do i = 1, lx * lx * lx
195  ux(i,1,1,e) = w3(i,1,1) &
196  * ( drdx(i,1,1,e) * ur(i,1,1) &
197  + dsdx(i,1,1,e) * us(i,1,1) &
198  + dtdx(i,1,1,e) * ut(i,1,1) )
199  uy(i,1,1,e) = w3(i,1,1) &
200  * ( dsdy(i,1,1,e) * us(i,1,1) &
201  + drdy(i,1,1,e) * ur(i,1,1) &
202  + dtdy(i,1,1,e) * ut(i,1,1) )
203  uz(i,1,1,e) = w3(i,1,1) &
204  * ( dtdz(i,1,1,e) * ut(i,1,1) &
205  + drdz(i,1,1,e) * ur(i,1,1) &
206  + dsdz(i,1,1,e) * us(i,1,1) )
207  end do
208  end do
209  end subroutine cpu_opgrad_lx18
210 
211  subroutine cpu_opgrad_lx17(ux, uy, uz, u, dx, dy, dz, &
212  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
213  integer, parameter :: lx = 17
214  integer, intent(in) :: n
215  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
216  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
217  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
218  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
219  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
220  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
221  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
222  real(kind=rp) :: ur(lx,lx,lx)
223  real(kind=rp) :: us(lx,lx,lx)
224  real(kind=rp) :: ut(lx,lx,lx)
225  integer :: e, i, j, k
226 
227  do e = 1, n
228  do j = 1, lx * lx
229  do i = 1, lx
230  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
231  + dx(i,2) * u(2,j,1,e) &
232  + dx(i,3) * u(3,j,1,e) &
233  + dx(i,4) * u(4,j,1,e) &
234  + dx(i,5) * u(5,j,1,e) &
235  + dx(i,6) * u(6,j,1,e) &
236  + dx(i,7) * u(7,j,1,e) &
237  + dx(i,8) * u(8,j,1,e) &
238  + dx(i,9) * u(9,j,1,e) &
239  + dx(i,10) * u(10,j,1,e) &
240  + dx(i,11) * u(11,j,1,e) &
241  + dx(i,12) * u(12,j,1,e) &
242  + dx(i,13) * u(13,j,1,e) &
243  + dx(i,14) * u(14,j,1,e) &
244  + dx(i,15) * u(15,j,1,e) &
245  + dx(i,16) * u(16,j,1,e) &
246  + dx(i,17) * u(17,j,1,e)
247  end do
248  end do
249 
250  do k = 1, lx
251  do j = 1, lx
252  do i = 1, lx
253  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
254  + dy(j,2) * u(i,2,k,e) &
255  + dy(j,3) * u(i,3,k,e) &
256  + dy(j,4) * u(i,4,k,e) &
257  + dy(j,5) * u(i,5,k,e) &
258  + dy(j,6) * u(i,6,k,e) &
259  + dy(j,7) * u(i,7,k,e) &
260  + dy(j,8) * u(i,8,k,e) &
261  + dy(j,9) * u(i,9,k,e) &
262  + dy(j,10) * u(i,10,k,e) &
263  + dy(j,11) * u(i,11,k,e) &
264  + dy(j,12) * u(i,12,k,e) &
265  + dy(j,13) * u(i,13,k,e) &
266  + dy(j,14) * u(i,14,k,e) &
267  + dy(j,15) * u(i,15,k,e) &
268  + dy(j,16) * u(i,16,k,e) &
269  + dy(j,17) * u(i,17,k,e)
270  end do
271  end do
272  end do
273 
274  do k = 1, lx
275  do i = 1, lx*lx
276  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
277  + dz(k,2) * u(i,1,2,e) &
278  + dz(k,3) * u(i,1,3,e) &
279  + dz(k,4) * u(i,1,4,e) &
280  + dz(k,5) * u(i,1,5,e) &
281  + dz(k,6) * u(i,1,6,e) &
282  + dz(k,7) * u(i,1,7,e) &
283  + dz(k,8) * u(i,1,8,e) &
284  + dz(k,9) * u(i,1,9,e) &
285  + dz(k,10) * u(i,1,10,e) &
286  + dz(k,11) * u(i,1,11,e) &
287  + dz(k,12) * u(i,1,12,e) &
288  + dz(k,13) * u(i,1,13,e) &
289  + dz(k,14) * u(i,1,14,e) &
290  + dz(k,15) * u(i,1,15,e) &
291  + dz(k,16) * u(i,1,16,e) &
292  + dz(k,17) * u(i,1,17,e)
293  end do
294  end do
295 
296  do i = 1, lx * lx * lx
297  ux(i,1,1,e) = w3(i,1,1) &
298  * ( drdx(i,1,1,e) * ur(i,1,1) &
299  + dsdx(i,1,1,e) * us(i,1,1) &
300  + dtdx(i,1,1,e) * ut(i,1,1) )
301  uy(i,1,1,e) = w3(i,1,1) &
302  * ( dsdy(i,1,1,e) * us(i,1,1) &
303  + drdy(i,1,1,e) * ur(i,1,1) &
304  + dtdy(i,1,1,e) * ut(i,1,1) )
305  uz(i,1,1,e) = w3(i,1,1) &
306  * ( dtdz(i,1,1,e) * ut(i,1,1) &
307  + drdz(i,1,1,e) * ur(i,1,1) &
308  + dsdz(i,1,1,e) * us(i,1,1) )
309  end do
310  end do
311  end subroutine cpu_opgrad_lx17
312 
313  subroutine cpu_opgrad_lx16(ux, uy, uz, u, dx, dy, dz, &
314  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
315  integer, parameter :: lx = 16
316  integer, intent(in) :: n
317  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
318  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
319  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
320  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
321  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
322  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
323  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
324  real(kind=rp) :: ur(lx,lx,lx)
325  real(kind=rp) :: us(lx,lx,lx)
326  real(kind=rp) :: ut(lx,lx,lx)
327  integer :: e, i, j, k
328 
329  do e = 1, n
330  do j = 1, lx * lx
331  do i = 1, lx
332  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
333  + dx(i,2) * u(2,j,1,e) &
334  + dx(i,3) * u(3,j,1,e) &
335  + dx(i,4) * u(4,j,1,e) &
336  + dx(i,5) * u(5,j,1,e) &
337  + dx(i,6) * u(6,j,1,e) &
338  + dx(i,7) * u(7,j,1,e) &
339  + dx(i,8) * u(8,j,1,e) &
340  + dx(i,9) * u(9,j,1,e) &
341  + dx(i,10) * u(10,j,1,e) &
342  + dx(i,11) * u(11,j,1,e) &
343  + dx(i,12) * u(12,j,1,e) &
344  + dx(i,13) * u(13,j,1,e) &
345  + dx(i,14) * u(14,j,1,e) &
346  + dx(i,15) * u(15,j,1,e) &
347  + dx(i,16) * u(16,j,1,e)
348  end do
349  end do
350 
351  do k = 1, lx
352  do j = 1, lx
353  do i = 1, lx
354  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
355  + dy(j,2) * u(i,2,k,e) &
356  + dy(j,3) * u(i,3,k,e) &
357  + dy(j,4) * u(i,4,k,e) &
358  + dy(j,5) * u(i,5,k,e) &
359  + dy(j,6) * u(i,6,k,e) &
360  + dy(j,7) * u(i,7,k,e) &
361  + dy(j,8) * u(i,8,k,e) &
362  + dy(j,9) * u(i,9,k,e) &
363  + dy(j,10) * u(i,10,k,e) &
364  + dy(j,11) * u(i,11,k,e) &
365  + dy(j,12) * u(i,12,k,e) &
366  + dy(j,13) * u(i,13,k,e) &
367  + dy(j,14) * u(i,14,k,e) &
368  + dy(j,15) * u(i,15,k,e) &
369  + dy(j,16) * u(i,16,k,e)
370  end do
371  end do
372  end do
373 
374  do k = 1, lx
375  do i = 1, lx*lx
376  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
377  + dz(k,2) * u(i,1,2,e) &
378  + dz(k,3) * u(i,1,3,e) &
379  + dz(k,4) * u(i,1,4,e) &
380  + dz(k,5) * u(i,1,5,e) &
381  + dz(k,6) * u(i,1,6,e) &
382  + dz(k,7) * u(i,1,7,e) &
383  + dz(k,8) * u(i,1,8,e) &
384  + dz(k,9) * u(i,1,9,e) &
385  + dz(k,10) * u(i,1,10,e) &
386  + dz(k,11) * u(i,1,11,e) &
387  + dz(k,12) * u(i,1,12,e) &
388  + dz(k,13) * u(i,1,13,e) &
389  + dz(k,14) * u(i,1,14,e) &
390  + dz(k,15) * u(i,1,15,e) &
391  + dz(k,16) * u(i,1,16,e)
392  end do
393  end do
394 
395  do i = 1, lx * lx * lx
396  ux(i,1,1,e) = w3(i,1,1) &
397  * ( drdx(i,1,1,e) * ur(i,1,1) &
398  + dsdx(i,1,1,e) * us(i,1,1) &
399  + dtdx(i,1,1,e) * ut(i,1,1) )
400  uy(i,1,1,e) = w3(i,1,1) &
401  * ( dsdy(i,1,1,e) * us(i,1,1) &
402  + drdy(i,1,1,e) * ur(i,1,1) &
403  + dtdy(i,1,1,e) * ut(i,1,1) )
404  uz(i,1,1,e) = w3(i,1,1) &
405  * ( dtdz(i,1,1,e) * ut(i,1,1) &
406  + drdz(i,1,1,e) * ur(i,1,1) &
407  + dsdz(i,1,1,e) * us(i,1,1) )
408  end do
409  end do
410  end subroutine cpu_opgrad_lx16
411 
412  subroutine cpu_opgrad_lx15(ux, uy, uz, u, dx, dy, dz, &
413  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
414  integer, parameter :: lx = 15
415  integer, intent(in) :: n
416  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
417  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
418  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
419  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
420  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
421  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
422  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
423  real(kind=rp) :: ur(lx,lx,lx)
424  real(kind=rp) :: us(lx,lx,lx)
425  real(kind=rp) :: ut(lx,lx,lx)
426  integer :: e, i, j, k
427 
428  do e = 1, n
429  do j = 1, lx * lx
430  do i = 1, lx
431  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
432  + dx(i,2) * u(2,j,1,e) &
433  + dx(i,3) * u(3,j,1,e) &
434  + dx(i,4) * u(4,j,1,e) &
435  + dx(i,5) * u(5,j,1,e) &
436  + dx(i,6) * u(6,j,1,e) &
437  + dx(i,7) * u(7,j,1,e) &
438  + dx(i,8) * u(8,j,1,e) &
439  + dx(i,9) * u(9,j,1,e) &
440  + dx(i,10) * u(10,j,1,e) &
441  + dx(i,11) * u(11,j,1,e) &
442  + dx(i,12) * u(12,j,1,e) &
443  + dx(i,13) * u(13,j,1,e) &
444  + dx(i,14) * u(14,j,1,e) &
445  + dx(i,15) * u(15,j,1,e)
446  end do
447  end do
448 
449  do k = 1, lx
450  do j = 1, lx
451  do i = 1, lx
452  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
453  + dy(j,2) * u(i,2,k,e) &
454  + dy(j,3) * u(i,3,k,e) &
455  + dy(j,4) * u(i,4,k,e) &
456  + dy(j,5) * u(i,5,k,e) &
457  + dy(j,6) * u(i,6,k,e) &
458  + dy(j,7) * u(i,7,k,e) &
459  + dy(j,8) * u(i,8,k,e) &
460  + dy(j,9) * u(i,9,k,e) &
461  + dy(j,10) * u(i,10,k,e) &
462  + dy(j,11) * u(i,11,k,e) &
463  + dy(j,12) * u(i,12,k,e) &
464  + dy(j,13) * u(i,13,k,e) &
465  + dy(j,14) * u(i,14,k,e) &
466  + dy(j,15) * u(i,15,k,e)
467  end do
468  end do
469  end do
470 
471  do k = 1, lx
472  do i = 1, lx*lx
473  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
474  + dz(k,2) * u(i,1,2,e) &
475  + dz(k,3) * u(i,1,3,e) &
476  + dz(k,4) * u(i,1,4,e) &
477  + dz(k,5) * u(i,1,5,e) &
478  + dz(k,6) * u(i,1,6,e) &
479  + dz(k,7) * u(i,1,7,e) &
480  + dz(k,8) * u(i,1,8,e) &
481  + dz(k,9) * u(i,1,9,e) &
482  + dz(k,10) * u(i,1,10,e) &
483  + dz(k,11) * u(i,1,11,e) &
484  + dz(k,12) * u(i,1,12,e) &
485  + dz(k,13) * u(i,1,13,e) &
486  + dz(k,14) * u(i,1,14,e) &
487  + dz(k,15) * u(i,1,15,e)
488  end do
489  end do
490 
491  do i = 1, lx * lx * lx
492  ux(i,1,1,e) = w3(i,1,1) &
493  * ( drdx(i,1,1,e) * ur(i,1,1) &
494  + dsdx(i,1,1,e) * us(i,1,1) &
495  + dtdx(i,1,1,e) * ut(i,1,1) )
496  uy(i,1,1,e) = w3(i,1,1) &
497  * ( dsdy(i,1,1,e) * us(i,1,1) &
498  + drdy(i,1,1,e) * ur(i,1,1) &
499  + dtdy(i,1,1,e) * ut(i,1,1) )
500  uz(i,1,1,e) = w3(i,1,1) &
501  * ( dtdz(i,1,1,e) * ut(i,1,1) &
502  + drdz(i,1,1,e) * ur(i,1,1) &
503  + dsdz(i,1,1,e) * us(i,1,1) )
504  end do
505  end do
506  end subroutine cpu_opgrad_lx15
507 
508  subroutine cpu_opgrad_lx14(ux, uy, uz, u, dx, dy, dz, &
509  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
510  integer, parameter :: lx = 14
511  integer, intent(in) :: n
512  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
513  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
514  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
515  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
516  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
517  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
518  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
519  real(kind=rp) :: ur(lx,lx,lx)
520  real(kind=rp) :: us(lx,lx,lx)
521  real(kind=rp) :: ut(lx,lx,lx)
522  integer :: e, i, j, k
523 
524  do e = 1, n
525  do j = 1, lx * lx
526  do i = 1, lx
527  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
528  + dx(i,2) * u(2,j,1,e) &
529  + dx(i,3) * u(3,j,1,e) &
530  + dx(i,4) * u(4,j,1,e) &
531  + dx(i,5) * u(5,j,1,e) &
532  + dx(i,6) * u(6,j,1,e) &
533  + dx(i,7) * u(7,j,1,e) &
534  + dx(i,8) * u(8,j,1,e) &
535  + dx(i,9) * u(9,j,1,e) &
536  + dx(i,10) * u(10,j,1,e) &
537  + dx(i,11) * u(11,j,1,e) &
538  + dx(i,12) * u(12,j,1,e) &
539  + dx(i,13) * u(13,j,1,e) &
540  + dx(i,14) * u(14,j,1,e)
541  end do
542  end do
543 
544  do k = 1, lx
545  do j = 1, lx
546  do i = 1, lx
547  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
548  + dy(j,2) * u(i,2,k,e) &
549  + dy(j,3) * u(i,3,k,e) &
550  + dy(j,4) * u(i,4,k,e) &
551  + dy(j,5) * u(i,5,k,e) &
552  + dy(j,6) * u(i,6,k,e) &
553  + dy(j,7) * u(i,7,k,e) &
554  + dy(j,8) * u(i,8,k,e) &
555  + dy(j,9) * u(i,9,k,e) &
556  + dy(j,10) * u(i,10,k,e) &
557  + dy(j,11) * u(i,11,k,e) &
558  + dy(j,12) * u(i,12,k,e) &
559  + dy(j,13) * u(i,13,k,e) &
560  + dy(j,14) * u(i,14,k,e)
561  end do
562  end do
563  end do
564 
565  do k = 1, lx
566  do i = 1, lx*lx
567  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
568  + dz(k,2) * u(i,1,2,e) &
569  + dz(k,3) * u(i,1,3,e) &
570  + dz(k,4) * u(i,1,4,e) &
571  + dz(k,5) * u(i,1,5,e) &
572  + dz(k,6) * u(i,1,6,e) &
573  + dz(k,7) * u(i,1,7,e) &
574  + dz(k,8) * u(i,1,8,e) &
575  + dz(k,9) * u(i,1,9,e) &
576  + dz(k,10) * u(i,1,10,e) &
577  + dz(k,11) * u(i,1,11,e) &
578  + dz(k,12) * u(i,1,12,e) &
579  + dz(k,13) * u(i,1,13,e) &
580  + dz(k,14) * u(i,1,14,e)
581  end do
582  end do
583 
584  do i = 1, lx * lx * lx
585  ux(i,1,1,e) = w3(i,1,1) &
586  * ( drdx(i,1,1,e) * ur(i,1,1) &
587  + dsdx(i,1,1,e) * us(i,1,1) &
588  + dtdx(i,1,1,e) * ut(i,1,1) )
589  uy(i,1,1,e) = w3(i,1,1) &
590  * ( dsdy(i,1,1,e) * us(i,1,1) &
591  + drdy(i,1,1,e) * ur(i,1,1) &
592  + dtdy(i,1,1,e) * ut(i,1,1) )
593  uz(i,1,1,e) = w3(i,1,1) &
594  * ( dtdz(i,1,1,e) * ut(i,1,1) &
595  + drdz(i,1,1,e) * ur(i,1,1) &
596  + dsdz(i,1,1,e) * us(i,1,1) )
597  end do
598  end do
599  end subroutine cpu_opgrad_lx14
600 
601  subroutine cpu_opgrad_lx13(ux, uy, uz, u, dx, dy, dz, &
602  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
603  integer, parameter :: lx = 13
604  integer, intent(in) :: n
605  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
606  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
607  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
608  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
609  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
610  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
611  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
612  real(kind=rp) :: ur(lx,lx,lx)
613  real(kind=rp) :: us(lx,lx,lx)
614  real(kind=rp) :: ut(lx,lx,lx)
615  integer :: e, i, j, k
616 
617  do e = 1, n
618  do j = 1, lx * lx
619  do i = 1, lx
620  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
621  + dx(i,2) * u(2,j,1,e) &
622  + dx(i,3) * u(3,j,1,e) &
623  + dx(i,4) * u(4,j,1,e) &
624  + dx(i,5) * u(5,j,1,e) &
625  + dx(i,6) * u(6,j,1,e) &
626  + dx(i,7) * u(7,j,1,e) &
627  + dx(i,8) * u(8,j,1,e) &
628  + dx(i,9) * u(9,j,1,e) &
629  + dx(i,10) * u(10,j,1,e) &
630  + dx(i,11) * u(11,j,1,e) &
631  + dx(i,12) * u(12,j,1,e) &
632  + dx(i,13) * u(13,j,1,e)
633  end do
634  end do
635 
636  do k = 1, lx
637  do j = 1, lx
638  do i = 1, lx
639  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
640  + dy(j,2) * u(i,2,k,e) &
641  + dy(j,3) * u(i,3,k,e) &
642  + dy(j,4) * u(i,4,k,e) &
643  + dy(j,5) * u(i,5,k,e) &
644  + dy(j,6) * u(i,6,k,e) &
645  + dy(j,7) * u(i,7,k,e) &
646  + dy(j,8) * u(i,8,k,e) &
647  + dy(j,9) * u(i,9,k,e) &
648  + dy(j,10) * u(i,10,k,e) &
649  + dy(j,11) * u(i,11,k,e) &
650  + dy(j,12) * u(i,12,k,e) &
651  + dy(j,13) * u(i,13,k,e)
652  end do
653  end do
654  end do
655 
656  do k = 1, lx
657  do i = 1, lx*lx
658  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
659  + dz(k,2) * u(i,1,2,e) &
660  + dz(k,3) * u(i,1,3,e) &
661  + dz(k,4) * u(i,1,4,e) &
662  + dz(k,5) * u(i,1,5,e) &
663  + dz(k,6) * u(i,1,6,e) &
664  + dz(k,7) * u(i,1,7,e) &
665  + dz(k,8) * u(i,1,8,e) &
666  + dz(k,9) * u(i,1,9,e) &
667  + dz(k,10) * u(i,1,10,e) &
668  + dz(k,11) * u(i,1,11,e) &
669  + dz(k,12) * u(i,1,12,e) &
670  + dz(k,13) * u(i,1,13,e)
671  end do
672  end do
673 
674  do i = 1, lx * lx * lx
675  ux(i,1,1,e) = w3(i,1,1) &
676  * ( drdx(i,1,1,e) * ur(i,1,1) &
677  + dsdx(i,1,1,e) * us(i,1,1) &
678  + dtdx(i,1,1,e) * ut(i,1,1) )
679  uy(i,1,1,e) = w3(i,1,1) &
680  * ( dsdy(i,1,1,e) * us(i,1,1) &
681  + drdy(i,1,1,e) * ur(i,1,1) &
682  + dtdy(i,1,1,e) * ut(i,1,1) )
683  uz(i,1,1,e) = w3(i,1,1) &
684  * ( dtdz(i,1,1,e) * ut(i,1,1) &
685  + drdz(i,1,1,e) * ur(i,1,1) &
686  + dsdz(i,1,1,e) * us(i,1,1) )
687  end do
688  end do
689  end subroutine cpu_opgrad_lx13
690 
691  subroutine cpu_opgrad_lx12(ux, uy, uz, u, dx, dy, dz, &
692  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
693  integer, parameter :: lx = 12
694  integer, intent(in) :: n
695  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
696  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
697  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
698  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
699  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
700  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
701  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
702  real(kind=rp) :: ur(lx,lx,lx)
703  real(kind=rp) :: us(lx,lx,lx)
704  real(kind=rp) :: ut(lx,lx,lx)
705  integer :: e, i, j, k
706 
707  do e = 1, n
708  do j = 1, lx * lx
709  do i = 1, lx
710  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
711  + dx(i,2) * u(2,j,1,e) &
712  + dx(i,3) * u(3,j,1,e) &
713  + dx(i,4) * u(4,j,1,e) &
714  + dx(i,5) * u(5,j,1,e) &
715  + dx(i,6) * u(6,j,1,e) &
716  + dx(i,7) * u(7,j,1,e) &
717  + dx(i,8) * u(8,j,1,e) &
718  + dx(i,9) * u(9,j,1,e) &
719  + dx(i,10) * u(10,j,1,e) &
720  + dx(i,11) * u(11,j,1,e) &
721  + dx(i,12) * u(12,j,1,e)
722  end do
723  end do
724 
725  do k = 1, lx
726  do j = 1, lx
727  do i = 1, lx
728  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
729  + dy(j,2) * u(i,2,k,e) &
730  + dy(j,3) * u(i,3,k,e) &
731  + dy(j,4) * u(i,4,k,e) &
732  + dy(j,5) * u(i,5,k,e) &
733  + dy(j,6) * u(i,6,k,e) &
734  + dy(j,7) * u(i,7,k,e) &
735  + dy(j,8) * u(i,8,k,e) &
736  + dy(j,9) * u(i,9,k,e) &
737  + dy(j,10) * u(i,10,k,e) &
738  + dy(j,11) * u(i,11,k,e) &
739  + dy(j,12) * u(i,12,k,e)
740  end do
741  end do
742  end do
743 
744  do k = 1, lx
745  do i = 1, lx*lx
746  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
747  + dz(k,2) * u(i,1,2,e) &
748  + dz(k,3) * u(i,1,3,e) &
749  + dz(k,4) * u(i,1,4,e) &
750  + dz(k,5) * u(i,1,5,e) &
751  + dz(k,6) * u(i,1,6,e) &
752  + dz(k,7) * u(i,1,7,e) &
753  + dz(k,8) * u(i,1,8,e) &
754  + dz(k,9) * u(i,1,9,e) &
755  + dz(k,10) * u(i,1,10,e) &
756  + dz(k,11) * u(i,1,11,e) &
757  + dz(k,12) * u(i,1,12,e)
758  end do
759  end do
760 
761  do i = 1, lx * lx * lx
762  ux(i,1,1,e) = w3(i,1,1) &
763  * ( drdx(i,1,1,e) * ur(i,1,1) &
764  + dsdx(i,1,1,e) * us(i,1,1) &
765  + dtdx(i,1,1,e) * ut(i,1,1) )
766  uy(i,1,1,e) = w3(i,1,1) &
767  * ( dsdy(i,1,1,e) * us(i,1,1) &
768  + drdy(i,1,1,e) * ur(i,1,1) &
769  + dtdy(i,1,1,e) * ut(i,1,1) )
770  uz(i,1,1,e) = w3(i,1,1) &
771  * ( dtdz(i,1,1,e) * ut(i,1,1) &
772  + drdz(i,1,1,e) * ur(i,1,1) &
773  + dsdz(i,1,1,e) * us(i,1,1) )
774  end do
775  end do
776  end subroutine cpu_opgrad_lx12
777 
778  subroutine cpu_opgrad_lx11(ux, uy, uz, u, dx, dy, dz, &
779  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
780  integer, parameter :: lx = 11
781  integer, intent(in) :: n
782  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
783  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
784  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
785  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
786  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
787  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
788  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
789  real(kind=rp) :: ur(lx,lx,lx)
790  real(kind=rp) :: us(lx,lx,lx)
791  real(kind=rp) :: ut(lx,lx,lx)
792  integer :: e, i, j, k
793 
794  do e = 1, n
795  do j = 1, lx * lx
796  do i = 1, lx
797  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
798  + dx(i,2) * u(2,j,1,e) &
799  + dx(i,3) * u(3,j,1,e) &
800  + dx(i,4) * u(4,j,1,e) &
801  + dx(i,5) * u(5,j,1,e) &
802  + dx(i,6) * u(6,j,1,e) &
803  + dx(i,7) * u(7,j,1,e) &
804  + dx(i,8) * u(8,j,1,e) &
805  + dx(i,9) * u(9,j,1,e) &
806  + dx(i,10) * u(10,j,1,e) &
807  + dx(i,11) * u(11,j,1,e)
808  end do
809  end do
810 
811  do k = 1, lx
812  do j = 1, lx
813  do i = 1, lx
814  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
815  + dy(j,2) * u(i,2,k,e) &
816  + dy(j,3) * u(i,3,k,e) &
817  + dy(j,4) * u(i,4,k,e) &
818  + dy(j,5) * u(i,5,k,e) &
819  + dy(j,6) * u(i,6,k,e) &
820  + dy(j,7) * u(i,7,k,e) &
821  + dy(j,8) * u(i,8,k,e) &
822  + dy(j,9) * u(i,9,k,e) &
823  + dy(j,10) * u(i,10,k,e) &
824  + dy(j,11) * u(i,11,k,e)
825  end do
826  end do
827  end do
828 
829  do k = 1, lx
830  do i = 1, lx*lx
831  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
832  + dz(k,2) * u(i,1,2,e) &
833  + dz(k,3) * u(i,1,3,e) &
834  + dz(k,4) * u(i,1,4,e) &
835  + dz(k,5) * u(i,1,5,e) &
836  + dz(k,6) * u(i,1,6,e) &
837  + dz(k,7) * u(i,1,7,e) &
838  + dz(k,8) * u(i,1,8,e) &
839  + dz(k,9) * u(i,1,9,e) &
840  + dz(k,10) * u(i,1,10,e) &
841  + dz(k,11) * u(i,1,11,e)
842  end do
843  end do
844 
845  do i = 1, lx * lx * lx
846  ux(i,1,1,e) = w3(i,1,1) &
847  * ( drdx(i,1,1,e) * ur(i,1,1) &
848  + dsdx(i,1,1,e) * us(i,1,1) &
849  + dtdx(i,1,1,e) * ut(i,1,1) )
850  uy(i,1,1,e) = w3(i,1,1) &
851  * ( dsdy(i,1,1,e) * us(i,1,1) &
852  + drdy(i,1,1,e) * ur(i,1,1) &
853  + dtdy(i,1,1,e) * ut(i,1,1) )
854  uz(i,1,1,e) = w3(i,1,1) &
855  * ( dtdz(i,1,1,e) * ut(i,1,1) &
856  + drdz(i,1,1,e) * ur(i,1,1) &
857  + dsdz(i,1,1,e) * us(i,1,1) )
858  end do
859  end do
860  end subroutine cpu_opgrad_lx11
861 
862  subroutine cpu_opgrad_lx10(ux, uy, uz, u, dx, dy, dz, &
863  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
864  integer, parameter :: lx = 10
865  integer, intent(in) :: n
866  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
867  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
868  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
869  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
870  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
871  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
872  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
873  real(kind=rp) :: ur(lx,lx,lx)
874  real(kind=rp) :: us(lx,lx,lx)
875  real(kind=rp) :: ut(lx,lx,lx)
876  integer :: e, i, j, k
877 
878  do e = 1, n
879  do j = 1, lx * lx
880  do i = 1, lx
881  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
882  + dx(i,2) * u(2,j,1,e) &
883  + dx(i,3) * u(3,j,1,e) &
884  + dx(i,4) * u(4,j,1,e) &
885  + dx(i,5) * u(5,j,1,e) &
886  + dx(i,6) * u(6,j,1,e) &
887  + dx(i,7) * u(7,j,1,e) &
888  + dx(i,8) * u(8,j,1,e) &
889  + dx(i,9) * u(9,j,1,e) &
890  + dx(i,10) * u(10,j,1,e)
891  end do
892  end do
893 
894  do k = 1, lx
895  do j = 1, lx
896  do i = 1, lx
897  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
898  + dy(j,2) * u(i,2,k,e) &
899  + dy(j,3) * u(i,3,k,e) &
900  + dy(j,4) * u(i,4,k,e) &
901  + dy(j,5) * u(i,5,k,e) &
902  + dy(j,6) * u(i,6,k,e) &
903  + dy(j,7) * u(i,7,k,e) &
904  + dy(j,8) * u(i,8,k,e) &
905  + dy(j,9) * u(i,9,k,e) &
906  + dy(j,10) * u(i,10,k,e)
907  end do
908  end do
909  end do
910 
911  do k = 1, lx
912  do i = 1, lx*lx
913  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
914  + dz(k,2) * u(i,1,2,e) &
915  + dz(k,3) * u(i,1,3,e) &
916  + dz(k,4) * u(i,1,4,e) &
917  + dz(k,5) * u(i,1,5,e) &
918  + dz(k,6) * u(i,1,6,e) &
919  + dz(k,7) * u(i,1,7,e) &
920  + dz(k,8) * u(i,1,8,e) &
921  + dz(k,9) * u(i,1,9,e) &
922  + dz(k,10) * u(i,1,10,e)
923  end do
924  end do
925 
926  do i = 1, lx * lx * lx
927  ux(i,1,1,e) = w3(i,1,1) &
928  * ( drdx(i,1,1,e) * ur(i,1,1) &
929  + dsdx(i,1,1,e) * us(i,1,1) &
930  + dtdx(i,1,1,e) * ut(i,1,1) )
931  uy(i,1,1,e) = w3(i,1,1) &
932  * ( dsdy(i,1,1,e) * us(i,1,1) &
933  + drdy(i,1,1,e) * ur(i,1,1) &
934  + dtdy(i,1,1,e) * ut(i,1,1) )
935  uz(i,1,1,e) = w3(i,1,1) &
936  * ( dtdz(i,1,1,e) * ut(i,1,1) &
937  + drdz(i,1,1,e) * ur(i,1,1) &
938  + dsdz(i,1,1,e) * us(i,1,1) )
939  end do
940  end do
941  end subroutine cpu_opgrad_lx10
942 
943  subroutine cpu_opgrad_lx9(ux, uy, uz, u, dx, dy, dz, &
944  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
945  integer, parameter :: lx = 9
946  integer, intent(in) :: n
947  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
948  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
949  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
950  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
951  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
952  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
953  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
954  real(kind=rp) :: ur(lx,lx,lx)
955  real(kind=rp) :: us(lx,lx,lx)
956  real(kind=rp) :: ut(lx,lx,lx)
957  integer :: e, i, j, k
958 
959  do e = 1, n
960  do j = 1, lx * lx
961  do i = 1, lx
962  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
963  + dx(i,2) * u(2,j,1,e) &
964  + dx(i,3) * u(3,j,1,e) &
965  + dx(i,4) * u(4,j,1,e) &
966  + dx(i,5) * u(5,j,1,e) &
967  + dx(i,6) * u(6,j,1,e) &
968  + dx(i,7) * u(7,j,1,e) &
969  + dx(i,8) * u(8,j,1,e) &
970  + dx(i,9) * u(9,j,1,e)
971  end do
972  end do
973 
974  do k = 1, lx
975  do j = 1, lx
976  do i = 1, lx
977  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
978  + dy(j,2) * u(i,2,k,e) &
979  + dy(j,3) * u(i,3,k,e) &
980  + dy(j,4) * u(i,4,k,e) &
981  + dy(j,5) * u(i,5,k,e) &
982  + dy(j,6) * u(i,6,k,e) &
983  + dy(j,7) * u(i,7,k,e) &
984  + dy(j,8) * u(i,8,k,e) &
985  + dy(j,9) * u(i,9,k,e)
986  end do
987  end do
988  end do
989 
990  do k = 1, lx
991  do i = 1, lx*lx
992  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
993  + dz(k,2) * u(i,1,2,e) &
994  + dz(k,3) * u(i,1,3,e) &
995  + dz(k,4) * u(i,1,4,e) &
996  + dz(k,5) * u(i,1,5,e) &
997  + dz(k,6) * u(i,1,6,e) &
998  + dz(k,7) * u(i,1,7,e) &
999  + dz(k,8) * u(i,1,8,e) &
1000  + dz(k,9) * u(i,1,9,e)
1001  end do
1002  end do
1003 
1004  do i = 1, lx * lx * lx
1005  ux(i,1,1,e) = w3(i,1,1) &
1006  * ( drdx(i,1,1,e) * ur(i,1,1) &
1007  + dsdx(i,1,1,e) * us(i,1,1) &
1008  + dtdx(i,1,1,e) * ut(i,1,1) )
1009  uy(i,1,1,e) = w3(i,1,1) &
1010  * ( dsdy(i,1,1,e) * us(i,1,1) &
1011  + drdy(i,1,1,e) * ur(i,1,1) &
1012  + dtdy(i,1,1,e) * ut(i,1,1) )
1013  uz(i,1,1,e) = w3(i,1,1) &
1014  * ( dtdz(i,1,1,e) * ut(i,1,1) &
1015  + drdz(i,1,1,e) * ur(i,1,1) &
1016  + dsdz(i,1,1,e) * us(i,1,1) )
1017  end do
1018  end do
1019  end subroutine cpu_opgrad_lx9
1020 
1021  subroutine cpu_opgrad_lx8(ux, uy, uz, u, dx, dy, dz, &
1022  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1023  integer, parameter :: lx = 8
1024  integer, intent(in) :: n
1025  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
1026  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
1027  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
1028  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
1029  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
1030  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
1031  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
1032  real(kind=rp) :: ur(lx,lx,lx)
1033  real(kind=rp) :: us(lx,lx,lx)
1034  real(kind=rp) :: ut(lx,lx,lx)
1035  integer :: e, i, j, k
1036 
1037  do e = 1, n
1038  do j = 1, lx * lx
1039  do i = 1, lx
1040  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1041  + dx(i,2) * u(2,j,1,e) &
1042  + dx(i,3) * u(3,j,1,e) &
1043  + dx(i,4) * u(4,j,1,e) &
1044  + dx(i,5) * u(5,j,1,e) &
1045  + dx(i,6) * u(6,j,1,e) &
1046  + dx(i,7) * u(7,j,1,e) &
1047  + dx(i,8) * u(8,j,1,e)
1048  end do
1049  end do
1050 
1051  do k = 1, lx
1052  do j = 1, lx
1053  do i = 1, lx
1054  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1055  + dy(j,2) * u(i,2,k,e) &
1056  + dy(j,3) * u(i,3,k,e) &
1057  + dy(j,4) * u(i,4,k,e) &
1058  + dy(j,5) * u(i,5,k,e) &
1059  + dy(j,6) * u(i,6,k,e) &
1060  + dy(j,7) * u(i,7,k,e) &
1061  + dy(j,8) * u(i,8,k,e)
1062  end do
1063  end do
1064  end do
1065 
1066  do k = 1, lx
1067  do i = 1, lx*lx
1068  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1069  + dz(k,2) * u(i,1,2,e) &
1070  + dz(k,3) * u(i,1,3,e) &
1071  + dz(k,4) * u(i,1,4,e) &
1072  + dz(k,5) * u(i,1,5,e) &
1073  + dz(k,6) * u(i,1,6,e) &
1074  + dz(k,7) * u(i,1,7,e) &
1075  + dz(k,8) * u(i,1,8,e)
1076  end do
1077  end do
1078 
1079  do i = 1, lx * lx * lx
1080  ux(i,1,1,e) = w3(i,1,1) &
1081  * ( drdx(i,1,1,e) * ur(i,1,1) &
1082  + dsdx(i,1,1,e) * us(i,1,1) &
1083  + dtdx(i,1,1,e) * ut(i,1,1) )
1084  uy(i,1,1,e) = w3(i,1,1) &
1085  * ( dsdy(i,1,1,e) * us(i,1,1) &
1086  + drdy(i,1,1,e) * ur(i,1,1) &
1087  + dtdy(i,1,1,e) * ut(i,1,1) )
1088  uz(i,1,1,e) = w3(i,1,1) &
1089  * ( dtdz(i,1,1,e) * ut(i,1,1) &
1090  + drdz(i,1,1,e) * ur(i,1,1) &
1091  + dsdz(i,1,1,e) * us(i,1,1) )
1092  end do
1093  end do
1094  end subroutine cpu_opgrad_lx8
1095 
1096  subroutine cpu_opgrad_lx7(ux, uy, uz, u, dx, dy, dz, &
1097  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1098  integer, parameter :: lx = 7
1099  integer, intent(in) :: n
1100  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
1101  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
1102  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
1103  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
1104  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
1105  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
1106  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
1107  real(kind=rp) :: ur(lx,lx,lx)
1108  real(kind=rp) :: us(lx,lx,lx)
1109  real(kind=rp) :: ut(lx,lx,lx)
1110  integer :: e, i, j, k
1111 
1112  do e = 1, n
1113  do j = 1, lx * lx
1114  do i = 1, lx
1115  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1116  + dx(i,2) * u(2,j,1,e) &
1117  + dx(i,3) * u(3,j,1,e) &
1118  + dx(i,4) * u(4,j,1,e) &
1119  + dx(i,5) * u(5,j,1,e) &
1120  + dx(i,6) * u(6,j,1,e) &
1121  + dx(i,7) * u(7,j,1,e)
1122  end do
1123  end do
1124 
1125  do k = 1, lx
1126  do j = 1, lx
1127  do i = 1, lx
1128  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1129  + dy(j,2) * u(i,2,k,e) &
1130  + dy(j,3) * u(i,3,k,e) &
1131  + dy(j,4) * u(i,4,k,e) &
1132  + dy(j,5) * u(i,5,k,e) &
1133  + dy(j,6) * u(i,6,k,e) &
1134  + dy(j,7) * u(i,7,k,e)
1135  end do
1136  end do
1137  end do
1138 
1139  do k = 1, lx
1140  do i = 1, lx*lx
1141  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1142  + dz(k,2) * u(i,1,2,e) &
1143  + dz(k,3) * u(i,1,3,e) &
1144  + dz(k,4) * u(i,1,4,e) &
1145  + dz(k,5) * u(i,1,5,e) &
1146  + dz(k,6) * u(i,1,6,e) &
1147  + dz(k,7) * u(i,1,7,e)
1148  end do
1149  end do
1150 
1151  do i = 1, lx * lx * lx
1152  ux(i,1,1,e) = w3(i,1,1) &
1153  * ( drdx(i,1,1,e) * ur(i,1,1) &
1154  + dsdx(i,1,1,e) * us(i,1,1) &
1155  + dtdx(i,1,1,e) * ut(i,1,1) )
1156  uy(i,1,1,e) = w3(i,1,1) &
1157  * ( dsdy(i,1,1,e) * us(i,1,1) &
1158  + drdy(i,1,1,e) * ur(i,1,1) &
1159  + dtdy(i,1,1,e) * ut(i,1,1) )
1160  uz(i,1,1,e) = w3(i,1,1) &
1161  * ( dtdz(i,1,1,e) * ut(i,1,1) &
1162  + drdz(i,1,1,e) * ur(i,1,1) &
1163  + dsdz(i,1,1,e) * us(i,1,1) )
1164  end do
1165  end do
1166  end subroutine cpu_opgrad_lx7
1167 
1168  subroutine cpu_opgrad_lx6(ux, uy, uz, u, dx, dy, dz, &
1169  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1170  integer, parameter :: lx = 6
1171  integer, intent(in) :: n
1172  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
1173  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
1174  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
1175  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
1176  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
1177  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
1178  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
1179  real(kind=rp) :: ur(lx,lx,lx)
1180  real(kind=rp) :: us(lx,lx,lx)
1181  real(kind=rp) :: ut(lx,lx,lx)
1182  integer :: e, i, j, k
1183 
1184  do e = 1, n
1185  do j = 1, lx * lx
1186  do i = 1, lx
1187  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1188  + dx(i,2) * u(2,j,1,e) &
1189  + dx(i,3) * u(3,j,1,e) &
1190  + dx(i,4) * u(4,j,1,e) &
1191  + dx(i,5) * u(5,j,1,e) &
1192  + dx(i,6) * u(6,j,1,e)
1193  end do
1194  end do
1195 
1196  do k = 1, lx
1197  do j = 1, lx
1198  do i = 1, lx
1199  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1200  + dy(j,2) * u(i,2,k,e) &
1201  + dy(j,3) * u(i,3,k,e) &
1202  + dy(j,4) * u(i,4,k,e) &
1203  + dy(j,5) * u(i,5,k,e) &
1204  + dy(j,6) * u(i,6,k,e)
1205  end do
1206  end do
1207  end do
1208 
1209  do k = 1, lx
1210  do i = 1, lx*lx
1211  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1212  + dz(k,2) * u(i,1,2,e) &
1213  + dz(k,3) * u(i,1,3,e) &
1214  + dz(k,4) * u(i,1,4,e) &
1215  + dz(k,5) * u(i,1,5,e) &
1216  + dz(k,6) * u(i,1,6,e)
1217  end do
1218  end do
1219 
1220  do i = 1, lx * lx * lx
1221  ux(i,1,1,e) = w3(i,1,1) &
1222  * ( drdx(i,1,1,e) * ur(i,1,1) &
1223  + dsdx(i,1,1,e) * us(i,1,1) &
1224  + dtdx(i,1,1,e) * ut(i,1,1) )
1225  uy(i,1,1,e) = w3(i,1,1) &
1226  * ( dsdy(i,1,1,e) * us(i,1,1) &
1227  + drdy(i,1,1,e) * ur(i,1,1) &
1228  + dtdy(i,1,1,e) * ut(i,1,1) )
1229  uz(i,1,1,e) = w3(i,1,1) &
1230  * ( dtdz(i,1,1,e) * ut(i,1,1) &
1231  + drdz(i,1,1,e) * ur(i,1,1) &
1232  + dsdz(i,1,1,e) * us(i,1,1) )
1233  end do
1234  end do
1235  end subroutine cpu_opgrad_lx6
1236 
1237  subroutine cpu_opgrad_lx5(ux, uy, uz, u, dx, dy, dz, &
1238  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1239  integer, parameter :: lx = 5
1240  integer, intent(in) :: n
1241  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
1242  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
1243  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
1244  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
1245  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
1246  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
1247  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
1248  real(kind=rp) :: ur(lx,lx,lx)
1249  real(kind=rp) :: us(lx,lx,lx)
1250  real(kind=rp) :: ut(lx,lx,lx)
1251  integer :: e, i, j, k
1252 
1253  do e = 1, n
1254  do j = 1, lx * lx
1255  do i = 1, lx
1256  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1257  + dx(i,2) * u(2,j,1,e) &
1258  + dx(i,3) * u(3,j,1,e) &
1259  + dx(i,4) * u(4,j,1,e) &
1260  + dx(i,5) * u(5,j,1,e)
1261  end do
1262  end do
1263 
1264  do k = 1, lx
1265  do j = 1, lx
1266  do i = 1, lx
1267  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1268  + dy(j,2) * u(i,2,k,e) &
1269  + dy(j,3) * u(i,3,k,e) &
1270  + dy(j,4) * u(i,4,k,e) &
1271  + dy(j,5) * u(i,5,k,e)
1272  end do
1273  end do
1274  end do
1275 
1276  do k = 1, lx
1277  do i = 1, lx*lx
1278  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1279  + dz(k,2) * u(i,1,2,e) &
1280  + dz(k,3) * u(i,1,3,e) &
1281  + dz(k,4) * u(i,1,4,e) &
1282  + dz(k,5) * u(i,1,5,e)
1283  end do
1284  end do
1285 
1286  do i = 1, lx * lx * lx
1287  ux(i,1,1,e) = w3(i,1,1) &
1288  * ( drdx(i,1,1,e) * ur(i,1,1) &
1289  + dsdx(i,1,1,e) * us(i,1,1) &
1290  + dtdx(i,1,1,e) * ut(i,1,1) )
1291  uy(i,1,1,e) = w3(i,1,1) &
1292  * ( dsdy(i,1,1,e) * us(i,1,1) &
1293  + drdy(i,1,1,e) * ur(i,1,1) &
1294  + dtdy(i,1,1,e) * ut(i,1,1) )
1295  uz(i,1,1,e) = w3(i,1,1) &
1296  * ( dtdz(i,1,1,e) * ut(i,1,1) &
1297  + drdz(i,1,1,e) * ur(i,1,1) &
1298  + dsdz(i,1,1,e) * us(i,1,1) )
1299  end do
1300  end do
1301  end subroutine cpu_opgrad_lx5
1302 
1303  subroutine cpu_opgrad_lx4(ux, uy, uz, u, dx, dy, dz, &
1304  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1305  integer, parameter :: lx = 4
1306  integer, intent(in) :: n
1307  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
1308  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
1309  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
1310  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
1311  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
1312  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
1313  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
1314  real(kind=rp) :: ur(lx,lx,lx)
1315  real(kind=rp) :: us(lx,lx,lx)
1316  real(kind=rp) :: ut(lx,lx,lx)
1317  integer :: e, i, j, k
1318 
1319  do e = 1, n
1320  do j = 1, lx * lx
1321  do i = 1, lx
1322  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1323  + dx(i,2) * u(2,j,1,e) &
1324  + dx(i,3) * u(3,j,1,e) &
1325  + dx(i,4) * u(4,j,1,e)
1326  end do
1327  end do
1328 
1329  do k = 1, lx
1330  do j = 1, lx
1331  do i = 1, lx
1332  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1333  + dy(j,2) * u(i,2,k,e) &
1334  + dy(j,3) * u(i,3,k,e) &
1335  + dy(j,4) * u(i,4,k,e)
1336  end do
1337  end do
1338  end do
1339 
1340  do k = 1, lx
1341  do i = 1, lx*lx
1342  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1343  + dz(k,2) * u(i,1,2,e) &
1344  + dz(k,3) * u(i,1,3,e) &
1345  + dz(k,4) * u(i,1,4,e)
1346  end do
1347  end do
1348 
1349  do i = 1, lx * lx * lx
1350  ux(i,1,1,e) = w3(i,1,1) &
1351  * ( drdx(i,1,1,e) * ur(i,1,1) &
1352  + dsdx(i,1,1,e) * us(i,1,1) &
1353  + dtdx(i,1,1,e) * ut(i,1,1) )
1354  uy(i,1,1,e) = w3(i,1,1) &
1355  * ( dsdy(i,1,1,e) * us(i,1,1) &
1356  + drdy(i,1,1,e) * ur(i,1,1) &
1357  + dtdy(i,1,1,e) * ut(i,1,1) )
1358  uz(i,1,1,e) = w3(i,1,1) &
1359  * ( dtdz(i,1,1,e) * ut(i,1,1) &
1360  + drdz(i,1,1,e) * ur(i,1,1) &
1361  + dsdz(i,1,1,e) * us(i,1,1) )
1362  end do
1363  end do
1364  end subroutine cpu_opgrad_lx4
1365 
1366  subroutine cpu_opgrad_lx3(ux, uy, uz, u, dx, dy, dz, &
1367  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1368  integer, parameter :: lx = 3
1369  integer, intent(in) :: n
1370  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
1371  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
1372  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
1373  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
1374  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
1375  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
1376  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
1377  real(kind=rp) :: ur(lx,lx,lx)
1378  real(kind=rp) :: us(lx,lx,lx)
1379  real(kind=rp) :: ut(lx,lx,lx)
1380  integer :: e, i, j, k
1381 
1382  do e = 1, n
1383  do j = 1, lx * lx
1384  do i = 1, lx
1385  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1386  + dx(i,2) * u(2,j,1,e) &
1387  + dx(i,3) * u(3,j,1,e)
1388  end do
1389  end do
1390 
1391  do k = 1, lx
1392  do j = 1, lx
1393  do i = 1, lx
1394  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1395  + dy(j,2) * u(i,2,k,e) &
1396  + dy(j,3) * u(i,3,k,e)
1397  end do
1398  end do
1399  end do
1400 
1401  do k = 1, lx
1402  do i = 1, lx*lx
1403  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1404  + dz(k,2) * u(i,1,2,e) &
1405  + dz(k,3) * u(i,1,3,e)
1406  end do
1407  end do
1408 
1409  do i = 1, lx * lx * lx
1410  ux(i,1,1,e) = w3(i,1,1) &
1411  * ( drdx(i,1,1,e) * ur(i,1,1) &
1412  + dsdx(i,1,1,e) * us(i,1,1) &
1413  + dtdx(i,1,1,e) * ut(i,1,1) )
1414  uy(i,1,1,e) = w3(i,1,1) &
1415  * ( dsdy(i,1,1,e) * us(i,1,1) &
1416  + drdy(i,1,1,e) * ur(i,1,1) &
1417  + dtdy(i,1,1,e) * ut(i,1,1) )
1418  uz(i,1,1,e) = w3(i,1,1) &
1419  * ( dtdz(i,1,1,e) * ut(i,1,1) &
1420  + drdz(i,1,1,e) * ur(i,1,1) &
1421  + dsdz(i,1,1,e) * us(i,1,1) )
1422  end do
1423  end do
1424  end subroutine cpu_opgrad_lx3
1425 
1426  subroutine cpu_opgrad_lx2(ux, uy, uz, u, dx, dy, dz, &
1427  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
1428  integer, parameter :: lx = 2
1429  integer, intent(in) :: n
1430  real(kind=rp), dimension(lx,lx,lx,n), intent(inout) :: ux, uy, uz
1431  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: u
1432  real(kind=rp), dimension(lx,lx), intent(in) :: dx, dy, dz
1433  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdx, dsdx, dtdx
1434  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdy, dsdy, dtdy
1435  real(kind=rp), dimension(lx,lx,lx,n), intent(in) :: drdz, dsdz, dtdz
1436  real(kind=rp), dimension(lx,lx,lx), intent(in) :: w3
1437  real(kind=rp) :: ur(lx,lx,lx)
1438  real(kind=rp) :: us(lx,lx,lx)
1439  real(kind=rp) :: ut(lx,lx,lx)
1440  integer :: e, i, j, k
1441 
1442  do e = 1, n
1443  do j = 1, lx * lx
1444  do i = 1, lx
1445  ur(i,j,1) = dx(i,1) * u(1,j,1,e) &
1446  + dx(i,2) * u(2,j,1,e)
1447  end do
1448  end do
1449 
1450  do k = 1, lx
1451  do j = 1, lx
1452  do i = 1, lx
1453  us(i,j,k) = dy(j,1) * u(i,1,k,e) &
1454  + dy(j,2) * u(i,2,k,e)
1455  end do
1456  end do
1457  end do
1458 
1459  do k = 1, lx
1460  do i = 1, lx*lx
1461  ut(i,1,k) = dz(k,1) * u(i,1,1,e) &
1462  + dz(k,2) * u(i,1,2,e)
1463  end do
1464  end do
1465 
1466  do i = 1, lx * lx * lx
1467  ux(i,1,1,e) = w3(i,1,1) &
1468  * ( drdx(i,1,1,e) * ur(i,1,1) &
1469  + dsdx(i,1,1,e) * us(i,1,1) &
1470  + dtdx(i,1,1,e) * ut(i,1,1) )
1471  uy(i,1,1,e) = w3(i,1,1) &
1472  * ( dsdy(i,1,1,e) * us(i,1,1) &
1473  + drdy(i,1,1,e) * ur(i,1,1) &
1474  + dtdy(i,1,1,e) * ut(i,1,1) )
1475  uz(i,1,1,e) = w3(i,1,1) &
1476  * ( dtdz(i,1,1,e) * ut(i,1,1) &
1477  + drdz(i,1,1,e) * ur(i,1,1) &
1478  + dsdz(i,1,1,e) * us(i,1,1) )
1479  end do
1480  end do
1481  end subroutine cpu_opgrad_lx2
1482 
1483 end module cpu_opgrad
Gradient kernels.
Definition: opgrad.f90:34
subroutine cpu_opgrad_lx16(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:315
subroutine cpu_opgrad_lx17(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:213
subroutine cpu_opgrad_lx11(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:780
subroutine cpu_opgrad_lx7(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:1098
subroutine cpu_opgrad_lx8(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:1023
subroutine cpu_opgrad_lx4(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:1305
subroutine cpu_opgrad_lx6(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:1170
subroutine cpu_opgrad_lx3(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:1368
subroutine cpu_opgrad_lx12(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:693
subroutine cpu_opgrad_lx15(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:414
subroutine cpu_opgrad_lx14(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:510
subroutine cpu_opgrad_lx2(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:1428
subroutine cpu_opgrad_lx(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n, lx)
Definition: opgrad.f90:42
subroutine cpu_opgrad_lx13(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:603
subroutine cpu_opgrad_lx18(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:108
subroutine cpu_opgrad_lx5(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:1239
subroutine cpu_opgrad_lx9(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:945
subroutine cpu_opgrad_lx10(ux, uy, uz, u, dx, dy, dz, drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
Definition: opgrad.f90:864
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12