Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
pnpn_res_sx.f90
Go to the documentation of this file.
1
3 use gather_scatter, only : gs_t, gs_op_add
4 use operators, only : opgrad, curl, cdtp
5 use field, only : field_t
6 use ax_product, only : ax_t
7 use coefs, only : coef_t
11 use mesh, only : mesh_t
12 use num_types, only : rp
13 use space, only : space_t
14 use math, only : copy, cmult2, invers2, rzero
15 use, intrinsic :: iso_c_binding, only : c_ptr
16 implicit none
17 private
18
19 type, public, extends(pnpn_prs_res_t) :: pnpn_prs_res_sx_t
20 contains
21 procedure, nopass :: compute => pnpn_prs_res_sx_compute
22 end type pnpn_prs_res_sx_t
23
24 type, public, extends(pnpn_vel_res_t) :: pnpn_vel_res_sx_t
25 contains
26 procedure, nopass :: compute => pnpn_vel_res_sx_compute
27 end type pnpn_vel_res_sx_t
28
29contains
30
31 subroutine pnpn_prs_res_sx_compute(p, p_res, u, v, w, u_e, v_e, w_e, f_x, &
32 f_y, f_z, c_Xh, gs_Xh, bc_prs_surface, bc_sym_surface, Ax, bd, dt, mu, &
33 rho, event)
34 type(field_t), intent(inout) :: p, u, v, w
35 type(field_t), intent(in) :: u_e, v_e, w_e
36 type(field_t), intent(inout) :: p_res
37 type(field_t), intent(in) :: f_x, f_y, f_z
38 type(coef_t), intent(inout) :: c_Xh
39 type(gs_t), intent(inout) :: gs_Xh
40 type(facet_normal_t), intent(in) :: bc_prs_surface
41 type(facet_normal_t), intent(in) :: bc_sym_surface
42 class(ax_t), intent(inout) :: Ax
43 real(kind=rp), intent(in) :: bd
44 real(kind=rp), intent(in) :: dt
45 type(field_t), intent(in) :: mu
46 type(field_t), intent(in) :: rho
47 type(c_ptr), intent(inout) :: event
48 real(kind=rp) :: dtbd, rho_val, mu_val
49 integer :: n
50 integer :: i
51 type(field_t), pointer :: ta1, ta2, ta3, wa1, wa2, wa3, work1, work2
52 integer :: temp_indices(8)
53
54 call neko_scratch_registry%request_field(ta1, temp_indices(1))
55 call neko_scratch_registry%request_field(ta2, temp_indices(2))
56 call neko_scratch_registry%request_field(ta3, temp_indices(3))
57 call neko_scratch_registry%request_field(wa1, temp_indices(4))
58 call neko_scratch_registry%request_field(wa2, temp_indices(5))
59 call neko_scratch_registry%request_field(wa3, temp_indices(6))
60 call neko_scratch_registry%request_field(work1, temp_indices(7))
61 call neko_scratch_registry%request_field(work2, temp_indices(8))
62
63 n = c_xh%dof%size()
64
65 ! We assume the material properties are constant
66 rho_val = rho%x(1,1,1,1)
67 mu_val = mu%x(1,1,1,1)
68 do i = 1, n
69 c_xh%h1(i,1,1,1) = 1.0_rp / rho_val
70 c_xh%h2(i,1,1,1) = 0.0_rp
71 end do
72 c_xh%ifh2 = .false.
73
74 call curl(ta1, ta2, ta3, u_e, v_e, w_e, work1, work2, c_xh)
75 call curl(wa1, wa2, wa3, ta1, ta2, ta3, work1, work2, c_xh)
76
77 do i = 1, n
78 wa1%x(i,1,1,1) = (wa1%x(i,1,1,1) * mu_val / rho_val) * c_xh%B(i,1,1,1)
79 wa2%x(i,1,1,1) = (wa2%x(i,1,1,1) * mu_val / rho_val) * c_xh%B(i,1,1,1)
80 wa3%x(i,1,1,1) = (wa3%x(i,1,1,1) * mu_val / rho_val) * c_xh%B(i,1,1,1)
81 end do
82
83 do i = 1, n
84 ta1%x(i,1,1,1) = f_x%x(i,1,1,1) / rho_val - wa1%x(i,1,1,1)
85 ta2%x(i,1,1,1) = f_y%x(i,1,1,1) / rho_val - wa2%x(i,1,1,1)
86 ta3%x(i,1,1,1) = f_z%x(i,1,1,1) / rho_val - wa3%x(i,1,1,1)
87 end do
88
89 call gs_xh%op(ta1, gs_op_add)
90 call gs_xh%op(ta2, gs_op_add)
91 call gs_xh%op(ta3, gs_op_add)
92
93 do i = 1, n
94 ta1%x(i,1,1,1) = ta1%x(i,1,1,1) * c_xh%Binv(i,1,1,1)
95 ta2%x(i,1,1,1) = ta2%x(i,1,1,1) * c_xh%Binv(i,1,1,1)
96 ta3%x(i,1,1,1) = ta3%x(i,1,1,1) * c_xh%Binv(i,1,1,1)
97 end do
98
99 call ax%compute(p_res%x, p%x, c_xh, p%msh, p%Xh)
100
101 call cdtp(wa1%x, ta1%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
102 call cdtp(wa2%x, ta2%x, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
103 call cdtp(wa3%x, ta3%x, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
104 do i = 1, n
105 p_res%x(i,1,1,1) = (-p_res%x(i,1,1,1)) &
106 + wa1%x(i,1,1,1) + wa2%x(i,1,1,1) + wa3%x(i,1,1,1)
107 end do
108
109 !
110 ! Surface velocity terms
111 !
112 do i = 1, n
113 wa1%x(i,1,1,1) = 0.0_rp
114 wa2%x(i,1,1,1) = 0.0_rp
115 wa3%x(i,1,1,1) = 0.0_rp
116 end do
117
118 call bc_sym_surface%apply_surfvec(wa1%x, wa2%x, wa3%x, ta1%x, ta2%x, ta3%x,&
119 n)
120
121 dtbd = bd / dt
122 do i = 1, n
123 ta1%x(i,1,1,1) = 0.0_rp
124 ta2%x(i,1,1,1) = 0.0_rp
125 ta3%x(i,1,1,1) = 0.0_rp
126 end do
127
128 call bc_prs_surface%apply_surfvec(ta1%x, ta2%x, ta3%x, u%x, v%x, w%x, n)
129
130 do i = 1, n
131 p_res%x(i,1,1,1) = p_res%x(i,1,1,1) &
132 - (dtbd * (ta1%x(i,1,1,1) + ta2%x(i,1,1,1) + ta3%x(i,1,1,1)))&
133 - (wa1%x(i,1,1,1) + wa2%x(i,1,1,1) + wa3%x(i,1,1,1))
134 end do
135
136 call neko_scratch_registry%relinquish_field(temp_indices)
137 end subroutine pnpn_prs_res_sx_compute
138
139 subroutine pnpn_vel_res_sx_compute(Ax, u, v, w, u_res, v_res, w_res, &
140 p, f_x, f_y, f_z, c_Xh, msh, Xh, mu, rho, bd, dt, n)
141 class(ax_t), intent(in) :: Ax
142 type(mesh_t), intent(inout) :: msh
143 type(space_t), intent(inout) :: Xh
144 type(field_t), intent(inout) :: p, u, v, w
145 type(field_t), intent(inout) :: u_res, v_res, w_res
146 type(field_t), intent(in) :: f_x, f_y, f_z
147 type(coef_t), intent(inout) :: c_Xh
148 type(field_t), intent(in) :: mu
149 type(field_t), intent(in) :: rho
150 real(kind=rp), intent(in) :: bd
151 real(kind=rp), intent(in) :: dt
152 integer, intent(in) :: n
153 integer :: temp_indices(3)
154 type(field_t), pointer :: ta1, ta2, ta3
155 integer :: i
156 real(kind=rp) :: rho_val, mu_val
157
158 rho_val = rho%x(1,1,1,1)
159 mu_val = mu%x(1,1,1,1)
160
161 do i = 1, n
162 c_xh%h1(i,1,1,1) = mu_val
163 c_xh%h2(i,1,1,1) = rho_val * bd / dt
164 end do
165
166 c_xh%ifh2 = .true.
167
168 call ax%compute(u_res%x, u%x, c_xh, msh, xh)
169 call ax%compute(v_res%x, v%x, c_xh, msh, xh)
170 call ax%compute(w_res%x, w%x, c_xh, msh, xh)
171
172 call neko_scratch_registry%request_field(ta1, temp_indices(1))
173 call neko_scratch_registry%request_field(ta2, temp_indices(2))
174 call neko_scratch_registry%request_field(ta3, temp_indices(3))
175
176 call opgrad(ta1%x, ta2%x, ta3%x, p%x, c_xh)
177
178 do i = 1, n
179 u_res%x(i,1,1,1) = (-u_res%x(i,1,1,1)) - ta1%x(i,1,1,1) + f_x%x(i,1,1,1)
180 v_res%x(i,1,1,1) = (-v_res%x(i,1,1,1)) - ta2%x(i,1,1,1) + f_y%x(i,1,1,1)
181 w_res%x(i,1,1,1) = (-w_res%x(i,1,1,1)) - ta3%x(i,1,1,1) + f_z%x(i,1,1,1)
182 end do
183
184 call neko_scratch_registry%relinquish_field(temp_indices)
185 end subroutine pnpn_vel_res_sx_compute
186
187end module pnpn_res_sx
Defines a Matrix-vector product.
Definition ax.f90:34
Coefficients.
Definition coef.f90:34
Dirichlet condition applied in the facet normal direction.
Defines a field.
Definition field.f90:34
Gather-scatter.
Definition math.f90:60
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:700
subroutine, public invers2(a, b, n)
Compute inverted vector .
Definition math.f90:499
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:238
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:194
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public opgrad(ux, uy, uz, u, coef, es, ee)
Compute the weak gradient of a scalar field, i.e. the gradient multiplied by the mass matrix.
subroutine, public curl(w1, w2, w3, u1, u2, u3, work1, work2, coef, event)
subroutine, public cdtp(dtx, x, dr, ds, dt, coef, es, ee)
Apply D^T to a scalar field, where D is the derivative matrix.
Residuals in the Pn-Pn formulation (SX version)
subroutine pnpn_vel_res_sx_compute(ax, u, v, w, u_res, v_res, w_res, p, f_x, f_y, f_z, c_xh, msh, xh, mu, rho, bd, dt, n)
subroutine pnpn_prs_res_sx_compute(p, p_res, u, v, w, u_e, v_e, w_e, f_x, f_y, f_z, c_xh, gs_xh, bc_prs_surface, bc_sym_surface, ax, bd, dt, mu, rho, event)
Defines Pressure and velocity residuals in the Pn-Pn formulation.
Definition pnpn_res.f90:34
Defines a registry for storing and requesting temporary fields This can be used when you have a funct...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Defines a function space.
Definition space.f90:34
Base type for a matrix-vector product providing .
Definition ax.f90:43
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Dirichlet condition in facet normal direction.
Abstract type to compute pressure residual.
Definition pnpn_res.f90:48
Abstract type to compute velocity residual.
Definition pnpn_res.f90:54
The function space for the SEM solution fields.
Definition space.f90:62