Skip to content

Commit aa9d6fa

Browse files
committed
Add visc_deriv subroutine to m_viscous
1 parent 31420aa commit aa9d6fa

File tree

2 files changed

+110
-105
lines changed

2 files changed

+110
-105
lines changed

src/simulation/m_rhs.f90

Lines changed: 6 additions & 105 deletions
Original file line numberDiff line numberDiff line change
@@ -896,19 +896,22 @@ subroutine s_compute_rhs(q_cons_vf, q_prim_vf, rhs_vf, t_step) ! -------
896896
dq_prim_dx_qp%vf(iv%beg:iv%end), &
897897
dqL_rsx_vf, dqL_rsy_vf, dqL_rsz_vf, &
898898
dqR_rsx_vf, dqR_rsy_vf, dqR_rsz_vf, &
899-
id, dqL_prim_dx_n(id)%vf(iv%beg:iv%end), dqR_prim_dx_n(id)%vf(iv%beg:iv%end))
899+
id, dqL_prim_dx_n(id)%vf(iv%beg:iv%end), dqR_prim_dx_n(id)%vf(iv%beg:iv%end), &
900+
ix, iy, iz)
900901
if (n > 0) then
901902
call s_reconstruct_cell_boundary_values_visc_deriv( &
902903
dq_prim_dy_qp%vf(iv%beg:iv%end), &
903904
dqL_rsx_vf, dqL_rsy_vf, dqL_rsz_vf, &
904905
dqR_rsx_vf, dqR_rsy_vf, dqR_rsz_vf, &
905-
id, dqL_prim_dy_n(id)%vf(iv%beg:iv%end), dqR_prim_dy_n(id)%vf(iv%beg:iv%end))
906+
id, dqL_prim_dy_n(id)%vf(iv%beg:iv%end), dqR_prim_dy_n(id)%vf(iv%beg:iv%end), &
907+
ix, iy, iz)
906908
if (p > 0) then
907909
call s_reconstruct_cell_boundary_values_visc_deriv( &
908910
dq_prim_dz_qp%vf(iv%beg:iv%end), &
909911
dqL_rsx_vf, dqL_rsy_vf, dqL_rsz_vf, &
910912
dqR_rsx_vf, dqR_rsy_vf, dqR_rsz_vf, &
911-
id, dqL_prim_dz_n(id)%vf(iv%beg:iv%end), dqR_prim_dz_n(id)%vf(iv%beg:iv%end))
913+
id, dqL_prim_dz_n(id)%vf(iv%beg:iv%end), dqR_prim_dz_n(id)%vf(iv%beg:iv%end), &
914+
ix, iy, iz)
912915
end if
913916
end if
914917
end if
@@ -2539,108 +2542,6 @@ subroutine s_reconstruct_cell_boundary_values(v_vf, vL_x, vL_y, vL_z, vR_x, vR_y
25392542
! ==================================================================
25402543
end subroutine s_reconstruct_cell_boundary_values ! --------------------
25412544

2542-
subroutine s_reconstruct_cell_boundary_values_visc_deriv(v_vf, vL_x, vL_y, vL_z, vR_x, vR_y, vR_z, & ! -
2543-
norm_dir, vL_prim_vf, vR_prim_vf)
2544-
2545-
type(scalar_field), dimension(iv%beg:iv%end), intent(IN) :: v_vf
2546-
type(scalar_field), dimension(iv%beg:iv%end), intent(INOUT) :: vL_prim_vf, vR_prim_vf
2547-
2548-
real(kind(0d0)), dimension(startx:, starty:, startz:, iv%beg:), intent(INOUT) :: vL_x, vL_y, vL_z, vR_x, vR_y, vR_z
2549-
2550-
integer, intent(IN) :: norm_dir
2551-
2552-
integer :: weno_dir !< Coordinate direction of the WENO reconstruction
2553-
2554-
integer :: i, j, k, l
2555-
! Reconstruction in s1-direction ===================================
2556-
2557-
if (norm_dir == 1) then
2558-
is1 = ix; is2 = iy; is3 = iz
2559-
weno_dir = 1; is1%beg = is1%beg + weno_polyn
2560-
is1%end = is1%end - weno_polyn
2561-
2562-
elseif (norm_dir == 2) then
2563-
is1 = iy; is2 = ix; is3 = iz
2564-
weno_dir = 2; is1%beg = is1%beg + weno_polyn
2565-
is1%end = is1%end - weno_polyn
2566-
2567-
else
2568-
is1 = iz; is2 = iy; is3 = ix
2569-
weno_dir = 3; is1%beg = is1%beg + weno_polyn
2570-
is1%end = is1%end - weno_polyn
2571-
2572-
end if
2573-
2574-
!$acc update device(is1, is2, is3, iv)
2575-
2576-
if (n > 0) then
2577-
if (p > 0) then
2578-
2579-
call s_weno(v_vf(iv%beg:iv%end), &
2580-
vL_x(:, :, :, iv%beg:iv%end), vL_y(:, :, :, iv%beg:iv%end), vL_z(:, :, :, iv%beg:iv%end), vR_x(:, :, :, iv%beg:iv%end), vR_y(:, :, :, iv%beg:iv%end), vR_z(:, :, :, iv%beg:iv%end), &
2581-
norm_dir, weno_dir, &
2582-
is1, is2, is3)
2583-
else
2584-
call s_weno(v_vf(iv%beg:iv%end), &
2585-
vL_x(:, :, :, iv%beg:iv%end), vL_y(:, :, :, iv%beg:iv%end), vL_z(:, :, :, :), vR_x(:, :, :, iv%beg:iv%end), vR_y(:, :, :, iv%beg:iv%end), vR_z(:, :, :, :), &
2586-
norm_dir, weno_dir, &
2587-
is1, is2, is3)
2588-
end if
2589-
else
2590-
2591-
call s_weno(v_vf(iv%beg:iv%end), &
2592-
vL_x(:, :, :, iv%beg:iv%end), vL_y(:, :, :, :), vL_z(:, :, :, :), vR_x(:, :, :, iv%beg:iv%end), vR_y(:, :, :, :), vR_z(:, :, :, :), &
2593-
norm_dir, weno_dir, &
2594-
is1, is2, is3)
2595-
end if
2596-
2597-
if (any(Re_size > 0)) then
2598-
if (weno_Re_flux) then
2599-
if (norm_dir == 2) then
2600-
!$acc parallel loop collapse(4) gang vector default(present)
2601-
do i = iv%beg, iv%end
2602-
do l = is3%beg, is3%end
2603-
do j = is1%beg, is1%end
2604-
do k = is2%beg, is2%end
2605-
vL_prim_vf(i)%sf(k, j, l) = vL_y(j, k, l, i)
2606-
vR_prim_vf(i)%sf(k, j, l) = vR_y(j, k, l, i)
2607-
end do
2608-
end do
2609-
end do
2610-
end do
2611-
elseif (norm_dir == 3) then
2612-
!$acc parallel loop collapse(4) gang vector default(present)
2613-
do i = iv%beg, iv%end
2614-
do j = is1%beg, is1%end
2615-
do k = is2%beg, is2%end
2616-
do l = is3%beg, is3%end
2617-
vL_prim_vf(i)%sf(l, k, j) = vL_z(j, k, l, i)
2618-
vR_prim_vf(i)%sf(l, k, j) = vR_z(j, k, l, i)
2619-
end do
2620-
end do
2621-
end do
2622-
end do
2623-
elseif (norm_dir == 1) then
2624-
!$acc parallel loop collapse(4) gang vector default(present)
2625-
do i = iv%beg, iv%end
2626-
do l = is3%beg, is3%end
2627-
do k = is2%beg, is2%end
2628-
do j = is1%beg, is1%end
2629-
vL_prim_vf(i)%sf(j, k, l) = vL_x(j, k, l, i)
2630-
vR_prim_vf(i)%sf(j, k, l) = vR_x(j, k, l, i)
2631-
end do
2632-
end do
2633-
end do
2634-
end do
2635-
end if
2636-
end if
2637-
end if
2638-
2639-
! ==================================================================
2640-
2641-
end subroutine s_reconstruct_cell_boundary_values_visc_deriv ! --------------------
2642-
2643-
26442545
!> Module deallocation and/or disassociation procedures
26452546
subroutine s_finalize_rhs_module() ! -----------------------------------
26462547

src/simulation/m_viscous.f90

Lines changed: 104 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -18,6 +18,7 @@ module m_viscous
1818
private; public s_get_viscous, &
1919
s_compute_viscous_stress_tensor, &
2020
s_initialize_viscous_module, &
21+
s_reconstruct_cell_boundary_values_visc_deriv, &
2122
s_finalize_viscous_module
2223

2324
type(int_bounds_info) :: iv
@@ -1359,6 +1360,109 @@ subroutine s_reconstruct_cell_boundary_values_visc(v_vf, vL_x, vL_y, vL_z, vR_x,
13591360

13601361
end subroutine s_reconstruct_cell_boundary_values_visc ! --------------------
13611362

1363+
subroutine s_reconstruct_cell_boundary_values_visc_deriv(v_vf, vL_x, vL_y, vL_z, vR_x, vR_y, vR_z, & ! -
1364+
norm_dir, vL_prim_vf, vR_prim_vf, ix, iy, iz)
1365+
1366+
type(scalar_field), dimension(iv%beg:iv%end), intent(IN) :: v_vf
1367+
type(scalar_field), dimension(iv%beg:iv%end), intent(INOUT) :: vL_prim_vf, vR_prim_vf
1368+
1369+
type(int_bounds_info) :: ix, iy, iz
1370+
1371+
real(kind(0d0)), dimension(startx:, starty:, startz:, iv%beg:), intent(INOUT) :: vL_x, vL_y, vL_z, vR_x, vR_y, vR_z
1372+
1373+
integer, intent(IN) :: norm_dir
1374+
1375+
integer :: weno_dir !< Coordinate direction of the WENO reconstruction
1376+
1377+
integer :: i, j, k, l
1378+
! Reconstruction in s1-direction ===================================
1379+
1380+
if (norm_dir == 1) then
1381+
is1 = ix; is2 = iy; is3 = iz
1382+
weno_dir = 1; is1%beg = is1%beg + weno_polyn
1383+
is1%end = is1%end - weno_polyn
1384+
1385+
elseif (norm_dir == 2) then
1386+
is1 = iy; is2 = ix; is3 = iz
1387+
weno_dir = 2; is1%beg = is1%beg + weno_polyn
1388+
is1%end = is1%end - weno_polyn
1389+
1390+
else
1391+
is1 = iz; is2 = iy; is3 = ix
1392+
weno_dir = 3; is1%beg = is1%beg + weno_polyn
1393+
is1%end = is1%end - weno_polyn
1394+
1395+
end if
1396+
1397+
!$acc update device(is1, is2, is3, iv)
1398+
1399+
if (n > 0) then
1400+
if (p > 0) then
1401+
1402+
call s_weno(v_vf(iv%beg:iv%end), &
1403+
vL_x(:, :, :, iv%beg:iv%end), vL_y(:, :, :, iv%beg:iv%end), vL_z(:, :, :, iv%beg:iv%end), vR_x(:, :, :, iv%beg:iv%end), vR_y(:, :, :, iv%beg:iv%end), vR_z(:, :, :, iv%beg:iv%end), &
1404+
norm_dir, weno_dir, &
1405+
is1, is2, is3)
1406+
else
1407+
call s_weno(v_vf(iv%beg:iv%end), &
1408+
vL_x(:, :, :, iv%beg:iv%end), vL_y(:, :, :, iv%beg:iv%end), vL_z(:, :, :, :), vR_x(:, :, :, iv%beg:iv%end), vR_y(:, :, :, iv%beg:iv%end), vR_z(:, :, :, :), &
1409+
norm_dir, weno_dir, &
1410+
is1, is2, is3)
1411+
end if
1412+
else
1413+
1414+
call s_weno(v_vf(iv%beg:iv%end), &
1415+
vL_x(:, :, :, iv%beg:iv%end), vL_y(:, :, :, :), vL_z(:, :, :, :), vR_x(:, :, :, iv%beg:iv%end), vR_y(:, :, :, :), vR_z(:, :, :, :), &
1416+
norm_dir, weno_dir, &
1417+
is1, is2, is3)
1418+
end if
1419+
1420+
if (any(Re_size > 0)) then
1421+
if (weno_Re_flux) then
1422+
if (norm_dir == 2) then
1423+
!$acc parallel loop collapse(4) gang vector default(present)
1424+
do i = iv%beg, iv%end
1425+
do l = is3%beg, is3%end
1426+
do j = is1%beg, is1%end
1427+
do k = is2%beg, is2%end
1428+
vL_prim_vf(i)%sf(k, j, l) = vL_y(j, k, l, i)
1429+
vR_prim_vf(i)%sf(k, j, l) = vR_y(j, k, l, i)
1430+
end do
1431+
end do
1432+
end do
1433+
end do
1434+
elseif (norm_dir == 3) then
1435+
!$acc parallel loop collapse(4) gang vector default(present)
1436+
do i = iv%beg, iv%end
1437+
do j = is1%beg, is1%end
1438+
do k = is2%beg, is2%end
1439+
do l = is3%beg, is3%end
1440+
vL_prim_vf(i)%sf(l, k, j) = vL_z(j, k, l, i)
1441+
vR_prim_vf(i)%sf(l, k, j) = vR_z(j, k, l, i)
1442+
end do
1443+
end do
1444+
end do
1445+
end do
1446+
elseif (norm_dir == 1) then
1447+
!$acc parallel loop collapse(4) gang vector default(present)
1448+
do i = iv%beg, iv%end
1449+
do l = is3%beg, is3%end
1450+
do k = is2%beg, is2%end
1451+
do j = is1%beg, is1%end
1452+
vL_prim_vf(i)%sf(j, k, l) = vL_x(j, k, l, i)
1453+
vR_prim_vf(i)%sf(j, k, l) = vR_x(j, k, l, i)
1454+
end do
1455+
end do
1456+
end do
1457+
end do
1458+
end if
1459+
end if
1460+
end if
1461+
1462+
! ==================================================================
1463+
1464+
end subroutine s_reconstruct_cell_boundary_values_visc_deriv ! --------------------
1465+
13621466
subroutine s_finalize_viscous_module()
13631467
deallocate (Res)
13641468
end subroutine s_finalize_viscous_module

0 commit comments

Comments
 (0)