Neko  0.9.0
A portable framework for high-order spectral element flow simulations
set_convect_rst.f90
Go to the documentation of this file.
1 ! Copyright (c) 2024, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 submodule(opr_cpu) cpu_set_convect_rst
35  use num_types, only : rp
36  implicit none
37 
38 contains
39 
40  module subroutine opr_cpu_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
41  type(space_t), intent(inout) :: Xh
42  type(coef_t), intent(inout) :: coef
43  real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
44  intent(inout) :: cr, cs, ct
45  real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
46  intent(in) :: cx, cy, cz
47  associate(drdx => coef%drdx, drdy => coef%drdy, drdz => coef%drdz, &
48  dsdx => coef%dsdx, dsdy => coef%dsdy, dsdz => coef%dsdz, &
49  dtdx => coef%dtdx, dtdy => coef%dtdy, dtdz => coef%dtdz, &
50  nelv => coef%msh%nelv, lx => xh%lx, w3 => xh%w3)
51 
52  select case (lx)
53  case (18)
54  call cpu_set_convect_rst_lx18(cr, cs, ct, cx, cy, cz, &
55  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
56  case (17)
57  call cpu_set_convect_rst_lx17(cr, cs, ct, cx, cy, cz, &
58  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
59  case (16)
60  call cpu_set_convect_rst_lx16(cr, cs, ct, cx, cy, cz, &
61  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
62  case (15)
63  call cpu_set_convect_rst_lx15(cr, cs, ct, cx, cy, cz, drdx, &
64  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
65  case (14)
66  call cpu_set_convect_rst_lx14(cr, cs, ct, cx, cy, cz, drdx, &
67  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
68  case (13)
69  call cpu_set_convect_rst_lx13(cr, cs, ct, cx, cy, cz, drdx, &
70  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
71  case (12)
72  call cpu_set_convect_rst_lx12(cr, cs, ct, cx, cy, cz, drdx, &
73  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
74  case (11)
75  call cpu_set_convect_rst_lx11(cr, cs, ct, cx, cy, cz, drdx, &
76  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
77  case (10)
78  call cpu_set_convect_rst_lx10(cr, cs, ct, cx, cy, cz, drdx, &
79  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
80  case (9)
81  call cpu_set_convect_rst_lx9(cr, cs, ct, cx, cy, cz, drdx, &
82  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
83  case (8)
84  call cpu_set_convect_rst_lx8(cr, cs, ct, cx, cy, cz, drdx, &
85  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
86  case (7)
87  call cpu_set_convect_rst_lx7(cr, cs, ct, cx, cy, cz, drdx, &
88  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
89  case (6)
90  call cpu_set_convect_rst_lx6(cr, cs, ct, cx, cy, cz, drdx, &
91  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
92  case (5)
93  call cpu_set_convect_rst_lx5(cr, cs, ct, cx, cy, cz, drdx, &
94  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
95  case (4)
96  call cpu_set_convect_rst_lx4(cr, cs, ct, cx, cy, cz, drdx, &
97  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
98  case (3)
99  call cpu_set_convect_rst_lx3(cr, cs, ct, cx, cy, cz, drdx, &
100  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
101  case (2)
102  call cpu_set_convect_rst_lx2(cr, cs, ct, cx, cy, cz, drdx, &
103  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
104  case default
105  call cpu_set_convect_rst_lx(cr, cs, ct, cx, cy, cz, drdx, &
106  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv, lx)
107  end select
108  end associate
109 
110  end subroutine opr_cpu_set_convect_rst
111 
112  subroutine cpu_set_convect_rst_lx(cr, cs, ct, cx, cy, cz, &
113  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n, lx)
114  integer, intent(in) :: n, lx
115  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
116  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
117  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
118  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
119  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
120  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
121  integer :: e, i
122 
123  do e = 1, n
124  do i = 1, lx * lx * lx
125  cr(i,1,1,e) = w3(i,1,1) &
126  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
127  + cy(i,1,1,e) * drdy(i,1,1,e) &
128  + cz(i,1,1,e) * drdz(i,1,1,e) )
129  cs(i,1,1,e) = w3(i,1,1) &
130  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
131  + cy(i,1,1,e) * dsdy(i,1,1,e) &
132  + cz(i,1,1,e) * dsdz(i,1,1,e))
133  ct(i,1,1,e) = w3(i,1,1) &
134  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
135  + cy(i,1,1,e) * dtdy(i,1,1,e) &
136  + cz(i,1,1,e) * dtdz(i,1,1,e))
137  end do
138  end do
139 
140  end subroutine cpu_set_convect_rst_lx
141 
142  subroutine cpu_set_convect_rst_lx18(cr, cs, ct, cx, cy, cz, &
143  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
144  integer, parameter :: lx = 18
145  integer, intent(in) :: n
146  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
147  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
148  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
149  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
150  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
151  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
152  integer :: e, i
153 
154  do e = 1, n
155  do i = 1, lx * lx * lx
156  cr(i,1,1,e) = w3(i,1,1) &
157  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
158  + cy(i,1,1,e) * drdy(i,1,1,e) &
159  + cz(i,1,1,e) * drdz(i,1,1,e) )
160  cs(i,1,1,e) = w3(i,1,1) &
161  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
162  + cy(i,1,1,e) * dsdy(i,1,1,e) &
163  + cz(i,1,1,e) * dsdz(i,1,1,e))
164  ct(i,1,1,e) = w3(i,1,1) &
165  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
166  + cy(i,1,1,e) * dtdy(i,1,1,e) &
167  + cz(i,1,1,e) * dtdz(i,1,1,e))
168  end do
169  end do
170 
171  end subroutine cpu_set_convect_rst_lx18
172 
173  subroutine cpu_set_convect_rst_lx17(cr, cs, ct, cx, cy, cz, &
174  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
175  integer, parameter :: lx = 17
176  integer, intent(in) :: n
177  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
178  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
179  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
180  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
181  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
182  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
183  integer :: e, i
184 
185  do e = 1, n
186  do i = 1, lx * lx * lx
187  cr(i,1,1,e) = w3(i,1,1) &
188  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
189  + cy(i,1,1,e) * drdy(i,1,1,e) &
190  + cz(i,1,1,e) * drdz(i,1,1,e) )
191  cs(i,1,1,e) = w3(i,1,1) &
192  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
193  + cy(i,1,1,e) * dsdy(i,1,1,e) &
194  + cz(i,1,1,e) * dsdz(i,1,1,e))
195  ct(i,1,1,e) = w3(i,1,1) &
196  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
197  + cy(i,1,1,e) * dtdy(i,1,1,e) &
198  + cz(i,1,1,e) * dtdz(i,1,1,e))
199  end do
200  end do
201 
202  end subroutine cpu_set_convect_rst_lx17
203 
204  subroutine cpu_set_convect_rst_lx16(cr, cs, ct, cx, cy, cz, &
205  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
206  integer, parameter :: lx = 16
207  integer, intent(in) :: n
208  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
209  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
210  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
211  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
212  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
213  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
214  integer :: e, i
215 
216  do e = 1, n
217  do i = 1, lx * lx * lx
218  cr(i,1,1,e) = w3(i,1,1) &
219  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
220  + cy(i,1,1,e) * drdy(i,1,1,e) &
221  + cz(i,1,1,e) * drdz(i,1,1,e) )
222  cs(i,1,1,e) = w3(i,1,1) &
223  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
224  + cy(i,1,1,e) * dsdy(i,1,1,e) &
225  + cz(i,1,1,e) * dsdz(i,1,1,e))
226  ct(i,1,1,e) = w3(i,1,1) &
227  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
228  + cy(i,1,1,e) * dtdy(i,1,1,e) &
229  + cz(i,1,1,e) * dtdz(i,1,1,e))
230  end do
231  end do
232 
233  end subroutine cpu_set_convect_rst_lx16
234 
235  subroutine cpu_set_convect_rst_lx15(cr, cs, ct, cx, cy, cz, &
236  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
237  integer, parameter :: lx = 15
238  integer, intent(in) :: n
239  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
240  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
241  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
242  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
243  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
244  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
245  integer :: e, i
246 
247  do e = 1, n
248  do i = 1, lx * lx * lx
249  cr(i,1,1,e) = w3(i,1,1) &
250  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
251  + cy(i,1,1,e) * drdy(i,1,1,e) &
252  + cz(i,1,1,e) * drdz(i,1,1,e) )
253  cs(i,1,1,e) = w3(i,1,1) &
254  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
255  + cy(i,1,1,e) * dsdy(i,1,1,e) &
256  + cz(i,1,1,e) * dsdz(i,1,1,e))
257  ct(i,1,1,e) = w3(i,1,1) &
258  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
259  + cy(i,1,1,e) * dtdy(i,1,1,e) &
260  + cz(i,1,1,e) * dtdz(i,1,1,e))
261  end do
262  end do
263 
264  end subroutine cpu_set_convect_rst_lx15
265 
266  subroutine cpu_set_convect_rst_lx14(cr, cs, ct, cx, cy, cz, &
267  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
268  integer, parameter :: lx = 14
269  integer, intent(in) :: n
270  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
271  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
272  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
273  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
274  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
275  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
276  integer :: e, i
277 
278  do e = 1, n
279  do i = 1, lx * lx * lx
280  cr(i,1,1,e) = w3(i,1,1) &
281  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
282  + cy(i,1,1,e) * drdy(i,1,1,e) &
283  + cz(i,1,1,e) * drdz(i,1,1,e) )
284  cs(i,1,1,e) = w3(i,1,1) &
285  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
286  + cy(i,1,1,e) * dsdy(i,1,1,e) &
287  + cz(i,1,1,e) * dsdz(i,1,1,e))
288  ct(i,1,1,e) = w3(i,1,1) &
289  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
290  + cy(i,1,1,e) * dtdy(i,1,1,e) &
291  + cz(i,1,1,e) * dtdz(i,1,1,e))
292  end do
293  end do
294 
295  end subroutine cpu_set_convect_rst_lx14
296 
297  subroutine cpu_set_convect_rst_lx13(cr, cs, ct, cx, cy, cz, &
298  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
299  integer, parameter :: lx = 13
300  integer, intent(in) :: n
301  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
302  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
303  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
304  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
305  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
306  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
307  integer :: e, i
308 
309  do e = 1, n
310  do i = 1, lx * lx * lx
311  cr(i,1,1,e) = w3(i,1,1) &
312  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
313  + cy(i,1,1,e) * drdy(i,1,1,e) &
314  + cz(i,1,1,e) * drdz(i,1,1,e) )
315  cs(i,1,1,e) = w3(i,1,1) &
316  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
317  + cy(i,1,1,e) * dsdy(i,1,1,e) &
318  + cz(i,1,1,e) * dsdz(i,1,1,e))
319  ct(i,1,1,e) = w3(i,1,1) &
320  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
321  + cy(i,1,1,e) * dtdy(i,1,1,e) &
322  + cz(i,1,1,e) * dtdz(i,1,1,e))
323  end do
324  end do
325 
326  end subroutine cpu_set_convect_rst_lx13
327 
328  subroutine cpu_set_convect_rst_lx12(cr, cs, ct, cx, cy, cz, &
329  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
330  integer, parameter :: lx = 12
331  integer, intent(in) :: n
332  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
333  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
334  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
335  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
336  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
337  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
338  integer :: e, i
339 
340  do e = 1, n
341  do i = 1, lx * lx * lx
342  cr(i,1,1,e) = w3(i,1,1) &
343  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
344  + cy(i,1,1,e) * drdy(i,1,1,e) &
345  + cz(i,1,1,e) * drdz(i,1,1,e) )
346  cs(i,1,1,e) = w3(i,1,1) &
347  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
348  + cy(i,1,1,e) * dsdy(i,1,1,e) &
349  + cz(i,1,1,e) * dsdz(i,1,1,e))
350  ct(i,1,1,e) = w3(i,1,1) &
351  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
352  + cy(i,1,1,e) * dtdy(i,1,1,e) &
353  + cz(i,1,1,e) * dtdz(i,1,1,e))
354  end do
355  end do
356 
357  end subroutine cpu_set_convect_rst_lx12
358 
359  subroutine cpu_set_convect_rst_lx11(cr, cs, ct, cx, cy, cz, &
360  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
361  integer, parameter :: lx = 11
362  integer, intent(in) :: n
363  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
364  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
365  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
366  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
367  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
368  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
369  integer :: e, i
370 
371  do e = 1, n
372  do i = 1, lx * lx * lx
373  cr(i,1,1,e) = w3(i,1,1) &
374  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
375  + cy(i,1,1,e) * drdy(i,1,1,e) &
376  + cz(i,1,1,e) * drdz(i,1,1,e) )
377  cs(i,1,1,e) = w3(i,1,1) &
378  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
379  + cy(i,1,1,e) * dsdy(i,1,1,e) &
380  + cz(i,1,1,e) * dsdz(i,1,1,e))
381  ct(i,1,1,e) = w3(i,1,1) &
382  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
383  + cy(i,1,1,e) * dtdy(i,1,1,e) &
384  + cz(i,1,1,e) * dtdz(i,1,1,e))
385  end do
386  end do
387 
388  end subroutine cpu_set_convect_rst_lx11
389 
390  subroutine cpu_set_convect_rst_lx10(cr, cs, ct, cx, cy, cz, &
391  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
392  integer, parameter :: lx = 10
393  integer, intent(in) :: n
394  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
395  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
396  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
397  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
398  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
399  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
400  integer :: e, i
401 
402  do e = 1, n
403  do i = 1, lx * lx * lx
404  cr(i,1,1,e) = w3(i,1,1) &
405  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
406  + cy(i,1,1,e) * drdy(i,1,1,e) &
407  + cz(i,1,1,e) * drdz(i,1,1,e) )
408  cs(i,1,1,e) = w3(i,1,1) &
409  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
410  + cy(i,1,1,e) * dsdy(i,1,1,e) &
411  + cz(i,1,1,e) * dsdz(i,1,1,e))
412  ct(i,1,1,e) = w3(i,1,1) &
413  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
414  + cy(i,1,1,e) * dtdy(i,1,1,e) &
415  + cz(i,1,1,e) * dtdz(i,1,1,e))
416  end do
417  end do
418 
419  end subroutine cpu_set_convect_rst_lx10
420 
421  subroutine cpu_set_convect_rst_lx9(cr, cs, ct, cx, cy, cz, &
422  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
423  integer, parameter :: lx = 9
424  integer, intent(in) :: n
425  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
426  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
427  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
428  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
429  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
430  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
431  integer :: e, i
432 
433  do e = 1, n
434  do i = 1, lx * lx * lx
435  cr(i,1,1,e) = w3(i,1,1) &
436  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
437  + cy(i,1,1,e) * drdy(i,1,1,e) &
438  + cz(i,1,1,e) * drdz(i,1,1,e) )
439  cs(i,1,1,e) = w3(i,1,1) &
440  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
441  + cy(i,1,1,e) * dsdy(i,1,1,e) &
442  + cz(i,1,1,e) * dsdz(i,1,1,e))
443  ct(i,1,1,e) = w3(i,1,1) &
444  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
445  + cy(i,1,1,e) * dtdy(i,1,1,e) &
446  + cz(i,1,1,e) * dtdz(i,1,1,e))
447  end do
448  end do
449 
450  end subroutine cpu_set_convect_rst_lx9
451 
452  subroutine cpu_set_convect_rst_lx8(cr, cs, ct, cx, cy, cz, &
453  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
454  integer, parameter :: lx = 8
455  integer, intent(in) :: n
456  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
457  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
458  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
459  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
460  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
461  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
462  integer :: e, i
463 
464  do e = 1, n
465  do i = 1, lx * lx * lx
466  cr(i,1,1,e) = w3(i,1,1) &
467  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
468  + cy(i,1,1,e) * drdy(i,1,1,e) &
469  + cz(i,1,1,e) * drdz(i,1,1,e) )
470  cs(i,1,1,e) = w3(i,1,1) &
471  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
472  + cy(i,1,1,e) * dsdy(i,1,1,e) &
473  + cz(i,1,1,e) * dsdz(i,1,1,e))
474  ct(i,1,1,e) = w3(i,1,1) &
475  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
476  + cy(i,1,1,e) * dtdy(i,1,1,e) &
477  + cz(i,1,1,e) * dtdz(i,1,1,e))
478  end do
479  end do
480 
481  end subroutine cpu_set_convect_rst_lx8
482 
483  subroutine cpu_set_convect_rst_lx7(cr, cs, ct, cx, cy, cz, &
484  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
485  integer, parameter :: lx = 7
486  integer, intent(in) :: n
487  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
488  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
489  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
490  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
491  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
492  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
493  integer :: e, i
494 
495  do e = 1, n
496  do i = 1, lx * lx * lx
497  cr(i,1,1,e) = w3(i,1,1) &
498  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
499  + cy(i,1,1,e) * drdy(i,1,1,e) &
500  + cz(i,1,1,e) * drdz(i,1,1,e) )
501  cs(i,1,1,e) = w3(i,1,1) &
502  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
503  + cy(i,1,1,e) * dsdy(i,1,1,e) &
504  + cz(i,1,1,e) * dsdz(i,1,1,e))
505  ct(i,1,1,e) = w3(i,1,1) &
506  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
507  + cy(i,1,1,e) * dtdy(i,1,1,e) &
508  + cz(i,1,1,e) * dtdz(i,1,1,e))
509  end do
510  end do
511 
512  end subroutine cpu_set_convect_rst_lx7
513 
514  subroutine cpu_set_convect_rst_lx6(cr, cs, ct, cx, cy, cz, &
515  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
516  integer, parameter :: lx = 6
517  integer, intent(in) :: n
518  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
519  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
520  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
521  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
522  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
523  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
524  integer :: e, i
525 
526  do e = 1, n
527  do i = 1, lx * lx * lx
528  cr(i,1,1,e) = w3(i,1,1) &
529  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
530  + cy(i,1,1,e) * drdy(i,1,1,e) &
531  + cz(i,1,1,e) * drdz(i,1,1,e) )
532  cs(i,1,1,e) = w3(i,1,1) &
533  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
534  + cy(i,1,1,e) * dsdy(i,1,1,e) &
535  + cz(i,1,1,e) * dsdz(i,1,1,e))
536  ct(i,1,1,e) = w3(i,1,1) &
537  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
538  + cy(i,1,1,e) * dtdy(i,1,1,e) &
539  + cz(i,1,1,e) * dtdz(i,1,1,e))
540  end do
541  end do
542 
543  end subroutine cpu_set_convect_rst_lx6
544 
545  subroutine cpu_set_convect_rst_lx5(cr, cs, ct, cx, cy, cz, &
546  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
547  integer, parameter :: lx = 5
548  integer, intent(in) :: n
549  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
550  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
551  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
552  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
553  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
554  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
555  integer :: e, i
556 
557  do e = 1, n
558  do i = 1, lx * lx * lx
559  cr(i,1,1,e) = w3(i,1,1) &
560  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
561  + cy(i,1,1,e) * drdy(i,1,1,e) &
562  + cz(i,1,1,e) * drdz(i,1,1,e) )
563  cs(i,1,1,e) = w3(i,1,1) &
564  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
565  + cy(i,1,1,e) * dsdy(i,1,1,e) &
566  + cz(i,1,1,e) * dsdz(i,1,1,e))
567  ct(i,1,1,e) = w3(i,1,1) &
568  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
569  + cy(i,1,1,e) * dtdy(i,1,1,e) &
570  + cz(i,1,1,e) * dtdz(i,1,1,e))
571  end do
572  end do
573 
574  end subroutine cpu_set_convect_rst_lx5
575 
576  subroutine cpu_set_convect_rst_lx4(cr, cs, ct, cx, cy, cz, &
577  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
578  integer, parameter :: lx = 4
579  integer, intent(in) :: n
580  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
581  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
582  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
583  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
584  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
585  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
586  integer :: e, i
587 
588  do e = 1, n
589  do i = 1, lx * lx * lx
590  cr(i,1,1,e) = w3(i,1,1) &
591  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
592  + cy(i,1,1,e) * drdy(i,1,1,e) &
593  + cz(i,1,1,e) * drdz(i,1,1,e) )
594  cs(i,1,1,e) = w3(i,1,1) &
595  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
596  + cy(i,1,1,e) * dsdy(i,1,1,e) &
597  + cz(i,1,1,e) * dsdz(i,1,1,e))
598  ct(i,1,1,e) = w3(i,1,1) &
599  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
600  + cy(i,1,1,e) * dtdy(i,1,1,e) &
601  + cz(i,1,1,e) * dtdz(i,1,1,e))
602  end do
603  end do
604 
605  end subroutine cpu_set_convect_rst_lx4
606 
607  subroutine cpu_set_convect_rst_lx3(cr, cs, ct, cx, cy, cz, &
608  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
609  integer, parameter :: lx = 3
610  integer, intent(in) :: n
611  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
612  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
613  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
614  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
615  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
616  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
617  integer :: e, i
618 
619  do e = 1, n
620  do i = 1, lx * lx * lx
621  cr(i,1,1,e) = w3(i,1,1) &
622  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
623  + cy(i,1,1,e) * drdy(i,1,1,e) &
624  + cz(i,1,1,e) * drdz(i,1,1,e) )
625  cs(i,1,1,e) = w3(i,1,1) &
626  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
627  + cy(i,1,1,e) * dsdy(i,1,1,e) &
628  + cz(i,1,1,e) * dsdz(i,1,1,e))
629  ct(i,1,1,e) = w3(i,1,1) &
630  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
631  + cy(i,1,1,e) * dtdy(i,1,1,e) &
632  + cz(i,1,1,e) * dtdz(i,1,1,e))
633  end do
634  end do
635 
636  end subroutine cpu_set_convect_rst_lx3
637 
638  subroutine cpu_set_convect_rst_lx2(cr, cs, ct, cx, cy, cz, &
639  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
640  integer, parameter :: lx = 2
641  integer, intent(in) :: n
642  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
643  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
644  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
645  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
646  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
647  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
648  integer :: e, i
649 
650  do e = 1, n
651  do i = 1, lx * lx * lx
652  cr(i,1,1,e) = w3(i,1,1) &
653  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
654  + cy(i,1,1,e) * drdy(i,1,1,e) &
655  + cz(i,1,1,e) * drdz(i,1,1,e) )
656  cs(i,1,1,e) = w3(i,1,1) &
657  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
658  + cy(i,1,1,e) * dsdy(i,1,1,e) &
659  + cz(i,1,1,e) * dsdz(i,1,1,e))
660  ct(i,1,1,e) = w3(i,1,1) &
661  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
662  + cy(i,1,1,e) * dtdy(i,1,1,e) &
663  + cz(i,1,1,e) * dtdz(i,1,1,e))
664  end do
665  end do
666 
667  end subroutine cpu_set_convect_rst_lx2
668 
669 end submodule cpu_set_convect_rst
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Operators CPU backend.
Definition: opr_cpu.f90:34