Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
pnpn_res_cpu.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, rotate_cyc
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_cpu_t
20 contains
21 procedure, nopass :: compute => pnpn_prs_res_cpu_compute
22 end type pnpn_prs_res_cpu_t
23
24 type, public, extends(pnpn_vel_res_t) :: pnpn_vel_res_cpu_t
25 contains
26 procedure, nopass :: compute => pnpn_vel_res_cpu_compute
27 end type pnpn_vel_res_cpu_t
28
29contains
30
31 subroutine pnpn_prs_res_cpu_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), .false.)
55 call neko_scratch_registry%request_field(ta2, temp_indices(2), .false.)
56 call neko_scratch_registry%request_field(ta3, temp_indices(3), .false.)
57 call neko_scratch_registry%request_field(wa1, temp_indices(4), .false.)
58 call neko_scratch_registry%request_field(wa2, temp_indices(5), .false.)
59 call neko_scratch_registry%request_field(wa3, temp_indices(6), .false.)
60 call neko_scratch_registry%request_field(work1, temp_indices(7), .false.)
61 call neko_scratch_registry%request_field(work2, temp_indices(8), .false.)
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 ! ta = f / rho - wa * mu / rho * B
78 do concurrent(i = 1:n)
79 ta1%x(i,1,1,1) = f_x%x(i,1,1,1) / rho_val &
80 - ((wa1%x(i,1,1,1) * (mu_val / rho_val)) * c_xh%B(i,1,1,1))
81 ta2%x(i,1,1,1) = f_y%x(i,1,1,1) / rho_val &
82 - ((wa2%x(i,1,1,1) * (mu_val / rho_val)) * c_xh%B(i,1,1,1))
83 ta3%x(i,1,1,1) = f_z%x(i,1,1,1) / rho_val &
84 - ((wa3%x(i,1,1,1) * (mu_val / rho_val)) * c_xh%B(i,1,1,1))
85 end do
86
87 call rotate_cyc(ta1%x, ta2%x, ta3%x, 1, c_xh)
88 call gs_xh%op(ta1, gs_op_add)
89 call gs_xh%op(ta2, gs_op_add)
90 call gs_xh%op(ta3, gs_op_add)
91 call rotate_cyc(ta1%x, ta2%x, ta3%x, 0, c_xh)
92
93 do concurrent(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 cdtp(wa1%x, ta1%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
100 call cdtp(wa2%x, ta2%x, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
101 call cdtp(wa3%x, ta3%x, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
102
103 call ax%compute(p_res%x, p%x, c_xh, p%msh, p%Xh)
104
105 do concurrent(i = 1:n)
106 p_res%x(i,1,1,1) = (-p_res%x(i,1,1,1)) &
107 + wa1%x(i,1,1,1) + wa2%x(i,1,1,1) + wa3%x(i,1,1,1)
108 end do
109
110 !
111 ! Surface velocity terms
112 !
113 do concurrent(i = 1:n)
114 wa1%x(i,1,1,1) = 0.0_rp
115 wa2%x(i,1,1,1) = 0.0_rp
116 wa3%x(i,1,1,1) = 0.0_rp
117 end do
118
119 call bc_sym_surface%apply_surfvec(wa1%x, wa2%x, wa3%x, ta1%x, ta2%x, ta3%x,&
120 n)
121
122 dtbd = bd / dt
123 do concurrent(i = 1:n)
124 ta1%x(i,1,1,1) = 0.0_rp
125 ta2%x(i,1,1,1) = 0.0_rp
126 ta3%x(i,1,1,1) = 0.0_rp
127 end do
128
129 call bc_prs_surface%apply_surfvec(ta1%x, ta2%x, ta3%x, u%x, v%x, w%x, n)
130
131 do concurrent(i = 1:n)
132 p_res%x(i,1,1,1) = p_res%x(i,1,1,1) &
133 - (dtbd * (ta1%x(i,1,1,1) + ta2%x(i,1,1,1) + ta3%x(i,1,1,1)))&
134 - (wa1%x(i,1,1,1) + wa2%x(i,1,1,1) + wa3%x(i,1,1,1))
135 end do
136
137 call neko_scratch_registry%relinquish_field(temp_indices)
138
139 end subroutine pnpn_prs_res_cpu_compute
140
141 subroutine pnpn_vel_res_cpu_compute(Ax, u, v, w, u_res, v_res, w_res, &
142 p, f_x, f_y, f_z, c_Xh, msh, Xh, mu, rho, bd, dt, n)
143 class(ax_t), intent(in) :: Ax
144 type(mesh_t), intent(inout) :: msh
145 type(space_t), intent(inout) :: Xh
146 type(field_t), intent(inout) :: p, u, v, w
147 type(field_t), intent(inout) :: u_res, v_res, w_res
148 type(field_t), intent(in) :: f_x, f_y, f_z
149 type(coef_t), intent(inout) :: c_Xh
150 type(field_t), intent(in) :: mu
151 type(field_t), intent(in) :: rho
152 real(kind=rp), intent(in) :: bd
153 real(kind=rp), intent(in) :: dt
154 real(kind=rp) :: rho_val, mu_val
155 integer :: temp_indices(3)
156 type(field_t), pointer :: ta1, ta2, ta3
157 integer, intent(in) :: n
158 integer :: i
159
160 ! We assume the material properties are constant
161 rho_val = rho%x(1,1,1,1)
162 mu_val = mu%x(1,1,1,1)
163
164 do concurrent(i = 1:n)
165 c_xh%h1(i,1,1,1) = mu_val
166 c_xh%h2(i,1,1,1) = rho_val * bd / dt
167 end do
168 c_xh%ifh2 = .true.
169
170 call ax%compute(u_res%x, u%x, c_xh, msh, xh)
171 call ax%compute(v_res%x, v%x, c_xh, msh, xh)
172 call ax%compute(w_res%x, w%x, c_xh, msh, xh)
173 call neko_scratch_registry%request_field(ta1, temp_indices(1), .false.)
174 call neko_scratch_registry%request_field(ta2, temp_indices(2), .false.)
175 call neko_scratch_registry%request_field(ta3, temp_indices(3), .false.)
176
177 call opgrad(ta1%x, ta2%x, ta3%x, p%x, c_xh)
178
179 do concurrent(i = 1:n)
180 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)
181 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)
182 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)
183 end do
184
185 call neko_scratch_registry%relinquish_field(temp_indices)
186
187 end subroutine pnpn_vel_res_cpu_compute
188
189end module pnpn_res_cpu
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:423
subroutine, public invers2(a, b, n)
Compute inverted vector .
Definition math.f90:639
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:249
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:205
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 (CPU version)
subroutine pnpn_vel_res_cpu_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_cpu_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 objects This can be used when you have a func...
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:56
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:63