Neko  0.8.99
A portable framework for high-order spectral element flow simulations
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
9  use facet_normal, only : facet_normal_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
16  use device_mathops, only : device_opcolv
19  implicit none
20  private
21 
24  type, public, extends(pnpn_prs_res_t) :: pnpn_prs_res_stress_device_t
25  contains
26  procedure, nopass :: compute => pnpn_prs_res_stress_device_compute
28 
31  type, public, extends(pnpn_vel_res_t) :: pnpn_vel_res_stress_device_t
32  contains
33  procedure, nopass :: compute => pnpn_vel_res_stress_device_compute
35 
36 #ifdef HAVE_HIP
37  interface
38  subroutine pnpn_prs_res_part1_hip(ta1_d, ta2_d, ta3_d, &
39  wa1_d, wa2_d, wa3_d, f_u_d, f_v_d, f_w_d, &
40  B_d, h1_d, mu, rho, n) &
41  bind(c, name = 'pnpn_prs_res_part1_hip')
42  use, intrinsic :: iso_c_binding
43  import c_rp
44  implicit none
45  type(c_ptr), value :: ta1_d, ta2_d, ta3_d
46  type(c_ptr), value :: wa1_d, wa2_d, wa3_d
47  type(c_ptr), value :: f_u_d, f_v_d, f_w_d
48  type(c_ptr), value :: B_d, h1_d
49  real(c_rp) :: mu, rho
50  integer(c_int) :: n
51  end subroutine pnpn_prs_res_part1_hip
52  end interface
53 
54  interface
55  subroutine pnpn_prs_res_part2_hip(p_res_d, wa1_d, wa2_d, wa3_d, n) &
56  bind(c, name = 'pnpn_prs_res_part2_hip')
57  use, intrinsic :: iso_c_binding
58  implicit none
59  type(c_ptr), value :: p_res_d, wa1_d, wa2_d, wa3_d
60  integer(c_int) :: n
61  end subroutine pnpn_prs_res_part2_hip
62  end interface
63 
64  interface
65  subroutine pnpn_prs_res_part3_hip(p_res_d, ta1_d, ta2_d, ta3_d, dtbd, n) &
66  bind(c, name = 'pnpn_prs_res_part3_hip')
67  use, intrinsic :: iso_c_binding
68  import c_rp
69  implicit none
70  type(c_ptr), value :: p_res_d, ta1_d, ta2_d, ta3_d
71  real(c_rp) :: dtbd
72  integer(c_int) :: n
73  end subroutine pnpn_prs_res_part3_hip
74  end interface
75 
76  interface
77  subroutine pnpn_vel_res_update_hip(u_res_d, v_res_d, w_res_d, &
78  ta1_d, ta2_d, ta3_d, f_u_d, f_v_d, f_w_d, n) &
79  bind(c, name = 'pnpn_vel_res_update_hip')
80  use, intrinsic :: iso_c_binding
81  implicit none
82  type(c_ptr), value :: u_res_d, v_res_d, w_res_d
83  type(c_ptr), value :: ta1_d, ta2_d, ta3_d
84  type(c_ptr), value :: f_u_d, f_v_d, f_w_d
85  integer(c_int) :: n
86  end subroutine pnpn_vel_res_update_hip
87  end interface
88 #elif HAVE_CUDA
89  interface
90  subroutine pnpn_prs_stress_res_part1_cuda(ta1_d, ta2_d, ta3_d, &
91  wa1_d, wa2_d, wa3_d, f_u_d, f_v_d, f_w_d, &
92  B_d, h1_d, rho_d, n) &
93  bind(c, name = 'pnpn_prs_stress_res_part1_cuda')
94  use, intrinsic :: iso_c_binding
95  import c_rp
96  implicit none
97  type(c_ptr), value :: ta1_d, ta2_d, ta3_d
98  type(c_ptr), value :: wa1_d, wa2_d, wa3_d
99  type(c_ptr), value :: f_u_d, f_v_d, f_w_d
100  type(c_ptr), value :: B_d, h1_d, rho_d
101  integer(c_int) :: n
102  end subroutine pnpn_prs_stress_res_part1_cuda
103  end interface
104 
105  interface
106  subroutine pnpn_prs_res_part2_cuda(p_res_d, wa1_d, wa2_d, wa3_d, n) &
107  bind(c, name = 'pnpn_prs_res_part2_cuda')
108  use, intrinsic :: iso_c_binding
109  implicit none
110  type(c_ptr), value :: p_res_d, wa1_d, wa2_d, wa3_d
111  integer(c_int) :: n
112  end subroutine pnpn_prs_res_part2_cuda
113  end interface
114 
115  interface
116  subroutine pnpn_prs_stress_res_part3_cuda(p_res_d, ta1_d, ta2_d, ta3_d, &
117  wa1_d, wa2_d, wa3_d, dtbd, n) &
118  bind(c, name = 'pnpn_prs_stress_res_part3_cuda')
119  use, intrinsic :: iso_c_binding
120  import c_rp
121  implicit none
122  type(c_ptr), value :: p_res_d, ta1_d, ta2_d, ta3_d
123  type(c_ptr), value :: wa1_d, wa2_d, wa3_d
124  real(c_rp) :: dtbd
125  integer(c_int) :: n
126  end subroutine pnpn_prs_stress_res_part3_cuda
127  end interface
128 
129  interface
130  subroutine pnpn_vel_res_update_cuda(u_res_d, v_res_d, w_res_d, &
131  ta1_d, ta2_d, ta3_d, f_u_d, f_v_d, f_w_d, n) &
132  bind(c, name = 'pnpn_vel_res_update_cuda')
133  use, intrinsic :: iso_c_binding
134  implicit none
135  type(c_ptr), value :: u_res_d, v_res_d, w_res_d
136  type(c_ptr), value :: ta1_d, ta2_d, ta3_d
137  type(c_ptr), value :: f_u_d, f_v_d, f_w_d
138  integer(c_int) :: n
139  end subroutine pnpn_vel_res_update_cuda
140  end interface
141 #elif HAVE_OPENCL
142  interface
143  subroutine pnpn_prs_res_part1_opencl(ta1_d, ta2_d, ta3_d, &
144  wa1_d, wa2_d, wa3_d, f_u_d, f_v_d, f_w_d, &
145  B_d, h1_d, mu, rho, n) &
146  bind(c, name = 'pnpn_prs_res_part1_opencl')
147  use, intrinsic :: iso_c_binding
148  import c_rp
149  implicit none
150  type(c_ptr), value :: ta1_d, ta2_d, ta3_d
151  type(c_ptr), value :: wa1_d, wa2_d, wa3_d
152  type(c_ptr), value :: f_u_d, f_v_d, f_w_d
153  type(c_ptr), value :: B_d, h1_d
154  real(c_rp) :: mu, rho
155  integer(c_int) :: n
156  end subroutine pnpn_prs_res_part1_opencl
157  end interface
158 
159  interface
160  subroutine pnpn_prs_res_part2_opencl(p_res_d, wa1_d, wa2_d, wa3_d, n) &
161  bind(c, name = 'pnpn_prs_res_part2_opencl')
162  use, intrinsic :: iso_c_binding
163  implicit none
164  type(c_ptr), value :: p_res_d, wa1_d, wa2_d, wa3_d
165  integer(c_int) :: n
166  end subroutine pnpn_prs_res_part2_opencl
167  end interface
168 
169  interface
170  subroutine pnpn_prs_res_part3_opencl(p_res_d, ta1_d, ta2_d, ta3_d, &
171  wa1_d, wa2_d, wa3_d, dtbd, n) &
172  bind(c, name = 'pnpn_prs_res_part3_opencl')
173  use, intrinsic :: iso_c_binding
174  import c_rp
175  implicit none
176  type(c_ptr), value :: p_res_d, ta1_d, ta2_d, ta3_d
177  type(c_ptr), value :: wa1_d, wa2_d, wa3_d
178  real(c_rp) :: dtbd
179  integer(c_int) :: n
180  end subroutine pnpn_prs_res_part3_opencl
181  end interface
182 
183  interface
184  subroutine pnpn_vel_res_update_opencl(u_res_d, v_res_d, w_res_d, &
185  ta1_d, ta2_d, ta3_d, f_u_d, f_v_d, f_w_d, n) &
186  bind(c, name = 'pnpn_vel_res_update_opencl')
187  use, intrinsic :: iso_c_binding
188  implicit none
189  type(c_ptr), value :: u_res_d, v_res_d, w_res_d
190  type(c_ptr), value :: ta1_d, ta2_d, ta3_d
191  type(c_ptr), value :: f_u_d, f_v_d, f_w_d
192  integer(c_int) :: n
193  end subroutine pnpn_vel_res_update_opencl
194  end interface
195 #endif
196 
197 contains
198 
199  subroutine pnpn_prs_res_stress_device_compute(p, p_res, u, v, w, u_e, v_e,&
200  w_e, f_x, f_y, f_z, c_Xh, gs_Xh, bc_prs_surface, bc_sym_surface, Ax, bd,&
201  dt, mu, rho)
202  type(field_t), intent(inout) :: p, u, v, w
203  type(field_t), intent(inout) :: u_e, v_e, w_e
204  type(field_t), intent(inout) :: p_res
205  type(field_t), intent(inout) :: f_x, f_y, f_z
206  type(coef_t), intent(inout) :: c_Xh
207  type(gs_t), intent(inout) :: gs_Xh
208  type(facet_normal_t), intent(inout) :: bc_prs_surface
209  type(facet_normal_t), intent(inout) :: bc_sym_surface
210  class(ax_t), intent(inout) :: Ax
211  real(kind=rp), intent(inout) :: bd
212  real(kind=rp), intent(in) :: dt
213  type(field_t), intent(in) :: mu
214  type(field_t), intent(in) :: rho
215  real(kind=rp) :: dtbd
216  integer :: n, nelv, lxyz, gdim
217  integer :: i, e
218  ! Work arrays
219  type(field_t), pointer :: ta1, ta2, ta3, wa1, wa2, wa3, work1, work2, work3
220  type(field_t), pointer :: s11, s22, s33, s12, s13, s23
221  integer :: temp_indices(15)
222 
223  ! Work arrays
224  call neko_scratch_registry%request_field(ta1, temp_indices(1))
225  call neko_scratch_registry%request_field(ta2, temp_indices(2))
226  call neko_scratch_registry%request_field(ta3, temp_indices(3))
227  call neko_scratch_registry%request_field(wa1, temp_indices(4))
228  call neko_scratch_registry%request_field(wa2, temp_indices(5))
229  call neko_scratch_registry%request_field(wa3, temp_indices(6))
230  call neko_scratch_registry%request_field(work1, temp_indices(7))
231  call neko_scratch_registry%request_field(work2, temp_indices(8))
232  call neko_scratch_registry%request_field(work3, temp_indices(9))
233 
234  ! Stress tensor
235  call neko_scratch_registry%request_field(s11, temp_indices(10))
236  call neko_scratch_registry%request_field(s22, temp_indices(11))
237  call neko_scratch_registry%request_field(s33, temp_indices(12))
238  call neko_scratch_registry%request_field(s12, temp_indices(13))
239  call neko_scratch_registry%request_field(s13, temp_indices(14))
240  call neko_scratch_registry%request_field(s23, temp_indices(15))
241 
242  n = c_xh%dof%size()
243  lxyz = c_xh%Xh%lxyz
244  nelv = c_xh%msh%nelv
245  gdim = c_xh%msh%gdim
246 
247  call device_copy(c_xh%h1_d, rho%x_d, n)
248  call device_invcol1(c_xh%h1_d, n)
249  call device_rzero(c_xh%h2_d, n)
250  c_xh%ifh2 = .false.
251 
252  ! mu times the double curl of the velocity
253  call curl(ta1, ta2, ta3, u_e, v_e, w_e, work1, work2, c_xh)
254  call curl(wa1, wa2, wa3, ta1, ta2, ta3, work1, work2, c_xh)
255 
256  call device_col2(wa1%x_d, mu%x_d, n)
257  call device_col2(wa2%x_d, mu%x_d, n)
258  call device_col2(wa3%x_d, mu%x_d, n)
259 
260 
261  ! The strain rate tensor
262  call strain_rate(s11%x, s22%x, s33%x, s12%x, s13%x, s23%x, &
263  u_e, v_e, w_e, c_xh)
264 
265 
266  ! Gradient of viscosity * 2
267  !call opgrad(ta1%x, ta2%x, ta3%x, mu%x, c_Xh)
268  call dudxyz(ta1%x, mu%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
269  call dudxyz(ta2%x, mu%x, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
270  call dudxyz(ta3%x, mu%x, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
271 
272  call device_cmult(ta1%x_d, 2.0_rp, n)
273  call device_cmult(ta2%x_d, 2.0_rp, n)
274  call device_cmult(ta3%x_d, 2.0_rp, n)
275 
276  ! S^T grad \mu
277  call device_vdot3 (work1%x_d, ta1%x_d, ta2%x_d, ta3%x_d, &
278  s11%x_d, s12%x_d, s13%x_d, n)
279 
280  call device_vdot3 (work2%x_d, ta1%x_d, ta2%x_d, ta3%x_d, &
281  s12%x_d, s22%x_d, s23%x_d, lxyz)
282 
283  call device_vdot3 (work3%x_d, ta1%x_d, ta2%x_d, ta3%x_d, &
284  s13%x_d, s23%x_d, s33%x_d, lxyz)
285 
286  ! Subtract the two terms of the viscous stress to get
287  ! \nabla x \nabla u - S^T \nabla \mu
288  ! The sign is consitent with the fact that we subtract the term
289  ! below.
290  call device_sub2(wa1%x_d, work1%x_d, n)
291  call device_sub2(wa2%x_d, work2%x_d, n)
292  call device_sub2(wa3%x_d, work3%x_d, n)
293 
294 #if HAVE_CUDA
295  call pnpn_prs_stress_res_part1_cuda(ta1%x_d, ta2%x_d, ta3%x_d, &
296  wa1%x_d, wa2%x_d, wa3%x_d, f_x%x_d, f_y%x_d, f_z%x_d, &
297  c_xh%B_d, c_xh%h1_d, rho%x_d, n)
298 #else
299  call neko_error('No device backend configured')
300 #endif
301 
302  call gs_xh%op(ta1, gs_op_add)
303  call gs_xh%op(ta2, gs_op_add)
304  call gs_xh%op(ta3, gs_op_add)
305 
306  call device_opcolv(ta1%x_d, ta2%x_d, ta3%x_d, c_xh%Binv_d, gdim, n)
307 
308  ! Compute the components of the divergence of the rhs
309  call cdtp(wa1%x, ta1%x, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
310  call cdtp(wa2%x, ta2%x, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
311  call cdtp(wa3%x, ta3%x, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
312 
313  ! The laplacian of the pressure
314  call ax%compute(p_res%x, p%x, c_xh, p%msh, p%Xh)
315 
316 #ifdef HAVE_HIP
317  call pnpn_prs_res_part2_hip(p_res%x_d, wa1%x_d, wa2%x_d, wa3%x_d, n)
318 #elif HAVE_CUDA
319  call pnpn_prs_res_part2_cuda(p_res%x_d, wa1%x_d, wa2%x_d, wa3%x_d, n)
320 #elif HAVE_OPENCL
321  call pnpn_prs_res_part2_opencl(p_res%x_d, wa1%x_d, wa2%x_d, wa3%x_d, n)
322 #endif
323 
324  !
325  ! Surface velocity terms
326  !
327  call device_rzero(wa1%x_d, n)
328  call device_rzero(wa2%x_d, n)
329  call device_rzero(wa3%x_d, n)
330 
331  call bc_sym_surface%apply_surfvec_dev(wa1%x_d, wa2%x_d, wa3%x_d, &
332  ta1%x_d , ta2%x_d, ta3%x_d)
333 
334  dtbd = bd / dt
335  call device_rzero(ta1%x_d, n)
336  call device_rzero(ta2%x_d, n)
337  call device_rzero(ta3%x_d, n)
338 
339  call bc_prs_surface%apply_surfvec_dev(ta1%x_d, ta2%x_d, ta3%x_d, &
340  u%x_D, v%x_d, w%x_d)
341 
342 #if HAVE_CUDA
343  call pnpn_prs_stress_res_part3_cuda(p_res%x_d, ta1%x_d, ta2%x_d, ta3%x_d, &
344  wa1%x_d, wa2%x_d, wa3%x_d, dtbd, n)
345 #else
346  call neko_error('No device backend configured')
347 #endif
348 
349  call neko_scratch_registry%relinquish_field(temp_indices)
350 
352 
353  subroutine pnpn_vel_res_stress_device_compute(Ax, u, v, w, u_res, v_res, &
354  w_res, p, f_x, f_y, f_z, c_Xh, msh, Xh, mu, rho, bd, dt, n)
355  class(ax_t), intent(in) :: Ax
356  type(mesh_t), intent(inout) :: msh
357  type(space_t), intent(inout) :: Xh
358  type(field_t), intent(inout) :: p, u, v, w
359  type(field_t), intent(inout) :: u_res, v_res, w_res
360  type(field_t), intent(inout) :: f_x, f_y, f_z
361  type(coef_t), intent(inout) :: c_Xh
362  type(field_t), intent(in) :: mu
363  type(field_t), intent(in) :: rho
364  real(kind=rp), intent(in) :: bd
365  real(kind=rp), intent(in) :: dt
366  real(kind=rp) :: bddt
367  integer :: temp_indices(3)
368  type(field_t), pointer :: ta1, ta2, ta3
369  integer, intent(in) :: n
370  integer :: i
371 
372  call device_copy(c_xh%h1_d, mu%x_d, n)
373  call device_copy(c_xh%h2_d, rho%x_d, n)
374 
375  bddt = bd / dt
376  call device_cmult(c_xh%h2_d, bddt, n)
377 
378  c_xh%ifh2 = .true.
379 
380  ! Viscous stresses
381  call ax%compute_vector(u_res%x, v_res%x, w_res%x, u%x, v%x, w%x, c_xh,&
382  msh, xh)
383 
384  call neko_scratch_registry%request_field(ta1, temp_indices(1))
385  call neko_scratch_registry%request_field(ta2, temp_indices(2))
386  call neko_scratch_registry%request_field(ta3, temp_indices(3))
387 
388  ! Pressure gradient
389  call opgrad(ta1%x, ta2%x, ta3%x, p%x, c_xh)
390 
391  ! Sum all the terms
392 #ifdef HAVE_HIP
393  call pnpn_vel_res_update_hip(u_res%x_d, v_res%x_d, w_res%x_d, &
394  ta1%x_d, ta2%x_d, ta3%x_d, f_x%x_d, f_y%x_d, f_z%x_d, n)
395 #elif HAVE_CUDA
396  call pnpn_vel_res_update_cuda(u_res%x_d, v_res%x_d, w_res%x_d, &
397  ta1%x_d, ta2%x_d, ta3%x_d, f_x%x_d, f_y%x_d, f_z%x_d, n)
398 #elif HAVE_OPENCL
399  call pnpn_vel_res_update_opencl(u_res%x_d, v_res%x_d, w_res%x_d, &
400  ta1%x_d, ta2%x_d, ta3%x_d, f_x%x_d, f_y%x_d, f_z%x_d, n)
401 #endif
402 
403  call neko_scratch_registry%relinquish_field(temp_indices)
404 
406 
407 end module pnpn_res_stress_device
Defines a Matrix-vector product.
Definition: ax.f90:34
Coefficients.
Definition: coef.f90:34
subroutine, public device_col2(a_d, b_d, n)
Vector multiplication .
subroutine, public device_rzero(a_d, n)
Zero a real vector.
subroutine, public device_invcol1(a_d, n)
Invert a vector .
subroutine, public device_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n)
Compute a dot product (3-d version) assuming vector components etc.
subroutine, public device_cmult(a_d, c, n)
Multiplication by constant c .
subroutine, public device_copy(a_d, b_d, n)
Copy a vector .
subroutine, public device_sub2(a_d, b_d, n)
Vector substraction .
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.
Definition: operators.f90:169
subroutine, public curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
Definition: operators.f90:360
subroutine, public strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
Compute the strain rate tensor, i.e 0.5 * du_i/dx_j + du_j/dx_i.
Definition: operators.f90:428
subroutine, public cdtp(dtx, x, dr, ds, dt, coef, es, ee)
Apply D^T to a scalar field, where D is the derivative matrix.
Definition: operators.f90:228
subroutine, public dudxyz(du, u, dr, ds, dt, coef)
Compute derivative of a scalar field along a single direction.
Definition: operators.f90:76
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)
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
Utilities.
Definition: utils.f90:35
void pnpn_prs_res_part3_opencl(void *p_res, void *ta1, void *ta2, void *ta3, real *dtbd, int *n)
Definition: pnpn_res.c:108
void pnpn_prs_res_part1_opencl(void *ta1, void *ta2, void *ta3, void *wa1, void *wa2, void *wa3, void *f_u, void *f_v, void *f_w, void *B, void *h1, real *mu, real *rho, int *n)
Definition: pnpn_res.c:44
void pnpn_prs_res_part2_opencl(void *p_res, void *wa1, void *wa2, void *wa3, int *n)
Definition: pnpn_res.c:82
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:135
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_cuda(void *ta1, void *ta2, void *ta3, void *wa1, void *wa2, void *wa3, 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:55
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:47
Abstract type to compute velocity residual.
Definition: pnpn_res.f90:53
The function space for the SEM solution fields.
Definition: space.f90:62