Neko  0.8.1
A portable framework for high-order spectral element flow simulations
opr_sx.f90
Go to the documentation of this file.
1 
2 module opr_sx
3  use sx_dudxyz
4  use sx_opgrad
5  use sx_conv1
6  use sx_cdtp
7  use sx_cfl
8  use sx_lambda2
10  use num_types, only : rp
11  use space, only : space_t
12  use coefs, only : coef_t
13  use math
14  use field, only : field_t
15  use mathops
16  implicit none
17  private
18 
21 
22 contains
23 
24  subroutine opr_sx_dudxyz(du, u, dr, ds, dt, coef)
25  type(coef_t), intent(in), target :: coef
26  real(kind=rp), dimension(coef%Xh%lx,coef%Xh%ly,coef%Xh%lz,coef%msh%nelv), intent(inout) :: du
27  real(kind=rp), dimension(coef%Xh%lx,coef%Xh%ly,coef%Xh%lz,coef%msh%nelv), intent(in) :: u, dr, ds, dt
28 
29  associate(xh => coef%Xh, msh => coef%msh, dof => coef%dof)
30  select case(coef%Xh%lx)
31  case(14)
32  call sx_dudxyz_lx14(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
33  coef%jacinv, msh%nelv, dof%size())
34  case(13)
35  call sx_dudxyz_lx13(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
36  coef%jacinv, msh%nelv, dof%size())
37  case(12)
38  call sx_dudxyz_lx12(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
39  coef%jacinv, msh%nelv, dof%size())
40  case(11)
41  call sx_dudxyz_lx11(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
42  coef%jacinv, msh%nelv, dof%size())
43  case(10)
44  call sx_dudxyz_lx10(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
45  coef%jacinv, msh%nelv, dof%size())
46  case(9)
47  call sx_dudxyz_lx9(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
48  coef%jacinv, msh%nelv, dof%size())
49  case(8)
50  call sx_dudxyz_lx8(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
51  coef%jacinv, msh%nelv, dof%size())
52  case(7)
53  call sx_dudxyz_lx7(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
54  coef%jacinv, msh%nelv, dof%size())
55  case(6)
56  call sx_dudxyz_lx6(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
57  coef%jacinv, msh%nelv, dof%size())
58  case(5)
59  call sx_dudxyz_lx5(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
60  coef%jacinv, msh%nelv, dof%size())
61  case(4)
62  call sx_dudxyz_lx4(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
63  coef%jacinv, msh%nelv, dof%size())
64  case(3)
65  call sx_dudxyz_lx3(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
66  coef%jacinv, msh%nelv, dof%size())
67  case(2)
68  call sx_dudxyz_lx2(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
69  coef%jacinv, msh%nelv, dof%size())
70  case default
71  call sx_dudxyz_lx(du, u, dr, ds, dt, xh%dx, xh%dy, xh%dz, &
72  coef%jacinv, msh%nelv, dof%size(), xh%lx)
73  end select
74  end associate
75 
76  end subroutine opr_sx_dudxyz
77 
78  subroutine opr_sx_opgrad(ux,uy,uz,u,coef)
79  type(coef_t), intent(in) :: coef
80  real(kind=rp), dimension(coef%Xh%lxyz,coef%msh%nelv), intent(inout) :: ux
81  real(kind=rp), dimension(coef%Xh%lxyz,coef%msh%nelv), intent(inout) :: uy
82  real(kind=rp), dimension(coef%Xh%lxyz,coef%msh%nelv), intent(inout) :: uz
83  real(kind=rp), dimension(coef%Xh%lxyz,coef%msh%nelv), intent(in) :: u
84 
85  associate(xh => coef%Xh, msh => coef%msh)
86  select case(xh%lx)
87  case(18)
88  call sx_opgrad_lx18(ux, uy, uz, u, &
89  xh%dx, xh%dy, xh%dz, &
90  coef%drdx, coef%dsdx, coef%dtdx, &
91  coef%drdy, coef%dsdy, coef%dtdy, &
92  coef%drdz, coef%dsdz, coef%dtdz, &
93  xh%w3, msh%nelv)
94  case(17)
95  call sx_opgrad_lx17(ux, uy, uz, u, &
96  xh%dx, xh%dy, xh%dz, &
97  coef%drdx, coef%dsdx, coef%dtdx, &
98  coef%drdy, coef%dsdy, coef%dtdy, &
99  coef%drdz, coef%dsdz, coef%dtdz, &
100  xh%w3, msh%nelv)
101  case(16)
102  call sx_opgrad_lx16(ux, uy, uz, u, &
103  xh%dx, xh%dy, xh%dz, &
104  coef%drdx, coef%dsdx, coef%dtdx, &
105  coef%drdy, coef%dsdy, coef%dtdy, &
106  coef%drdz, coef%dsdz, coef%dtdz, &
107  xh%w3, msh%nelv)
108  case(15)
109  call sx_opgrad_lx15(ux, uy, uz, u, &
110  xh%dx, xh%dy, xh%dz, &
111  coef%drdx, coef%dsdx, coef%dtdx, &
112  coef%drdy, coef%dsdy, coef%dtdy, &
113  coef%drdz, coef%dsdz, coef%dtdz, &
114  xh%w3, msh%nelv)
115  case(14)
116  call sx_opgrad_lx14(ux, uy, uz, u, &
117  xh%dx, xh%dy, xh%dz, &
118  coef%drdx, coef%dsdx, coef%dtdx, &
119  coef%drdy, coef%dsdy, coef%dtdy, &
120  coef%drdz, coef%dsdz, coef%dtdz, &
121  xh%w3, msh%nelv)
122  case(13)
123  call sx_opgrad_lx13(ux, uy, uz, u, &
124  xh%dx, xh%dy, xh%dz, &
125  coef%drdx, coef%dsdx, coef%dtdx, &
126  coef%drdy, coef%dsdy, coef%dtdy, &
127  coef%drdz, coef%dsdz, coef%dtdz, &
128  xh%w3, msh%nelv)
129  case(12)
130  call sx_opgrad_lx12(ux, uy, uz, u, &
131  xh%dx, xh%dy, xh%dz, &
132  coef%drdx, coef%dsdx, coef%dtdx, &
133  coef%drdy, coef%dsdy, coef%dtdy, &
134  coef%drdz, coef%dsdz, coef%dtdz, &
135  xh%w3, msh%nelv)
136  case(11)
137  call sx_opgrad_lx11(ux, uy, uz, u, &
138  xh%dx, xh%dy, xh%dz, &
139  coef%drdx, coef%dsdx, coef%dtdx, &
140  coef%drdy, coef%dsdy, coef%dtdy, &
141  coef%drdz, coef%dsdz, coef%dtdz, &
142  xh%w3, msh%nelv)
143  case(10)
144  call sx_opgrad_lx10(ux, uy, uz, u, &
145  xh%dx, xh%dy, xh%dz, &
146  coef%drdx, coef%dsdx, coef%dtdx, &
147  coef%drdy, coef%dsdy, coef%dtdy, &
148  coef%drdz, coef%dsdz, coef%dtdz, &
149  xh%w3, msh%nelv)
150  case(9)
151  call sx_opgrad_lx9(ux, uy, uz, u, &
152  xh%dx, xh%dy, xh%dz, &
153  coef%drdx, coef%dsdx, coef%dtdx, &
154  coef%drdy, coef%dsdy, coef%dtdy, &
155  coef%drdz, coef%dsdz, coef%dtdz, &
156  xh%w3, msh%nelv)
157  case(8)
158  call sx_opgrad_lx8(ux, uy, uz, u, &
159  xh%dx, xh%dy, xh%dz, &
160  coef%drdx, coef%dsdx, coef%dtdx, &
161  coef%drdy, coef%dsdy, coef%dtdy, &
162  coef%drdz, coef%dsdz, coef%dtdz, &
163  xh%w3, msh%nelv)
164  case(7)
165  call sx_opgrad_lx7(ux, uy, uz, u, &
166  xh%dx, xh%dy, xh%dz, &
167  coef%drdx, coef%dsdx, coef%dtdx, &
168  coef%drdy, coef%dsdy, coef%dtdy, &
169  coef%drdz, coef%dsdz, coef%dtdz, &
170  xh%w3, msh%nelv)
171  case(6)
172  call sx_opgrad_lx6(ux, uy, uz, u, &
173  xh%dx, xh%dy, xh%dz, &
174  coef%drdx, coef%dsdx, coef%dtdx, &
175  coef%drdy, coef%dsdy, coef%dtdy, &
176  coef%drdz, coef%dsdz, coef%dtdz, &
177  xh%w3, msh%nelv)
178  case(5)
179  call sx_opgrad_lx5(ux, uy, uz, u, &
180  xh%dx, xh%dy, xh%dz, &
181  coef%drdx, coef%dsdx, coef%dtdx, &
182  coef%drdy, coef%dsdy, coef%dtdy, &
183  coef%drdz, coef%dsdz, coef%dtdz, &
184  xh%w3, msh%nelv)
185  case(4)
186  call sx_opgrad_lx4(ux, uy, uz, u, &
187  xh%dx, xh%dy, xh%dz, &
188  coef%drdx, coef%dsdx, coef%dtdx, &
189  coef%drdy, coef%dsdy, coef%dtdy, &
190  coef%drdz, coef%dsdz, coef%dtdz, &
191  xh%w3, msh%nelv)
192  case(3)
193  call sx_opgrad_lx3(ux, uy, uz, u, &
194  xh%dx, xh%dy, xh%dz, &
195  coef%drdx, coef%dsdx, coef%dtdx, &
196  coef%drdy, coef%dsdy, coef%dtdy, &
197  coef%drdz, coef%dsdz, coef%dtdz, &
198  xh%w3, msh%nelv)
199  case(2)
200  call sx_opgrad_lx2(ux, uy, uz, u, &
201  xh%dx, xh%dy, xh%dz, &
202  coef%drdx, coef%dsdx, coef%dtdx, &
203  coef%drdy, coef%dsdy, coef%dtdy, &
204  coef%drdz, coef%dsdz, coef%dtdz, &
205  xh%w3, msh%nelv)
206  case default
207  call sx_opgrad_lx(ux, uy, uz, u, &
208  xh%dx, xh%dy, xh%dz, &
209  coef%drdx, coef%dsdx, coef%dtdx, &
210  coef%drdy, coef%dsdy, coef%dtdy, &
211  coef%drdz, coef%dsdz, coef%dtdz, &
212  xh%w3, msh%nelv, xh%lx)
213  end select
214  end associate
215 
216  end subroutine opr_sx_opgrad
217 
218  subroutine opr_sx_cdtp(dtx,x,dr,ds,dt, coef)
219  type(coef_t), intent(in) :: coef
220  real(kind=rp), dimension(coef%Xh%lxyz,coef%msh%nelv), intent(inout) :: dtx
221  real(kind=rp), dimension(coef%Xh%lxyz,coef%msh%nelv), intent(inout) :: x
222  real(kind=rp), dimension(coef%Xh%lxyz,coef%msh%nelv), intent(in) :: dr
223  real(kind=rp), dimension(coef%Xh%lxyz,coef%msh%nelv), intent(in) :: ds
224  real(kind=rp), dimension(coef%Xh%lxyz,coef%msh%nelv), intent(in) :: dt
225 
226  associate(xh => coef%Xh, msh => coef%msh, dof => coef%dof)
227  select case(xh%lx)
228  case(14)
229  call sx_cdtp_lx14(dtx, x, dr, ds, dt, &
230  xh%dxt, xh%dyt, xh%dzt, &
231  coef%B, coef%jac, msh%nelv, dof%size())
232  case(13)
233  call sx_cdtp_lx13(dtx, x, dr, ds, dt, &
234  xh%dxt, xh%dyt, xh%dzt, &
235  coef%B, coef%jac, msh%nelv, dof%size())
236  case(12)
237  call sx_cdtp_lx12(dtx, x, dr, ds, dt, &
238  xh%dxt, xh%dyt, xh%dzt, &
239  coef%B, coef%jac, msh%nelv, dof%size())
240  case(11)
241  call sx_cdtp_lx11(dtx, x, dr, ds, dt, &
242  xh%dxt, xh%dyt, xh%dzt, &
243  coef%B, coef%jac, msh%nelv, dof%size())
244  case(10)
245  call sx_cdtp_lx10(dtx, x, dr, ds, dt, &
246  xh%dxt, xh%dyt, xh%dzt, &
247  coef%B, coef%jac, msh%nelv, dof%size())
248  case(9)
249  call sx_cdtp_lx9(dtx, x, dr, ds, dt, &
250  xh%dxt, xh%dyt, xh%dzt, &
251  coef%B, coef%jac, msh%nelv, dof%size())
252  case(8)
253  call sx_cdtp_lx8(dtx, x, dr, ds, dt, &
254  xh%dxt, xh%dyt, xh%dzt, &
255  coef%B, coef%jac, msh%nelv, dof%size())
256  case(7)
257  call sx_cdtp_lx7(dtx, x, dr, ds, dt, &
258  xh%dxt, xh%dyt, xh%dzt, &
259  coef%B, coef%jac, msh%nelv, dof%size())
260  case(6)
261  call sx_cdtp_lx6(dtx, x, dr, ds, dt, &
262  xh%dxt, xh%dyt, xh%dzt, &
263  coef%B, coef%jac, msh%nelv, dof%size())
264  case(5)
265  call sx_cdtp_lx5(dtx, x, dr, ds, dt, &
266  xh%dxt, xh%dyt, xh%dzt, &
267  coef%B, coef%jac, msh%nelv, dof%size())
268  case(4)
269  call sx_cdtp_lx4(dtx, x, dr, ds, dt, &
270  xh%dxt, xh%dyt, xh%dzt, &
271  coef%B, coef%jac, msh%nelv, dof%size())
272  case(3)
273  call sx_cdtp_lx3(dtx, x, dr, ds, dt, &
274  xh%dxt, xh%dyt, xh%dzt, &
275  coef%B, coef%jac, msh%nelv, dof%size())
276  case(2)
277  call sx_cdtp_lx2(dtx, x, dr, ds, dt, &
278  xh%dxt, xh%dyt, xh%dzt, &
279  coef%B, coef%jac, msh%nelv, dof%size())
280  case default
281  call sx_cdtp_lx(dtx, x, dr, ds, dt, &
282  xh%dxt, xh%dyt, xh%dzt, &
283  coef%B, coef%jac, msh%nelv, dof%size(), xh%lx)
284  end select
285  end associate
286 
287  end subroutine opr_sx_cdtp
288 
289  subroutine opr_sx_conv1(du,u, vx, vy, vz, Xh, coef, nelv, gdim)
290  type(space_t), intent(inout) :: xh
291  type(coef_t), intent(inout) :: coef
292  integer, intent(in) :: nelv, gdim
293  real(kind=rp), intent(inout) :: du(xh%lxyz,nelv)
294  real(kind=rp), intent(inout), dimension(Xh%lx,Xh%ly,Xh%lz,nelv) :: u
295  real(kind=rp), intent(inout), dimension(Xh%lx,Xh%ly,Xh%lz,nelv) :: vx
296  real(kind=rp), intent(inout), dimension(Xh%lx,Xh%ly,Xh%lz,nelv) :: vy
297  real(kind=rp), intent(inout), dimension(Xh%lx,Xh%ly,Xh%lz,nelv) :: vz
298 
299  select case(xh%lx)
300  case(14)
301  call sx_conv1_lx14(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
302  coef%drdx, coef%dsdx, coef%dtdx, &
303  coef%drdy, coef%dsdy, coef%dtdy, &
304  coef%drdz, coef%dsdz, coef%dtdz, &
305  coef%jacinv, nelv, gdim)
306  case(13)
307  call sx_conv1_lx13(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
308  coef%drdx, coef%dsdx, coef%dtdx, &
309  coef%drdy, coef%dsdy, coef%dtdy, &
310  coef%drdz, coef%dsdz, coef%dtdz, &
311  coef%jacinv, nelv, gdim)
312  case(12)
313  call sx_conv1_lx12(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
314  coef%drdx, coef%dsdx, coef%dtdx, &
315  coef%drdy, coef%dsdy, coef%dtdy, &
316  coef%drdz, coef%dsdz, coef%dtdz, &
317  coef%jacinv, nelv, gdim)
318  case(11)
319  call sx_conv1_lx11(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
320  coef%drdx, coef%dsdx, coef%dtdx, &
321  coef%drdy, coef%dsdy, coef%dtdy, &
322  coef%drdz, coef%dsdz, coef%dtdz, &
323  coef%jacinv, nelv, gdim)
324  case(10)
325  call sx_conv1_lx10(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
326  coef%drdx, coef%dsdx, coef%dtdx, &
327  coef%drdy, coef%dsdy, coef%dtdy, &
328  coef%drdz, coef%dsdz, coef%dtdz, &
329  coef%jacinv, nelv, gdim)
330  case(9)
331  call sx_conv1_lx9(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
332  coef%drdx, coef%dsdx, coef%dtdx, &
333  coef%drdy, coef%dsdy, coef%dtdy, &
334  coef%drdz, coef%dsdz, coef%dtdz, &
335  coef%jacinv, nelv, gdim)
336  case(8)
337  call sx_conv1_lx8(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
338  coef%drdx, coef%dsdx, coef%dtdx, &
339  coef%drdy, coef%dsdy, coef%dtdy, &
340  coef%drdz, coef%dsdz, coef%dtdz, &
341  coef%jacinv, nelv, gdim)
342  case(7)
343  call sx_conv1_lx7(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
344  coef%drdx, coef%dsdx, coef%dtdx, &
345  coef%drdy, coef%dsdy, coef%dtdy, &
346  coef%drdz, coef%dsdz, coef%dtdz, &
347  coef%jacinv, nelv, gdim)
348  case(6)
349  call sx_conv1_lx6(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
350  coef%drdx, coef%dsdx, coef%dtdx, &
351  coef%drdy, coef%dsdy, coef%dtdy, &
352  coef%drdz, coef%dsdz, coef%dtdz, &
353  coef%jacinv, nelv, gdim)
354  case(5)
355  call sx_conv1_lx5(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
356  coef%drdx, coef%dsdx, coef%dtdx, &
357  coef%drdy, coef%dsdy, coef%dtdy, &
358  coef%drdz, coef%dsdz, coef%dtdz, &
359  coef%jacinv, nelv, gdim)
360  case(4)
361  call sx_conv1_lx4(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
362  coef%drdx, coef%dsdx, coef%dtdx, &
363  coef%drdy, coef%dsdy, coef%dtdy, &
364  coef%drdz, coef%dsdz, coef%dtdz, &
365  coef%jacinv, nelv, gdim)
366  case(3)
367  call sx_conv1_lx3(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
368  coef%drdx, coef%dsdx, coef%dtdx, &
369  coef%drdy, coef%dsdy, coef%dtdy, &
370  coef%drdz, coef%dsdz, coef%dtdz, &
371  coef%jacinv, nelv, gdim)
372  case(2)
373  call sx_conv1_lx2(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
374  coef%drdx, coef%dsdx, coef%dtdx, &
375  coef%drdy, coef%dsdy, coef%dtdy, &
376  coef%drdz, coef%dsdz, coef%dtdz, &
377  coef%jacinv, nelv, gdim)
378  case default
379  call sx_conv1_lx(du, u, vx, vy, vz, xh%dx, xh%dy, xh%dz, &
380  coef%drdx, coef%dsdx, coef%dtdx, &
381  coef%drdy, coef%dsdy, coef%dtdy, &
382  coef%drdz, coef%dsdz, coef%dtdz, &
383  coef%jacinv, nelv, gdim, xh%lx)
384  end select
385 
386  end subroutine opr_sx_conv1
387 
388  subroutine opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_Xh)
389  type(field_t), intent(inout) :: w1
390  type(field_t), intent(inout) :: w2
391  type(field_t), intent(inout) :: w3
392  type(field_t), intent(inout) :: u1
393  type(field_t), intent(inout) :: u2
394  type(field_t), intent(inout) :: u3
395  type(field_t), intent(inout) :: work1
396  type(field_t), intent(inout) :: work2
397  type(coef_t), intent(in) :: c_xh
398  integer :: gdim, n
399 
400  n = w1%dof%size()
401  gdim = c_xh%msh%gdim
402 
403  ! this%work1=dw/dy ; this%work2=dv/dz
404  call opr_sx_dudxyz(work1%x, u3%x, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
405  if (gdim .eq. 3) then
406  call opr_sx_dudxyz(work2%x, u2%x, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
407  call sub3(w1%x, work1%x, work2%x, n)
408  else
409  call copy(w1%x, work1%x, n)
410  end if
411  ! this%work1=du/dz ; this%work2=dw/dx
412  if (gdim .eq. 3) then
413  call opr_sx_dudxyz(work1%x, u1%x, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
414  call opr_sx_dudxyz(work2%x, u3%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
415  call sub3(w2%x, work1%x, work2%x, n)
416  else
417  call rzero (work1%x, n)
418  call opr_sx_dudxyz(work2%x, u3%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
419  call sub3(w2%x, work1%x, work2%x, n)
420  end if
421  ! this%work1=dv/dx ; this%work2=du/dy
422  call opr_sx_dudxyz(work1%x, u2%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
423  call opr_sx_dudxyz(work2%x, u1%x, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
424  call sub3(w3%x, work1%x, work2%x, n)
425  !! BC dependent, Needs to change if cyclic
426 
427  call opcolv(w1%x,w2%x,w3%x,c_xh%B, gdim, n)
428  call c_xh%gs_h%op(w1, gs_op_add)
429  call c_xh%gs_h%op(w2, gs_op_add)
430  call c_xh%gs_h%op(w3, gs_op_add)
431  call opcolv(w1%x,w2%x,w3%x,c_xh%Binv, gdim, n)
432 
433  end subroutine opr_sx_curl
434 
435  function opr_sx_cfl(dt, u, v, w, Xh, coef, nelv, gdim) result(cfl)
436  type(space_t) :: xh
437  type(coef_t) :: coef
438  integer :: nelv, gdim
439  real(kind=rp) :: dt
440  real(kind=rp), dimension(Xh%lx,Xh%ly,Xh%lz,nelv) :: u, v, w
441  real(kind=rp) :: cfl
442 
443  select case(xh%lx)
444  case (14)
445  cfl = sx_cfl_lx14(dt, u, v, w, &
446  coef%drdx, coef%dsdx, coef%dtdx, &
447  coef%drdy, coef%dsdy, coef%dtdy, &
448  coef%drdz, coef%dsdz, coef%dtdz, &
449  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
450  coef%jacinv, nelv, gdim)
451  case (13)
452  cfl = sx_cfl_lx13(dt, u, v, w, &
453  coef%drdx, coef%dsdx, coef%dtdx, &
454  coef%drdy, coef%dsdy, coef%dtdy, &
455  coef%drdz, coef%dsdz, coef%dtdz, &
456  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
457  coef%jacinv, nelv, gdim)
458  case (12)
459  cfl = sx_cfl_lx12(dt, u, v, w, &
460  coef%drdx, coef%dsdx, coef%dtdx, &
461  coef%drdy, coef%dsdy, coef%dtdy, &
462  coef%drdz, coef%dsdz, coef%dtdz, &
463  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
464  coef%jacinv, nelv, gdim)
465  case (11)
466  cfl = sx_cfl_lx11(dt, u, v, w, &
467  coef%drdx, coef%dsdx, coef%dtdx, &
468  coef%drdy, coef%dsdy, coef%dtdy, &
469  coef%drdz, coef%dsdz, coef%dtdz, &
470  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
471  coef%jacinv, nelv, gdim)
472  case (10)
473  cfl = sx_cfl_lx10(dt, u, v, w, &
474  coef%drdx, coef%dsdx, coef%dtdx, &
475  coef%drdy, coef%dsdy, coef%dtdy, &
476  coef%drdz, coef%dsdz, coef%dtdz, &
477  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
478  coef%jacinv, nelv, gdim)
479  case (9)
480  cfl = sx_cfl_lx9(dt, u, v, w, &
481  coef%drdx, coef%dsdx, coef%dtdx, &
482  coef%drdy, coef%dsdy, coef%dtdy, &
483  coef%drdz, coef%dsdz, coef%dtdz, &
484  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
485  coef%jacinv, nelv, gdim)
486  case (8)
487  cfl = sx_cfl_lx8(dt, u, v, w, &
488  coef%drdx, coef%dsdx, coef%dtdx, &
489  coef%drdy, coef%dsdy, coef%dtdy, &
490  coef%drdz, coef%dsdz, coef%dtdz, &
491  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
492  coef%jacinv, nelv, gdim)
493  case (7)
494  cfl = sx_cfl_lx7(dt, u, v, w, &
495  coef%drdx, coef%dsdx, coef%dtdx, &
496  coef%drdy, coef%dsdy, coef%dtdy, &
497  coef%drdz, coef%dsdz, coef%dtdz, &
498  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
499  coef%jacinv, nelv, gdim)
500  case (6)
501  cfl = sx_cfl_lx6(dt, u, v, w, &
502  coef%drdx, coef%dsdx, coef%dtdx, &
503  coef%drdy, coef%dsdy, coef%dtdy, &
504  coef%drdz, coef%dsdz, coef%dtdz, &
505  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
506  coef%jacinv, nelv, gdim)
507  case (5)
508  cfl = sx_cfl_lx5(dt, u, v, w, &
509  coef%drdx, coef%dsdx, coef%dtdx, &
510  coef%drdy, coef%dsdy, coef%dtdy, &
511  coef%drdz, coef%dsdz, coef%dtdz, &
512  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
513  coef%jacinv, nelv, gdim)
514  case (4)
515  cfl = sx_cfl_lx4(dt, u, v, w, &
516  coef%drdx, coef%dsdx, coef%dtdx, &
517  coef%drdy, coef%dsdy, coef%dtdy, &
518  coef%drdz, coef%dsdz, coef%dtdz, &
519  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
520  coef%jacinv, nelv, gdim)
521  case (3)
522  cfl = sx_cfl_lx3(dt, u, v, w, &
523  coef%drdx, coef%dsdx, coef%dtdx, &
524  coef%drdy, coef%dsdy, coef%dtdy, &
525  coef%drdz, coef%dsdz, coef%dtdz, &
526  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
527  coef%jacinv, nelv, gdim)
528  case (2)
529  cfl = sx_cfl_lx2(dt, u, v, w, &
530  coef%drdx, coef%dsdx, coef%dtdx, &
531  coef%drdy, coef%dsdy, coef%dtdy, &
532  coef%drdz, coef%dsdz, coef%dtdz, &
533  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
534  coef%jacinv, nelv, gdim)
535  case default
536  cfl = sx_cfl_lx(dt, u, v, w, &
537  coef%drdx, coef%dsdx, coef%dtdx, &
538  coef%drdy, coef%dsdy, coef%dtdy, &
539  coef%drdz, coef%dsdz, coef%dtdz, &
540  xh%dr_inv, xh%ds_inv, xh%dt_inv, &
541  coef%jacinv, nelv, gdim, xh%lx)
542  end select
543 
544  end function opr_sx_cfl
545 
546  subroutine opr_sx_lambda2(lambda2, u, v, w, coef)
547  type(coef_t), intent(in) :: coef
548  type(field_t), intent(inout) :: lambda2
549  type(field_t), intent(in) :: u, v, w
550 
551  associate(xh => coef%Xh, msh => coef%msh)
552  select case(xh%lx)
553  case (18)
554  call sx_lambda2_lx18(lambda2%x, u%x, v%x, w%x, &
555  xh%dx, xh%dy, xh%dz, &
556  coef%drdx, coef%dsdx, coef%dtdx, &
557  coef%drdy, coef%dsdy, coef%dtdy, &
558  coef%drdz, coef%dsdz, coef%dtdz, &
559  xh%w3, coef%B, msh%nelv)
560  case (17)
561  call sx_lambda2_lx17(lambda2%x, u%x, v%x, w%x, &
562  xh%dx, xh%dy, xh%dz, &
563  coef%drdx, coef%dsdx, coef%dtdx, &
564  coef%drdy, coef%dsdy, coef%dtdy, &
565  coef%drdz, coef%dsdz, coef%dtdz, &
566  xh%w3, coef%B, msh%nelv)
567  case (16)
568  call sx_lambda2_lx16(lambda2%x, u%x, v%x, w%x, &
569  xh%dx, xh%dy, xh%dz, &
570  coef%drdx, coef%dsdx, coef%dtdx, &
571  coef%drdy, coef%dsdy, coef%dtdy, &
572  coef%drdz, coef%dsdz, coef%dtdz, &
573  xh%w3, coef%B, msh%nelv)
574  case (15)
575  call sx_lambda2_lx15(lambda2%x, u%x, v%x, w%x, &
576  xh%dx, xh%dy, xh%dz, &
577  coef%drdx, coef%dsdx, coef%dtdx, &
578  coef%drdy, coef%dsdy, coef%dtdy, &
579  coef%drdz, coef%dsdz, coef%dtdz, &
580  xh%w3, coef%B, msh%nelv)
581  case (14)
582  call sx_lambda2_lx14(lambda2%x, u%x, v%x, w%x, &
583  xh%dx, xh%dy, xh%dz, &
584  coef%drdx, coef%dsdx, coef%dtdx, &
585  coef%drdy, coef%dsdy, coef%dtdy, &
586  coef%drdz, coef%dsdz, coef%dtdz, &
587  xh%w3, coef%B, msh%nelv)
588  case (13)
589  call sx_lambda2_lx13(lambda2%x, u%x, v%x, w%x, &
590  xh%dx, xh%dy, xh%dz, &
591  coef%drdx, coef%dsdx, coef%dtdx, &
592  coef%drdy, coef%dsdy, coef%dtdy, &
593  coef%drdz, coef%dsdz, coef%dtdz, &
594  xh%w3, coef%B, msh%nelv)
595  case (12)
596  call sx_lambda2_lx12(lambda2%x, u%x, v%x, w%x, &
597  xh%dx, xh%dy, xh%dz, &
598  coef%drdx, coef%dsdx, coef%dtdx, &
599  coef%drdy, coef%dsdy, coef%dtdy, &
600  coef%drdz, coef%dsdz, coef%dtdz, &
601  xh%w3, coef%B, msh%nelv)
602  case (11)
603  call sx_lambda2_lx11(lambda2%x, u%x, v%x, w%x, &
604  xh%dx, xh%dy, xh%dz, &
605  coef%drdx, coef%dsdx, coef%dtdx, &
606  coef%drdy, coef%dsdy, coef%dtdy, &
607  coef%drdz, coef%dsdz, coef%dtdz, &
608  xh%w3, coef%B, msh%nelv)
609  case (10)
610  call sx_lambda2_lx10(lambda2%x, u%x, v%x, w%x, &
611  xh%dx, xh%dy, xh%dz, &
612  coef%drdx, coef%dsdx, coef%dtdx, &
613  coef%drdy, coef%dsdy, coef%dtdy, &
614  coef%drdz, coef%dsdz, coef%dtdz, &
615  xh%w3, coef%B, msh%nelv)
616  case (9)
617  call sx_lambda2_lx9(lambda2%x, u%x, v%x, w%x, &
618  xh%dx, xh%dy, xh%dz, &
619  coef%drdx, coef%dsdx, coef%dtdx, &
620  coef%drdy, coef%dsdy, coef%dtdy, &
621  coef%drdz, coef%dsdz, coef%dtdz, &
622  xh%w3, coef%B, msh%nelv)
623  case (8)
624  call sx_lambda2_lx8(lambda2%x, u%x, v%x, w%x, &
625  xh%dx, xh%dy, xh%dz, &
626  coef%drdx, coef%dsdx, coef%dtdx, &
627  coef%drdy, coef%dsdy, coef%dtdy, &
628  coef%drdz, coef%dsdz, coef%dtdz, &
629  xh%w3, coef%B, msh%nelv)
630  case (7)
631  call sx_lambda2_lx7(lambda2%x, u%x, v%x, w%x, &
632  xh%dx, xh%dy, xh%dz, &
633  coef%drdx, coef%dsdx, coef%dtdx, &
634  coef%drdy, coef%dsdy, coef%dtdy, &
635  coef%drdz, coef%dsdz, coef%dtdz, &
636  xh%w3, coef%B, msh%nelv)
637  case (6)
638  call sx_lambda2_lx6(lambda2%x, u%x, v%x, w%x, &
639  xh%dx, xh%dy, xh%dz, &
640  coef%drdx, coef%dsdx, coef%dtdx, &
641  coef%drdy, coef%dsdy, coef%dtdy, &
642  coef%drdz, coef%dsdz, coef%dtdz, &
643  xh%w3, coef%B, msh%nelv)
644  case (5)
645  call sx_lambda2_lx5(lambda2%x, u%x, v%x, w%x, &
646  xh%dx, xh%dy, xh%dz, &
647  coef%drdx, coef%dsdx, coef%dtdx, &
648  coef%drdy, coef%dsdy, coef%dtdy, &
649  coef%drdz, coef%dsdz, coef%dtdz, &
650  xh%w3, coef%B, msh%nelv)
651  case (4)
652  call sx_lambda2_lx4(lambda2%x, u%x, v%x, w%x, &
653  xh%dx, xh%dy, xh%dz, &
654  coef%drdx, coef%dsdx, coef%dtdx, &
655  coef%drdy, coef%dsdy, coef%dtdy, &
656  coef%drdz, coef%dsdz, coef%dtdz, &
657  xh%w3, coef%B, msh%nelv)
658  case (3)
659  call sx_lambda2_lx3(lambda2%x, u%x, v%x, w%x, &
660  xh%dx, xh%dy, xh%dz, &
661  coef%drdx, coef%dsdx, coef%dtdx, &
662  coef%drdy, coef%dsdy, coef%dtdy, &
663  coef%drdz, coef%dsdz, coef%dtdz, &
664  xh%w3, coef%B, msh%nelv)
665  case (2)
666  call sx_lambda2_lx2(lambda2%x, u%x, v%x, w%x, &
667  xh%dx, xh%dy, xh%dz, &
668  coef%drdx, coef%dsdx, coef%dtdx, &
669  coef%drdy, coef%dsdy, coef%dtdy, &
670  coef%drdz, coef%dsdz, coef%dtdz, &
671  xh%w3, coef%B, msh%nelv)
672  case default
673  call sx_lambda2_lx(lambda2%x, u%x, v%x, w%x, &
674  xh%dx, xh%dy, xh%dz, &
675  coef%drdx, coef%dsdx, coef%dtdx, &
676  coef%drdy, coef%dsdy, coef%dtdy, &
677  coef%drdz, coef%dsdz, coef%dtdz, &
678  xh%w3, coef%B, msh%nelv, xh%lx)
679  end select
680  end associate
681 
682  end subroutine opr_sx_lambda2
683 
684 end module opr_sx
Coefficients.
Definition: coef.f90:34
Defines a field.
Definition: field.f90:34
Gather-scatter.
A simulation component that computes lambda2 The values are stored in the field registry under the na...
Definition: lambda2.f90:37
Definition: math.f90:60
Collection of vector field operations operating on and . Note that in general the indices and ....
Definition: mathops.f90:65
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Operators SX-Aurora backend.
Definition: opr_sx.f90:2
subroutine, public opr_sx_conv1(du, u, vx, vy, vz, Xh, coef, nelv, gdim)
Definition: opr_sx.f90:290
subroutine, public opr_sx_opgrad(ux, uy, uz, u, coef)
Definition: opr_sx.f90:79
subroutine, public opr_sx_dudxyz(du, u, dr, ds, dt, coef)
Definition: opr_sx.f90:25
real(kind=rp) function, public opr_sx_cfl(dt, u, v, w, Xh, coef, nelv, gdim)
Definition: opr_sx.f90:436
subroutine, public opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_Xh)
Definition: opr_sx.f90:389
subroutine, public opr_sx_cdtp(dtx, x, dr, ds, dt, coef)
Definition: opr_sx.f90:219
subroutine, public opr_sx_lambda2(lambda2, u, v, w, coef)
Definition: opr_sx.f90:547
Defines a function space.
Definition: space.f90:34
DT*X kernels for SX-Aurora.
Definition: sx_cdtp.f90:34
CFL SX-Aurora kernels.
Definition: sx_cfl.f90:34
conv1 SX-Aurora kernels
Definition: sx_conv1.f90:34
Derivative kernels for SX-Aurora.
Definition: sx_dudxyz.f90:34
subroutine, public sx_dudxyz_lx8(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:467
subroutine, public sx_dudxyz_lx4(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:707
subroutine, public sx_dudxyz_lx12(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:227
subroutine, public sx_dudxyz_lx(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd, lx)
Definition: sx_dudxyz.f90:48
subroutine, public sx_dudxyz_lx3(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:767
subroutine, public sx_dudxyz_lx6(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:587
subroutine, public sx_dudxyz_lx7(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:527
subroutine, public sx_dudxyz_lx13(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:167
subroutine, public sx_dudxyz_lx2(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:827
subroutine, public sx_dudxyz_lx14(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:107
subroutine, public sx_dudxyz_lx10(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:347
subroutine, public sx_dudxyz_lx11(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:287
subroutine, public sx_dudxyz_lx5(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:647
subroutine, public sx_dudxyz_lx9(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, nd)
Definition: sx_dudxyz.f90:407
Lambda2 kernels for SX-Aurora.
Definition: sx_lambda2.f90:34
Gradient kernels for SX-Aurora.
Definition: sx_opgrad.f90:34
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:54
The function space for the SEM solution fields.
Definition: space.f90:62