Neko 1.99.4
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
pnpn_res_stress_device.F90
Go to the documentation of this file.
1
3 use gather_scatter, only : gs_t, gs_op_add
4 use utils, only : neko_error
6 use field, only : field_t
7 use ax_product, only : ax_t
8 use coefs, only : coef_t
12 use mesh, only : mesh_t
13 use num_types, only : rp, c_rp
14 use space, only : space_t
15 use, intrinsic :: iso_c_binding, only : c_ptr, c_int
19 implicit none
20 private
21
25 contains
26 procedure, nopass :: compute => pnpn_prs_res_stress_device_compute
28
32 contains
33 procedure, nopass :: compute => pnpn_vel_res_stress_device_compute
35
36#ifdef HAVE_HIP
37 interface
38 subroutine pnpn_prs_stress_res_part1_hip(ta1_d, ta2_d, ta3_d, &
39 wa1_d, wa2_d, wa3_d, s11_d, s22_d, s33_d, &
40 s12_d, s13_d, s23_d, f_u_d, f_v_d, f_w_d, &
41 B_d, h1_d, rho_d, n) &
42 bind(c, name = 'pnpn_prs_stress_res_part1_hip')
43 use, intrinsic :: iso_c_binding
44 import c_rp
45 implicit none
46 type(c_ptr), value :: ta1_d, ta2_d, ta3_d
47 type(c_ptr), value :: wa1_d, wa2_d, wa3_d
48 type(c_ptr), value :: s11_d, s22_d, s33_d, s12_d, s13_d, s23_d
49 type(c_ptr), value :: f_u_d, f_v_d, f_w_d
50 type(c_ptr), value :: B_d, h1_d, rho_d
51 integer(c_int) :: n
53 end interface
54
55 interface
56 subroutine pnpn_prs_res_part2_hip(p_res_d, wa1_d, wa2_d, wa3_d, n) &
57 bind(c, name = 'pnpn_prs_res_part2_hip')
58 use, intrinsic :: iso_c_binding
59 implicit none
60 type(c_ptr), value :: p_res_d, wa1_d, wa2_d, wa3_d
61 integer(c_int) :: n
62 end subroutine pnpn_prs_res_part2_hip
63 end interface
64
65 interface
66 subroutine pnpn_prs_stress_res_part3_hip(p_res_d, ta1_d, ta2_d, ta3_d, &
67 wa1_d, wa2_d, wa3_d, dtbd, n) &
68 bind(c, name = 'pnpn_prs_stress_res_part3_hip')
69 use, intrinsic :: iso_c_binding
70 import c_rp
71 implicit none
72 type(c_ptr), value :: p_res_d, ta1_d, ta2_d, ta3_d
73 type(c_ptr), value :: wa1_d, wa2_d, wa3_d
74 real(c_rp) :: dtbd
75 integer(c_int) :: n
77 end interface
78
79 interface
80 subroutine pnpn_vel_res_update_hip(u_res_d, v_res_d, w_res_d, &
81 ta1_d, ta2_d, ta3_d, f_u_d, f_v_d, f_w_d, n) &
82 bind(c, name = 'pnpn_vel_res_update_hip')
83 use, intrinsic :: iso_c_binding
84 implicit none
85 type(c_ptr), value :: u_res_d, v_res_d, w_res_d
86 type(c_ptr), value :: ta1_d, ta2_d, ta3_d
87 type(c_ptr), value :: f_u_d, f_v_d, f_w_d
88 integer(c_int) :: n
89 end subroutine pnpn_vel_res_update_hip
90 end interface
91#elif HAVE_CUDA
92 interface
93 subroutine pnpn_prs_stress_res_part1_cuda(ta1_d, ta2_d, ta3_d, &
94 wa1_d, wa2_d, wa3_d, s11_d, s22_d, s33_d, &
95 s12_d, s13_d, s23_d, f_u_d, f_v_d, f_w_d, &
96 B_d, h1_d, rho_d, n) &
97 bind(c, name = 'pnpn_prs_stress_res_part1_cuda')
98 use, intrinsic :: iso_c_binding
99 import c_rp
100 implicit none
101 type(c_ptr), value :: ta1_d, ta2_d, ta3_d
102 type(c_ptr), value :: wa1_d, wa2_d, wa3_d
103 type(c_ptr), value :: s11_d, s22_d, s33_d, s12_d, s13_d, s23_d
104 type(c_ptr), value :: f_u_d, f_v_d, f_w_d
105 type(c_ptr), value :: B_d, h1_d, rho_d
106 integer(c_int) :: n
107 end subroutine pnpn_prs_stress_res_part1_cuda
108 end interface
109
110 interface
111 subroutine pnpn_prs_res_part2_cuda(p_res_d, wa1_d, wa2_d, wa3_d, n) &
112 bind(c, name = 'pnpn_prs_res_part2_cuda')
113 use, intrinsic :: iso_c_binding
114 implicit none
115 type(c_ptr), value :: p_res_d, wa1_d, wa2_d, wa3_d
116 integer(c_int) :: n
117 end subroutine pnpn_prs_res_part2_cuda
118 end interface
119
120 interface
121 subroutine pnpn_prs_stress_res_part3_cuda(p_res_d, ta1_d, ta2_d, ta3_d, &
122 wa1_d, wa2_d, wa3_d, dtbd, n) &
123 bind(c, name = 'pnpn_prs_stress_res_part3_cuda')
124 use, intrinsic :: iso_c_binding
125 import c_rp
126 implicit none
127 type(c_ptr), value :: p_res_d, ta1_d, ta2_d, ta3_d
128 type(c_ptr), value :: wa1_d, wa2_d, wa3_d
129 real(c_rp) :: dtbd
130 integer(c_int) :: n
131 end subroutine pnpn_prs_stress_res_part3_cuda
132 end interface
133
134 interface
135 subroutine pnpn_vel_res_update_cuda(u_res_d, v_res_d, w_res_d, &
136 ta1_d, ta2_d, ta3_d, f_u_d, f_v_d, f_w_d, n) &
137 bind(c, name = 'pnpn_vel_res_update_cuda')
138 use, intrinsic :: iso_c_binding
139 implicit none
140 type(c_ptr), value :: u_res_d, v_res_d, w_res_d
141 type(c_ptr), value :: ta1_d, ta2_d, ta3_d
142 type(c_ptr), value :: f_u_d, f_v_d, f_w_d
143 integer(c_int) :: n
144 end subroutine pnpn_vel_res_update_cuda
145 end interface
146#elif HAVE_OPENCL
147 interface
148 subroutine pnpn_prs_stress_res_part1_opencl(ta1_d, ta2_d, ta3_d, &
149 wa1_d, wa2_d, wa3_d, s11_d, s22_d, s33_d, &
150 s12_d, s13_d, s23_d, f_u_d, f_v_d, f_w_d, &
151 B_d, h1_d, rho_d, n) &
152 bind(c, name = 'pnpn_prs_stress_res_part1_opencl')
153 use, intrinsic :: iso_c_binding
154 import c_rp
155 implicit none
156 type(c_ptr), value :: ta1_d, ta2_d, ta3_d
157 type(c_ptr), value :: wa1_d, wa2_d, wa3_d
158 type(c_ptr), value :: s11_d, s22_d, s33_d, s12_d, s13_d, s23_d
159 type(c_ptr), value :: f_u_d, f_v_d, f_w_d
160 type(c_ptr), value :: B_d, h1_d, rho_d
161 integer(c_int) :: n
163 end interface
164
165 interface
166 subroutine pnpn_prs_res_part2_opencl(p_res_d, wa1_d, wa2_d, wa3_d, n) &
167 bind(c, name = 'pnpn_prs_res_part2_opencl')
168 use, intrinsic :: iso_c_binding
169 implicit none
170 type(c_ptr), value :: p_res_d, wa1_d, wa2_d, wa3_d
171 integer(c_int) :: n
172 end subroutine pnpn_prs_res_part2_opencl
173 end interface
174
175 interface
176 subroutine pnpn_prs_stress_res_part3_opencl(p_res_d, ta1_d, ta2_d, ta3_d, &
177 wa1_d, wa2_d, wa3_d, dtbd, n) &
178 bind(c, name = 'pnpn_prs_stress_res_part3_opencl')
179 use, intrinsic :: iso_c_binding
180 import c_rp
181 implicit none
182 type(c_ptr), value :: p_res_d, ta1_d, ta2_d, ta3_d
183 type(c_ptr), value :: wa1_d, wa2_d, wa3_d
184 real(c_rp) :: dtbd
185 integer(c_int) :: n
187 end interface
188
189
190 interface
191 subroutine pnpn_vel_res_update_opencl(u_res_d, v_res_d, w_res_d, &
192 ta1_d, ta2_d, ta3_d, f_u_d, f_v_d, f_w_d, n) &
193 bind(c, name = 'pnpn_vel_res_update_opencl')
194 use, intrinsic :: iso_c_binding
195 implicit none
196 type(c_ptr), value :: u_res_d, v_res_d, w_res_d
197 type(c_ptr), value :: ta1_d, ta2_d, ta3_d
198 type(c_ptr), value :: f_u_d, f_v_d, f_w_d
199 integer(c_int) :: n
200 end subroutine pnpn_vel_res_update_opencl
201 end interface
202#elif HAVE_METAL
203 interface
204 subroutine pnpn_prs_res_part2_metal(p_res_d, wa1_d, wa2_d, wa3_d, n) &
205 bind(c, name = 'pnpn_prs_res_part2_metal')
206 use, intrinsic :: iso_c_binding
207 implicit none
208 type(c_ptr), value :: p_res_d, wa1_d, wa2_d, wa3_d
209 integer(c_int) :: n
210 end subroutine pnpn_prs_res_part2_metal
211 end interface
212
213
214 interface
215 subroutine pnpn_vel_res_update_metal(u_res_d, v_res_d, w_res_d, &
216 ta1_d, ta2_d, ta3_d, f_u_d, f_v_d, f_w_d, n) &
217 bind(c, name = 'pnpn_vel_res_update_metal')
218 use, intrinsic :: iso_c_binding
219 implicit none
220 type(c_ptr), value :: u_res_d, v_res_d, w_res_d
221 type(c_ptr), value :: ta1_d, ta2_d, ta3_d
222 type(c_ptr), value :: f_u_d, f_v_d, f_w_d
223 integer(c_int) :: n
224 end subroutine pnpn_vel_res_update_metal
225 end interface
226#endif
227
228contains
229
230 subroutine pnpn_prs_res_stress_device_compute(p, p_res, u, v, w, u_e, v_e,&
231 w_e, f_x, f_y, f_z, c_Xh, gs_Xh, bc_prs_surface, bc_sym_surface, Ax, bd,&
232 dt, mu, rho, event)
233 type(field_t), intent(inout) :: p, u, v, w
234 type(field_t), intent(in) :: u_e, v_e, w_e
235 type(field_t), intent(inout) :: p_res
236 type(field_t), intent(in) :: f_x, f_y, f_z
237 type(coef_t), intent(inout) :: c_Xh
238 type(gs_t), intent(inout) :: gs_Xh
239 type(facet_normal_t), intent(in) :: bc_prs_surface
240 type(facet_normal_t), intent(in) :: bc_sym_surface
241 class(ax_t), intent(inout) :: Ax
242 real(kind=rp), intent(in) :: bd
243 real(kind=rp), intent(in) :: dt
244 type(field_t), intent(in) :: mu
245 type(field_t), intent(in) :: rho
246 type(c_ptr), intent(inout) :: event
247 real(kind=rp) :: dtbd
248 integer :: n, nelv, lxyz, gdim
249 integer :: i, e
250 ! Work arrays
251 type(field_t), pointer :: ta1, ta2, ta3, wa1, wa2, wa3, work1, work2
252 type(field_t), pointer :: s11, s22, s33, s12, s13, s23
253 integer :: temp_indices(14)
254
255 ! Work arrays
256 call neko_scratch_registry%request_field(ta1, temp_indices(1), .false.)
257 call neko_scratch_registry%request_field(ta2, temp_indices(2), .false.)
258 call neko_scratch_registry%request_field(ta3, temp_indices(3), .false.)
259 call neko_scratch_registry%request_field(wa1, temp_indices(4), .false.)
260 call neko_scratch_registry%request_field(wa2, temp_indices(5), .false.)
261 call neko_scratch_registry%request_field(wa3, temp_indices(6), .false.)
262 call neko_scratch_registry%request_field(work1, temp_indices(7), .false.)
263 call neko_scratch_registry%request_field(work2, temp_indices(8), .false.)
264
265 ! Stress tensor
266 call neko_scratch_registry%request_field(s11, temp_indices(9), .false.)
267 call neko_scratch_registry%request_field(s22, temp_indices(10), .false.)
268 call neko_scratch_registry%request_field(s33, temp_indices(11), .false.)
269 call neko_scratch_registry%request_field(s12, temp_indices(12), .false.)
270 call neko_scratch_registry%request_field(s13, temp_indices(13), .false.)
271 call neko_scratch_registry%request_field(s23, temp_indices(14), .false.)
272
273 n = c_xh%dof%size()
274 lxyz = c_xh%Xh%lxyz
275 nelv = c_xh%msh%nelv
276 gdim = c_xh%msh%gdim
277
278 call device_copy(c_xh%h1_d, rho%x_d, n)
279 call device_invcol1(c_xh%h1_d, n)
280 call device_rzero(c_xh%h2_d, n)
281 c_xh%ifh2 = .false.
282
283 ! mu times the double curl of the velocity
284 call curl(ta1, ta2, ta3, u_e, v_e, w_e, work1, work2, c_xh, event)
285 call curl(wa1, wa2, wa3, ta1, ta2, ta3, work1, work2, c_xh, event)
286
287 call device_col2(wa1%x_d, mu%x_d, n)
288 call device_col2(wa2%x_d, mu%x_d, n)
289 call device_col2(wa3%x_d, mu%x_d, n)
290
291
292 ! The strain rate tensor
293 call strain_rate(s11%x_d, s22%x_d, s33%x_d, s12%x_d, s13%x_d, s23%x_d, &
294 u_e%x_d, v_e%x_d, w_e%x_d, c_xh)
295
296
297 ! Gradient of viscosity * 2
298 !call opgrad(ta1%x, ta2%x, ta3%x, mu%x, c_Xh)
299 call dudxyz(ta1%x_d, mu%x_d, c_xh%drdx_d, c_xh%dsdx_d, c_xh%dtdx_d, c_xh)
300 call dudxyz(ta2%x_d, mu%x_d, c_xh%drdy_d, c_xh%dsdy_d, c_xh%dtdy_d, c_xh)
301 call dudxyz(ta3%x_d, mu%x_d, c_xh%drdz_d, c_xh%dsdz_d, c_xh%dtdz_d, c_xh)
302
303#ifdef HAVE_HIP
304 call pnpn_prs_stress_res_part1_hip(ta1%x_d, ta2%x_d, ta3%x_d, &
305 wa1%x_d, wa2%x_d, wa3%x_d, &
306 s11%x_d, s22%x_d, s33%x_d, s12%x_d, s13%x_d, s23%x_d, &
307 f_x%x_d, f_y%x_d, f_z%x_d, &
308 c_xh%B_d, c_xh%h1_d, rho%x_d, n)
309#elif HAVE_CUDA
310 call pnpn_prs_stress_res_part1_cuda(ta1%x_d, ta2%x_d, ta3%x_d, &
311 wa1%x_d, wa2%x_d, wa3%x_d, &
312 s11%x_d, s22%x_d, s33%x_d, s12%x_d, s13%x_d, s23%x_d, &
313 f_x%x_d, f_y%x_d, f_z%x_d, &
314 c_xh%B_d, c_xh%h1_d, rho%x_d, n)
315#elif HAVE_OPENCL
316 call pnpn_prs_stress_res_part1_opencl(ta1%x_d, ta2%x_d, ta3%x_d, &
317 wa1%x_d, wa2%x_d, wa3%x_d, &
318 s11%x_d, s22%x_d, s33%x_d, s12%x_d, s13%x_d, s23%x_d, &
319 f_x%x_d, f_y%x_d, f_z%x_d, &
320 c_xh%B_d, c_xh%h1_d, rho%x_d, n)
321#elif HAVE_METAL
322 call neko_error('Stress formulation prs res is not implemented for Metal')
323#endif
324
325 call rotate_cyc(ta1%x_d, ta2%x_d, ta3%x_d, 1, c_xh)
326 call gs_xh%op(ta1, gs_op_add)
327 call gs_xh%op(ta2, gs_op_add)
328 call gs_xh%op(ta3, gs_op_add)
329 call rotate_cyc(ta1%x_d, ta2%x_d, ta3%x_d, 0, c_xh)
330
331 call device_opcolv(ta1%x_d, ta2%x_d, ta3%x_d, c_xh%Binv_d, gdim, n)
332
333 ! Compute the components of the divergence of the rhs
334 call cdtp(wa1%x, ta1%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
335 call cdtp(wa2%x, ta2%x, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
336 call cdtp(wa3%x, ta3%x, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
337
338 ! The laplacian of the pressure
339 call ax%compute(p_res%x, p%x, c_xh, p%msh, p%Xh)
340
341#ifdef HAVE_HIP
342 call pnpn_prs_res_part2_hip(p_res%x_d, wa1%x_d, wa2%x_d, wa3%x_d, n)
343#elif HAVE_CUDA
344 call pnpn_prs_res_part2_cuda(p_res%x_d, wa1%x_d, wa2%x_d, wa3%x_d, n)
345#elif HAVE_OPENCL
346 call pnpn_prs_res_part2_opencl(p_res%x_d, wa1%x_d, wa2%x_d, wa3%x_d, n)
347#elif HAVE_METAL
348 call pnpn_prs_res_part2_metal(p_res%x_d, wa1%x_d, wa2%x_d, wa3%x_d, n)
349#endif
350
351 !
352 ! Surface velocity terms
353 !
354 call device_rzero(wa1%x_d, n)
355 call device_rzero(wa2%x_d, n)
356 call device_rzero(wa3%x_d, n)
357
358 call bc_sym_surface%apply_surfvec_dev(wa1%x_d, wa2%x_d, wa3%x_d, &
359 ta1%x_d , ta2%x_d, ta3%x_d)
360
361 dtbd = bd / dt
362 call device_rzero(ta1%x_d, n)
363 call device_rzero(ta2%x_d, n)
364 call device_rzero(ta3%x_d, n)
365
366 call bc_prs_surface%apply_surfvec_dev(ta1%x_d, ta2%x_d, ta3%x_d, &
367 u%x_D, v%x_d, w%x_d)
368
369#ifdef HAVE_HIP
370 call pnpn_prs_stress_res_part3_hip(p_res%x_d, ta1%x_d, ta2%x_d, ta3%x_d, &
371 wa1%x_d, wa2%x_d, wa3%x_d, dtbd, n)
372#elif HAVE_CUDA
373 call pnpn_prs_stress_res_part3_cuda(p_res%x_d, ta1%x_d, ta2%x_d, ta3%x_d, &
374 wa1%x_d, wa2%x_d, wa3%x_d, dtbd, n)
375#elif HAVE_OPENCL
376 call pnpn_prs_stress_res_part3_opencl(p_res%x_d, ta1%x_d, ta2%x_d, ta3%x_d,&
377 wa1%x_d, wa2%x_d, wa3%x_d, dtbd, n)
378#elif HAVE_METAL
379 call neko_error('Stress formulation prs ress is not implemented for Metal')
380#endif
381
382 call neko_scratch_registry%relinquish_field(temp_indices)
383
385
386 subroutine pnpn_vel_res_stress_device_compute(Ax, u, v, w, u_res, v_res, &
387 w_res, p, f_x, f_y, f_z, c_Xh, msh, Xh, mu, rho, bd, dt, n)
388 class(ax_t), intent(in) :: Ax
389 type(mesh_t), intent(inout) :: msh
390 type(space_t), intent(inout) :: Xh
391 type(field_t), intent(inout) :: p, u, v, w
392 type(field_t), intent(inout) :: u_res, v_res, w_res
393 type(field_t), intent(in) :: f_x, f_y, f_z
394 type(coef_t), intent(inout) :: c_Xh
395 type(field_t), intent(in) :: mu
396 type(field_t), intent(in) :: rho
397 real(kind=rp), intent(in) :: bd
398 real(kind=rp), intent(in) :: dt
399 real(kind=rp) :: bddt
400 integer :: temp_indices(3)
401 type(field_t), pointer :: ta1, ta2, ta3
402 integer, intent(in) :: n
403 integer :: i
404
405 call device_copy(c_xh%h1_d, mu%x_d, n)
406 call device_copy(c_xh%h2_d, rho%x_d, n)
407
408 bddt = bd / dt
409 call device_cmult(c_xh%h2_d, bddt, n)
410
411 c_xh%ifh2 = .true.
412
413 ! Viscous stresses
414 call ax%compute_vector(u_res%x, v_res%x, w_res%x, u%x, v%x, w%x, c_xh,&
415 msh, xh)
416
417 call neko_scratch_registry%request_field(ta1, temp_indices(1), .false.)
418 call neko_scratch_registry%request_field(ta2, temp_indices(2), .false.)
419 call neko_scratch_registry%request_field(ta3, temp_indices(3), .false.)
420
421 ! Pressure gradient
422 call opgrad(ta1%x, ta2%x, ta3%x, p%x, c_xh)
423
424 ! Sum all the terms
425#ifdef HAVE_HIP
426 call pnpn_vel_res_update_hip(u_res%x_d, v_res%x_d, w_res%x_d, &
427 ta1%x_d, ta2%x_d, ta3%x_d, f_x%x_d, f_y%x_d, f_z%x_d, n)
428#elif HAVE_CUDA
429 call pnpn_vel_res_update_cuda(u_res%x_d, v_res%x_d, w_res%x_d, &
430 ta1%x_d, ta2%x_d, ta3%x_d, f_x%x_d, f_y%x_d, f_z%x_d, n)
431#elif HAVE_OPENCL
432 call pnpn_vel_res_update_opencl(u_res%x_d, v_res%x_d, w_res%x_d, &
433 ta1%x_d, ta2%x_d, ta3%x_d, f_x%x_d, f_y%x_d, f_z%x_d, n)
434#elif HAVE_METAL
435 call pnpn_vel_res_update_metal(u_res%x_d, v_res%x_d, w_res%x_d, &
436 ta1%x_d, ta2%x_d, ta3%x_d, f_x%x_d, f_y%x_d, f_z%x_d, n)
437#endif
438
439 call neko_scratch_registry%relinquish_field(temp_indices)
440
442
443end module pnpn_res_stress_device
Compute derivative of a scalar field along a single direction.
Definition operators.f90:78
Apply cyclic boundary condition to a vector field.
Compute the strain rate tensor of a vector field.
Defines a Matrix-vector product.
Definition ax.f90:34
Coefficients.
Definition coef.f90:34
subroutine, public device_invcol1(a_d, n, strm)
Invert a vector .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_cmult(a_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
Dirichlet condition applied in the facet normal direction.
Defines a field.
Definition field.f90:34
Gather-scatter.
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public c_rp
Definition num_types.f90:13
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 (device version)
subroutine pnpn_vel_res_stress_device_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_stress_device_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
Utilities.
Definition utils.f90:35
void pnpn_prs_res_part2_opencl(void *p_res, void *wa1, void *wa2, void *wa3, int *n)
Definition pnpn_res.c:88
void pnpn_vel_res_update_opencl(void *u_res, void *v_res, void *w_res, void *ta1, void *ta2, void *ta3, void *f_u, void *f_v, void *f_w, int *n)
Definition pnpn_res.c:143
void pnpn_prs_res_part2_cuda(void *p_res, void *wa1, void *wa2, void *wa3, int *n)
Definition pnpn_res.cu:63
void pnpn_vel_res_update_cuda(void *u_res, void *v_res, void *w_res, void *ta1, void *ta2, void *ta3, void *f_u, void *f_v, void *f_w, int *n)
Definition pnpn_res.cu:92
void pnpn_prs_stress_res_part1_opencl(void *ta1, void *ta2, void *ta3, void *wa1, void *wa2, void *wa3, void *s11, void *s22, void *s33, void *s12, void *s13, void *s23, void *f_u, void *f_v, void *f_w, void *B, void *h1, void *rho, int *n)
void pnpn_prs_stress_res_part3_opencl(void *p_res, void *ta1, void *ta2, void *ta3, void *wa1, void *wa2, void *wa3, real *dtbd, int *n)
void pnpn_prs_stress_res_part1_cuda(void *ta1, void *ta2, void *ta3, void *wa1, void *wa2, void *wa3, void *s11, void *s22, void *s33, void *s12, void *s13, void *s23, void *f_u, void *f_v, void *f_w, void *B, void *h1, void *rho, int *n)
void pnpn_prs_stress_res_part3_cuda(void *p_res, void *ta1, void *ta2, void *ta3, void *wa1, void *wa2, void *wa3, real *dtbd, int *n)
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:63
Dirichlet condition in facet normal direction.
Device implementation of the pressure residual for the PnPn fluid with full viscous stress formulatio...
Device implementation of the velocity residual for the PnPn fluid with full viscous stress formulatio...
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