Neko  0.9.0
A portable framework for high-order spectral element flow simulations
sx_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 !
33 submodule(opr_sx) sx_set_convect_rst
34  use num_types, only : rp
35  implicit none
36 
37 contains
38 
39  module subroutine opr_sx_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
40  type(space_t), intent(inout) :: Xh
41  type(coef_t), intent(inout) :: coef
42  real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
43  intent(inout) :: cr, cs, ct
44  real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
45  intent(in) :: cx, cy, cz
46  associate(drdx => coef%drdx, drdy => coef%drdy, drdz => coef%drdz, &
47  dsdx => coef%dsdx, dsdy => coef%dsdy, dsdz => coef%dsdz, &
48  dtdx => coef%dtdx, dtdy => coef%dtdy, dtdz => coef%dtdz, &
49  nelv => coef%msh%nelv, lx => xh%lx, w3 => xh%w3)
50 
51  select case (lx)
52  case (18)
53  call sx_set_convect_rst_lx18(cr, cs, ct, cx, cy, cz, drdx, &
54  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
55  case (17)
56  call sx_set_convect_rst_lx17(cr, cs, ct, cx, cy, cz, drdx, &
57  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
58  case (16)
59  call sx_set_convect_rst_lx16(cr, cs, ct, cx, cy, cz, drdx, &
60  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
61  case (15)
62  call sx_set_convect_rst_lx15(cr, cs, ct, cx, cy, cz, drdx, &
63  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
64  case (14)
65  call sx_set_convect_rst_lx14(cr, cs, ct, cx, cy, cz, drdx, &
66  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
67  case (13)
68  call sx_set_convect_rst_lx13(cr, cs, ct, cx, cy, cz, drdx, &
69  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
70  case (12)
71  call sx_set_convect_rst_lx12(cr, cs, ct, cx, cy, cz, drdx, &
72  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
73  case (11)
74  call sx_set_convect_rst_lx11(cr, cs, ct, cx, cy, cz, drdx, &
75  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
76  case (10)
77  call sx_set_convect_rst_lx10(cr, cs, ct, cx, cy, cz, drdx, &
78  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
79  case (9)
80  call sx_set_convect_rst_lx9(cr, cs, ct, cx, cy, cz, drdx, &
81  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
82  case (8)
83  call sx_set_convect_rst_lx8(cr, cs, ct, cx, cy, cz, drdx, &
84  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
85  case (7)
86  call sx_set_convect_rst_lx7(cr, cs, ct, cx, cy, cz, drdx, &
87  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
88  case (6)
89  call sx_set_convect_rst_lx6(cr, cs, ct, cx, cy, cz, drdx, &
90  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
91  case (5)
92  call sx_set_convect_rst_lx5(cr, cs, ct, cx, cy, cz, drdx, &
93  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
94  case (4)
95  call sx_set_convect_rst_lx4(cr, cs, ct, cx, cy, cz, drdx, &
96  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
97  case (3)
98  call sx_set_convect_rst_lx3(cr, cs, ct, cx, cy, cz, drdx, &
99  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
100  case (2)
101  call sx_set_convect_rst_lx2(cr, cs, ct, cx, cy, cz, drdx, &
102  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv)
103  case default
104  call sx_set_convect_rst_lx(cr, cs, ct, cx, cy, cz, drdx, &
105  dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, nelv, lx)
106  end select
107  end associate
108  end subroutine opr_sx_set_convect_rst
109 
110  subroutine sx_set_convect_rst_lx(cr, cs, ct, cx, cy, cz, &
111  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n, lx)
112  integer, intent(in) :: n, lx
113  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
114  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
115  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
116  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
117  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
118  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
119  integer :: e, i
120 
121  do i = 1, lx * lx * lx
122  do e = 1, n
123  cr(i,1,1,e) = w3(i,1,1) &
124  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
125  + cy(i,1,1,e) * drdy(i,1,1,e) &
126  + cz(i,1,1,e) * drdz(i,1,1,e) )
127  cs(i,1,1,e) = w3(i,1,1) &
128  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
129  + cy(i,1,1,e) * dsdy(i,1,1,e) &
130  + cz(i,1,1,e) * dsdz(i,1,1,e))
131  ct(i,1,1,e) = w3(i,1,1) &
132  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
133  + cy(i,1,1,e) * dtdy(i,1,1,e) &
134  + cz(i,1,1,e) * dtdz(i,1,1,e))
135  end do
136  end do
137 
138  end subroutine sx_set_convect_rst_lx
139 
140  subroutine sx_set_convect_rst_lx18(cr, cs, ct, cx, cy, cz, &
141  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
142  integer, parameter :: lx = 18
143  integer, intent(in) :: n
144  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
145  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
146  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
147  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
148  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
149  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
150  integer :: e, i
151 
152  do i = 1, lx * lx * lx
153  do e = 1, n
154  cr(i,1,1,e) = w3(i,1,1) &
155  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
156  + cy(i,1,1,e) * drdy(i,1,1,e) &
157  + cz(i,1,1,e) * drdz(i,1,1,e) )
158  cs(i,1,1,e) = w3(i,1,1) &
159  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
160  + cy(i,1,1,e) * dsdy(i,1,1,e) &
161  + cz(i,1,1,e) * dsdz(i,1,1,e))
162  ct(i,1,1,e) = w3(i,1,1) &
163  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
164  + cy(i,1,1,e) * dtdy(i,1,1,e) &
165  + cz(i,1,1,e) * dtdz(i,1,1,e))
166  end do
167  end do
168 
169  end subroutine sx_set_convect_rst_lx18
170 
171  subroutine sx_set_convect_rst_lx17(cr, cs, ct, cx, cy, cz, &
172  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
173  integer, parameter :: lx = 17
174  integer, intent(in) :: n
175  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
176  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
177  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
178  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
179  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
180  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
181  integer :: e, i
182 
183  do i = 1, lx * lx * lx
184  do e = 1, n
185  cr(i,1,1,e) = w3(i,1,1) &
186  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
187  + cy(i,1,1,e) * drdy(i,1,1,e) &
188  + cz(i,1,1,e) * drdz(i,1,1,e) )
189  cs(i,1,1,e) = w3(i,1,1) &
190  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
191  + cy(i,1,1,e) * dsdy(i,1,1,e) &
192  + cz(i,1,1,e) * dsdz(i,1,1,e))
193  ct(i,1,1,e) = w3(i,1,1) &
194  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
195  + cy(i,1,1,e) * dtdy(i,1,1,e) &
196  + cz(i,1,1,e) * dtdz(i,1,1,e))
197  end do
198  end do
199 
200  end subroutine sx_set_convect_rst_lx17
201 
202  subroutine sx_set_convect_rst_lx16(cr, cs, ct, cx, cy, cz, &
203  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
204  integer, parameter :: lx = 16
205  integer, intent(in) :: n
206  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
207  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
208  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
209  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
210  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
211  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
212  integer :: e, i
213 
214  do i = 1, lx * lx * lx
215  do e = 1, n
216  cr(i,1,1,e) = w3(i,1,1) &
217  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
218  + cy(i,1,1,e) * drdy(i,1,1,e) &
219  + cz(i,1,1,e) * drdz(i,1,1,e) )
220  cs(i,1,1,e) = w3(i,1,1) &
221  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
222  + cy(i,1,1,e) * dsdy(i,1,1,e) &
223  + cz(i,1,1,e) * dsdz(i,1,1,e))
224  ct(i,1,1,e) = w3(i,1,1) &
225  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
226  + cy(i,1,1,e) * dtdy(i,1,1,e) &
227  + cz(i,1,1,e) * dtdz(i,1,1,e))
228  end do
229  end do
230 
231  end subroutine sx_set_convect_rst_lx16
232 
233  subroutine sx_set_convect_rst_lx15(cr, cs, ct, cx, cy, cz, &
234  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
235  integer, parameter :: lx = 15
236  integer, intent(in) :: n
237  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
238  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
239  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
240  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
241  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
242  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
243  integer :: e, i
244 
245  do i = 1, lx * lx * lx
246  do e = 1, n
247  cr(i,1,1,e) = w3(i,1,1) &
248  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
249  + cy(i,1,1,e) * drdy(i,1,1,e) &
250  + cz(i,1,1,e) * drdz(i,1,1,e) )
251  cs(i,1,1,e) = w3(i,1,1) &
252  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
253  + cy(i,1,1,e) * dsdy(i,1,1,e) &
254  + cz(i,1,1,e) * dsdz(i,1,1,e))
255  ct(i,1,1,e) = w3(i,1,1) &
256  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
257  + cy(i,1,1,e) * dtdy(i,1,1,e) &
258  + cz(i,1,1,e) * dtdz(i,1,1,e))
259  end do
260  end do
261 
262  end subroutine sx_set_convect_rst_lx15
263 
264  subroutine sx_set_convect_rst_lx14(cr, cs, ct, cx, cy, cz, &
265  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
266  integer, parameter :: lx = 14
267  integer, intent(in) :: n
268  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
269  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
270  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
271  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
272  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
273  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
274  integer :: e, i
275 
276  do i = 1, lx * lx * lx
277  do e = 1, n
278  cr(i,1,1,e) = w3(i,1,1) &
279  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
280  + cy(i,1,1,e) * drdy(i,1,1,e) &
281  + cz(i,1,1,e) * drdz(i,1,1,e) )
282  cs(i,1,1,e) = w3(i,1,1) &
283  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
284  + cy(i,1,1,e) * dsdy(i,1,1,e) &
285  + cz(i,1,1,e) * dsdz(i,1,1,e))
286  ct(i,1,1,e) = w3(i,1,1) &
287  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
288  + cy(i,1,1,e) * dtdy(i,1,1,e) &
289  + cz(i,1,1,e) * dtdz(i,1,1,e))
290  end do
291  end do
292 
293  end subroutine sx_set_convect_rst_lx14
294 
295  subroutine sx_set_convect_rst_lx13(cr, cs, ct, cx, cy, cz, &
296  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
297  integer, parameter :: lx = 13
298  integer, intent(in) :: n
299  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
300  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
301  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
302  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
303  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
304  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
305  integer :: e, i
306 
307  do i = 1, lx * lx * lx
308  do e = 1, n
309  cr(i,1,1,e) = w3(i,1,1) &
310  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
311  + cy(i,1,1,e) * drdy(i,1,1,e) &
312  + cz(i,1,1,e) * drdz(i,1,1,e) )
313  cs(i,1,1,e) = w3(i,1,1) &
314  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
315  + cy(i,1,1,e) * dsdy(i,1,1,e) &
316  + cz(i,1,1,e) * dsdz(i,1,1,e))
317  ct(i,1,1,e) = w3(i,1,1) &
318  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
319  + cy(i,1,1,e) * dtdy(i,1,1,e) &
320  + cz(i,1,1,e) * dtdz(i,1,1,e))
321  end do
322  end do
323 
324  end subroutine sx_set_convect_rst_lx13
325 
326  subroutine sx_set_convect_rst_lx12(cr, cs, ct, cx, cy, cz, &
327  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
328  integer, parameter :: lx = 12
329  integer, intent(in) :: n
330  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
331  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
332  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
333  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
334  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
335  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
336  integer :: e, i
337 
338  do i = 1, lx * lx * lx
339  do e = 1, n
340  cr(i,1,1,e) = w3(i,1,1) &
341  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
342  + cy(i,1,1,e) * drdy(i,1,1,e) &
343  + cz(i,1,1,e) * drdz(i,1,1,e) )
344  cs(i,1,1,e) = w3(i,1,1) &
345  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
346  + cy(i,1,1,e) * dsdy(i,1,1,e) &
347  + cz(i,1,1,e) * dsdz(i,1,1,e))
348  ct(i,1,1,e) = w3(i,1,1) &
349  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
350  + cy(i,1,1,e) * dtdy(i,1,1,e) &
351  + cz(i,1,1,e) * dtdz(i,1,1,e))
352  end do
353  end do
354 
355  end subroutine sx_set_convect_rst_lx12
356 
357  subroutine sx_set_convect_rst_lx11(cr, cs, ct, cx, cy, cz, &
358  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
359  integer, parameter :: lx = 11
360  integer, intent(in) :: n
361  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
362  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
363  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
364  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
365  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
366  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
367  integer :: e, i
368 
369  do i = 1, lx * lx * lx
370  do e = 1, n
371  cr(i,1,1,e) = w3(i,1,1) &
372  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
373  + cy(i,1,1,e) * drdy(i,1,1,e) &
374  + cz(i,1,1,e) * drdz(i,1,1,e) )
375  cs(i,1,1,e) = w3(i,1,1) &
376  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
377  + cy(i,1,1,e) * dsdy(i,1,1,e) &
378  + cz(i,1,1,e) * dsdz(i,1,1,e))
379  ct(i,1,1,e) = w3(i,1,1) &
380  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
381  + cy(i,1,1,e) * dtdy(i,1,1,e) &
382  + cz(i,1,1,e) * dtdz(i,1,1,e))
383  end do
384  end do
385 
386  end subroutine sx_set_convect_rst_lx11
387 
388  subroutine sx_set_convect_rst_lx10(cr, cs, ct, cx, cy, cz, &
389  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
390  integer, parameter :: lx = 10
391  integer, intent(in) :: n
392  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
393  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
394  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
395  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
396  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
397  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
398  integer :: e, i
399 
400  do i = 1, lx * lx * lx
401  do e = 1, n
402  cr(i,1,1,e) = w3(i,1,1) &
403  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
404  + cy(i,1,1,e) * drdy(i,1,1,e) &
405  + cz(i,1,1,e) * drdz(i,1,1,e) )
406  cs(i,1,1,e) = w3(i,1,1) &
407  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
408  + cy(i,1,1,e) * dsdy(i,1,1,e) &
409  + cz(i,1,1,e) * dsdz(i,1,1,e))
410  ct(i,1,1,e) = w3(i,1,1) &
411  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
412  + cy(i,1,1,e) * dtdy(i,1,1,e) &
413  + cz(i,1,1,e) * dtdz(i,1,1,e))
414  end do
415  end do
416 
417  end subroutine sx_set_convect_rst_lx10
418 
419  subroutine sx_set_convect_rst_lx9(cr, cs, ct, cx, cy, cz, &
420  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
421  integer, parameter :: lx = 9
422  integer, intent(in) :: n
423  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
424  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
425  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
426  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
427  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
428  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
429  integer :: e, i
430 
431  do i = 1, lx * lx * lx
432  do e = 1, n
433  cr(i,1,1,e) = w3(i,1,1) &
434  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
435  + cy(i,1,1,e) * drdy(i,1,1,e) &
436  + cz(i,1,1,e) * drdz(i,1,1,e) )
437  cs(i,1,1,e) = w3(i,1,1) &
438  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
439  + cy(i,1,1,e) * dsdy(i,1,1,e) &
440  + cz(i,1,1,e) * dsdz(i,1,1,e))
441  ct(i,1,1,e) = w3(i,1,1) &
442  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
443  + cy(i,1,1,e) * dtdy(i,1,1,e) &
444  + cz(i,1,1,e) * dtdz(i,1,1,e))
445  end do
446  end do
447 
448  end subroutine sx_set_convect_rst_lx9
449 
450  subroutine sx_set_convect_rst_lx8(cr, cs, ct, cx, cy, cz, &
451  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
452  integer, parameter :: lx = 8
453  integer, intent(in) :: n
454  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
455  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
456  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
457  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
458  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
459  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
460  integer :: e, i
461 
462  do i = 1, lx * lx * lx
463  do e = 1, n
464  cr(i,1,1,e) = w3(i,1,1) &
465  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
466  + cy(i,1,1,e) * drdy(i,1,1,e) &
467  + cz(i,1,1,e) * drdz(i,1,1,e) )
468  cs(i,1,1,e) = w3(i,1,1) &
469  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
470  + cy(i,1,1,e) * dsdy(i,1,1,e) &
471  + cz(i,1,1,e) * dsdz(i,1,1,e))
472  ct(i,1,1,e) = w3(i,1,1) &
473  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
474  + cy(i,1,1,e) * dtdy(i,1,1,e) &
475  + cz(i,1,1,e) * dtdz(i,1,1,e))
476  end do
477  end do
478 
479  end subroutine sx_set_convect_rst_lx8
480 
481  subroutine sx_set_convect_rst_lx7(cr, cs, ct, cx, cy, cz, &
482  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
483  integer, parameter :: lx = 7
484  integer, intent(in) :: n
485  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
486  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
487  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
488  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
489  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
490  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
491  integer :: e, i
492 
493  do i = 1, lx * lx * lx
494  do e = 1, n
495  cr(i,1,1,e) = w3(i,1,1) &
496  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
497  + cy(i,1,1,e) * drdy(i,1,1,e) &
498  + cz(i,1,1,e) * drdz(i,1,1,e) )
499  cs(i,1,1,e) = w3(i,1,1) &
500  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
501  + cy(i,1,1,e) * dsdy(i,1,1,e) &
502  + cz(i,1,1,e) * dsdz(i,1,1,e))
503  ct(i,1,1,e) = w3(i,1,1) &
504  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
505  + cy(i,1,1,e) * dtdy(i,1,1,e) &
506  + cz(i,1,1,e) * dtdz(i,1,1,e))
507  end do
508  end do
509 
510  end subroutine sx_set_convect_rst_lx7
511 
512  subroutine sx_set_convect_rst_lx6(cr, cs, ct, cx, cy, cz, &
513  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
514  integer, parameter :: lx = 6
515  integer, intent(in) :: n
516  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
517  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
518  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
519  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
520  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
521  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
522  integer :: e, i
523 
524  do i = 1, lx * lx * lx
525  do e = 1, n
526  cr(i,1,1,e) = w3(i,1,1) &
527  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
528  + cy(i,1,1,e) * drdy(i,1,1,e) &
529  + cz(i,1,1,e) * drdz(i,1,1,e) )
530  cs(i,1,1,e) = w3(i,1,1) &
531  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
532  + cy(i,1,1,e) * dsdy(i,1,1,e) &
533  + cz(i,1,1,e) * dsdz(i,1,1,e))
534  ct(i,1,1,e) = w3(i,1,1) &
535  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
536  + cy(i,1,1,e) * dtdy(i,1,1,e) &
537  + cz(i,1,1,e) * dtdz(i,1,1,e))
538  end do
539  end do
540 
541 end subroutine sx_set_convect_rst_lx6
542 
543 subroutine sx_set_convect_rst_lx5(cr, cs, ct, cx, cy, cz, &
544  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
545  integer, parameter :: lx = 5
546  integer, intent(in) :: n
547  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
548  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
549  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
550  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
551  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
552  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
553  integer :: e, i
554 
555  do i = 1, lx * lx * lx
556  do e = 1, n
557  cr(i,1,1,e) = w3(i,1,1) &
558  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
559  + cy(i,1,1,e) * drdy(i,1,1,e) &
560  + cz(i,1,1,e) * drdz(i,1,1,e) )
561  cs(i,1,1,e) = w3(i,1,1) &
562  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
563  + cy(i,1,1,e) * dsdy(i,1,1,e) &
564  + cz(i,1,1,e) * dsdz(i,1,1,e))
565  ct(i,1,1,e) = w3(i,1,1) &
566  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
567  + cy(i,1,1,e) * dtdy(i,1,1,e) &
568  + cz(i,1,1,e) * dtdz(i,1,1,e))
569  end do
570  end do
571 
572  end subroutine sx_set_convect_rst_lx5
573 
574  subroutine sx_set_convect_rst_lx4(cr, cs, ct, cx, cy, cz, &
575  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
576  integer, parameter :: lx = 4
577  integer, intent(in) :: n
578  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
579  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
580  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
581  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
582  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
583  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
584  integer :: e, i
585 
586  do i = 1, lx * lx * lx
587  do e = 1, n
588  cr(i,1,1,e) = w3(i,1,1) &
589  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
590  + cy(i,1,1,e) * drdy(i,1,1,e) &
591  + cz(i,1,1,e) * drdz(i,1,1,e) )
592  cs(i,1,1,e) = w3(i,1,1) &
593  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
594  + cy(i,1,1,e) * dsdy(i,1,1,e) &
595  + cz(i,1,1,e) * dsdz(i,1,1,e))
596  ct(i,1,1,e) = w3(i,1,1) &
597  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
598  + cy(i,1,1,e) * dtdy(i,1,1,e) &
599  + cz(i,1,1,e) * dtdz(i,1,1,e))
600  end do
601  end do
602 
603  end subroutine sx_set_convect_rst_lx4
604 
605  subroutine sx_set_convect_rst_lx3(cr, cs, ct, cx, cy, cz, &
606  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
607  integer, parameter :: lx = 3
608  integer, intent(in) :: n
609  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
610  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
611  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
612  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
613  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
614  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
615  integer :: e, i
616 
617  do i = 1, lx * lx * lx
618  do e = 1, n
619  cr(i,1,1,e) = w3(i,1,1) &
620  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
621  + cy(i,1,1,e) * drdy(i,1,1,e) &
622  + cz(i,1,1,e) * drdz(i,1,1,e) )
623  cs(i,1,1,e) = w3(i,1,1) &
624  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
625  + cy(i,1,1,e) * dsdy(i,1,1,e) &
626  + cz(i,1,1,e) * dsdz(i,1,1,e))
627  ct(i,1,1,e) = w3(i,1,1) &
628  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
629  + cy(i,1,1,e) * dtdy(i,1,1,e) &
630  + cz(i,1,1,e) * dtdz(i,1,1,e))
631  end do
632  end do
633 
634  end subroutine sx_set_convect_rst_lx3
635 
636  subroutine sx_set_convect_rst_lx2(cr, cs, ct, cx, cy, cz, &
637  drdx, dsdx, dtdx, drdy, dsdy, dtdy, drdz, dsdz, dtdz, w3, n)
638  integer, parameter :: lx = 2
639  integer, intent(in) :: n
640  real(kind=rp), dimension(lx, lx, lx, n), intent(inout) :: cr, cs, ct
641  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: cx, cy, cz
642  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdx, dsdx, dtdx
643  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdy, dsdy, dtdy
644  real(kind=rp), dimension(lx, lx, lx, n), intent(in) :: drdz, dsdz, dtdz
645  real(kind=rp), dimension(lx, lx, lx), intent(in) :: w3
646  integer :: e, i
647 
648  do i = 1, lx * lx * lx
649  do e = 1, n
650  cr(i,1,1,e) = w3(i,1,1) &
651  * ( cx(i,1,1,e) * drdx(i,1,1,e) &
652  + cy(i,1,1,e) * drdy(i,1,1,e) &
653  + cz(i,1,1,e) * drdz(i,1,1,e) )
654  cs(i,1,1,e) = w3(i,1,1) &
655  * ( cx(i,1,1,e) * dsdx(i,1,1,e) &
656  + cy(i,1,1,e) * dsdy(i,1,1,e) &
657  + cz(i,1,1,e) * dsdz(i,1,1,e))
658  ct(i,1,1,e) = w3(i,1,1) &
659  * ( cx(i,1,1,e) * dtdx(i,1,1,e) &
660  + cy(i,1,1,e) * dtdy(i,1,1,e) &
661  + cz(i,1,1,e) * dtdz(i,1,1,e))
662  end do
663  end do
664 
665  end subroutine sx_set_convect_rst_lx2
666 
667 end submodule sx_set_convect_rst
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Operators SX-Aurora backend.
Definition: opr_sx.f90:2