Neko  0.9.99
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 gather_scatter, only : gs_t, gs_op_add
4  use interpolation, only : interpolator_t
5  use num_types, only : rp
6  use space, only : space_t
7  use coefs, only : coef_t
8  use math, only : sub3, copy, rzero
9  use field, only : field_t
10  use mathops, only : opcolv
11  implicit none
12  private
13 
14  public :: opr_sx_dudxyz, opr_sx_opgrad, opr_sx_cdtp, opr_sx_conv1, &
15  opr_sx_curl, opr_sx_cfl, opr_sx_lambda2, opr_sx_convect_scalar, &
16  opr_sx_set_convect_rst
17 
18 
19  interface
20  module subroutine opr_sx_dudxyz(du, u, dr, ds, dt, coef)
21  type(coef_t), intent(in), target :: coef
22  real(kind=rp), intent(inout), &
23  dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, coef%msh%nelv) :: du
24  real(kind=rp), intent(in), &
25  dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, coef%msh%nelv) :: &
26  u, dr, ds, dt
27  end subroutine opr_sx_dudxyz
28 
29  module subroutine opr_sx_opgrad(ux, uy, uz, u, coef)
30  type(coef_t), intent(in) :: coef
31  real(kind=rp), intent(inout) :: ux(coef%Xh%lxyz, coef%msh%nelv)
32  real(kind=rp), intent(inout) :: uy(coef%Xh%lxyz, coef%msh%nelv)
33  real(kind=rp), intent(inout) :: uz(coef%Xh%lxyz, coef%msh%nelv)
34  real(kind=rp), intent(in) :: u(coef%Xh%lxyz, coef%msh%nelv)
35  end subroutine opr_sx_opgrad
36 
37  module subroutine opr_sx_cdtp(dtx, x, dr, ds, dt, coef)
38  type(coef_t), intent(in) :: coef
39  real(kind=rp), intent(inout) :: dtx(coef%Xh%lxyz, coef%msh%nelv)
40  real(kind=rp), intent(inout) :: x(coef%Xh%lxyz, coef%msh%nelv)
41  real(kind=rp), intent(in) :: dr(coef%Xh%lxyz, coef%msh%nelv)
42  real(kind=rp), intent(in) :: ds(coef%Xh%lxyz, coef%msh%nelv)
43  real(kind=rp), intent(in) :: dt(coef%Xh%lxyz, coef%msh%nelv)
44  end subroutine opr_sx_cdtp
45 
46  module subroutine opr_sx_conv1(du, u, vx, vy, vz, xh, coef, nelv)
47  type(space_t), intent(inout) :: Xh
48  type(coef_t), intent(inout) :: coef
49  integer, intent(in) :: nelv
50  real(kind=rp), intent(inout) :: du(xh%lxyz, nelv)
51  real(kind=rp), intent(inout) :: u(xh%lx, xh%ly, xh%lz, nelv)
52  real(kind=rp), intent(inout) :: vx(xh%lx, xh%ly, xh%lz, nelv)
53  real(kind=rp), intent(inout) :: vy(xh%lx, xh%ly, xh%lz, nelv)
54  real(kind=rp), intent(inout) :: vz(xh%lx, xh%ly, xh%lz, nelv)
55  end subroutine opr_sx_conv1
56 
57  module subroutine opr_sx_convect_scalar(du, u, c, xh_gll, xh_gl, &
58  coef_gll, coef_gl, gll_to_gl)
59  type(space_t), intent(in) :: Xh_GL
60  type(space_t), intent(in) :: Xh_GLL
61  type(coef_t), intent(in) :: coef_GLL
62  type(coef_t), intent(in) :: coef_GL
63  type(interpolator_t), intent(inout) :: GLL_to_GL
64  real(kind=rp), intent(inout) :: &
65  du(xh_gll%lx, xh_gll%ly, xh_gll%lz, coef_gl%msh%nelv)
66  real(kind=rp), intent(inout) :: &
67  u(xh_gl%lx, xh_gl%lx, xh_gl%lx, coef_gl%msh%nelv)
68  real(kind=rp), intent(inout) :: c(xh_gl%lxyz, coef_gl%msh%nelv, 3)
69 
70  end subroutine opr_sx_convect_scalar
71 
72  module function opr_sx_cfl(dt, u, v, w, xh, coef, nelv) result(cfl)
73  type(space_t), intent(in) :: Xh
74  type(coef_t), intent(in) :: coef
75  integer, intent(in) :: nelv
76  real(kind=rp), intent(in) :: dt
77  real(kind=rp), dimension(Xh%lx, Xh%ly, Xh%lz, nelv) :: u, v, w
78  real(kind=rp) :: cfl
79  end function opr_sx_cfl
80 
81  module subroutine opr_sx_lambda2(lambda2, u, v, w, coef)
82  type(coef_t), intent(in) :: coef
83  type(field_t), intent(inout) :: lambda2
84  type(field_t), intent(in) :: u, v, w
85  end subroutine opr_sx_lambda2
86 
87  module subroutine opr_sx_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
88  type(space_t), intent(inout) :: Xh
89  type(coef_t), intent(inout) :: coef
90  real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
91  intent(inout) :: cr, cs, ct
92  real(kind=rp), dimension(Xh%lxyz, coef%msh%nelv), &
93  intent(in) :: cx, cy, cz
94  end subroutine opr_sx_set_convect_rst
95  end interface
96 
97 contains
98 
99  subroutine opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_Xh)
100  type(field_t), intent(inout) :: w1
101  type(field_t), intent(inout) :: w2
102  type(field_t), intent(inout) :: w3
103  type(field_t), intent(inout) :: u1
104  type(field_t), intent(inout) :: u2
105  type(field_t), intent(inout) :: u3
106  type(field_t), intent(inout) :: work1
107  type(field_t), intent(inout) :: work2
108  type(coef_t), intent(in) :: c_xh
109  integer :: gdim, n
110 
111  n = w1%dof%size()
112  gdim = c_xh%msh%gdim
113 
114  ! this%work1=dw/dy ; this%work2=dv/dz
115  call opr_sx_dudxyz(work1%x, u3%x, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
116  if (gdim .eq. 3) then
117  call opr_sx_dudxyz(work2%x, u2%x, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
118  call sub3(w1%x, work1%x, work2%x, n)
119  else
120  call copy(w1%x, work1%x, n)
121  end if
122  ! this%work1=du/dz ; this%work2=dw/dx
123  if (gdim .eq. 3) then
124  call opr_sx_dudxyz(work1%x, u1%x, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
125  call opr_sx_dudxyz(work2%x, u3%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
126  call sub3(w2%x, work1%x, work2%x, n)
127  else
128  call rzero (work1%x, n)
129  call opr_sx_dudxyz(work2%x, u3%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
130  call sub3(w2%x, work1%x, work2%x, n)
131  end if
132  ! this%work1=dv/dx ; this%work2=du/dy
133  call opr_sx_dudxyz(work1%x, u2%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
134  call opr_sx_dudxyz(work2%x, u1%x, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
135  call sub3(w3%x, work1%x, work2%x, n)
136  !! BC dependent, Needs to change if cyclic
137 
138  call opcolv(w1%x, w2%x, w3%x, c_xh%B, gdim, n)
139  call c_xh%gs_h%op(w1, gs_op_add)
140  call c_xh%gs_h%op(w2, gs_op_add)
141  call c_xh%gs_h%op(w3, gs_op_add)
142  call opcolv(w1%x, w2%x, w3%x, c_xh%Binv, gdim, n)
143 
144  end subroutine opr_sx_curl
145 
146 
147 
148 
149 
150 end module opr_sx
Coefficients.
Definition: coef.f90:34
Defines a field.
Definition: field.f90:34
Gather-scatter.
Routines to interpolate between different spaces.
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
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition: math.f90:642
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:239
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:195
Collection of vector field operations operating on and . Note that in general the indices and ....
Definition: mathops.f90:65
subroutine, public opcolv(a1, a2, a3, c, gdim, n)
Definition: mathops.f90:97
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_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_Xh)
Definition: opr_sx.f90:100
Defines a function space.
Definition: space.f90:34
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:55
Interpolation between two space::space_t.
The function space for the SEM solution fields.
Definition: space.f90:62