Skip to content

Commit c97c5be

Browse files
committed
refactor: consolidate MUSCL bounds triplet into muscl_bounds_t derived type
1 parent a0d3cce commit c97c5be

1 file changed

Lines changed: 47 additions & 42 deletions

File tree

src/simulation/m_muscl.fpp

Lines changed: 47 additions & 42 deletions
Original file line numberDiff line numberDiff line change
@@ -21,11 +21,15 @@ module m_muscl
2121

2222
private; public :: s_initialize_muscl_module, s_muscl, s_finalize_muscl_module
2323

24+
type, private :: muscl_bounds_t
25+
type(int_bounds_info) :: x, y, z
26+
end type muscl_bounds_t
27+
2428
integer :: v_size
2529
$:GPU_DECLARE(create='[v_size]')
2630

27-
type(int_bounds_info) :: is1_muscl, is2_muscl, is3_muscl
28-
$:GPU_DECLARE(create='[is1_muscl, is2_muscl, is3_muscl]')
31+
type(muscl_bounds_t) :: is_muscl
32+
$:GPU_DECLARE(create='[is_muscl]')
2933

3034
!> @name The cell-average variables that will be MUSCL-reconstructed, unpacked into an array for performance
3135
!> @{
@@ -39,73 +43,74 @@ contains
3943
subroutine s_initialize_muscl_module()
4044

4145
! Initializing in x-direction
42-
is1_muscl%beg = -buff_size; is1_muscl%end = m - is1_muscl%beg
46+
is_muscl%x%beg = -buff_size; is_muscl%x%end = m - is_muscl%x%beg
4347
if (n == 0) then
44-
is2_muscl%beg = 0
48+
is_muscl%y%beg = 0
4549
else
46-
is2_muscl%beg = -buff_size
50+
is_muscl%y%beg = -buff_size
4751
end if
4852

49-
is2_muscl%end = n - is2_muscl%beg
53+
is_muscl%y%end = n - is_muscl%y%beg
5054

5155
if (p == 0) then
52-
is3_muscl%beg = 0
56+
is_muscl%z%beg = 0
5357
else
54-
is3_muscl%beg = -buff_size
58+
is_muscl%z%beg = -buff_size
5559
end if
5660

57-
is3_muscl%end = p - is3_muscl%beg
61+
is_muscl%z%end = p - is_muscl%z%beg
5862

59-
@:ALLOCATE(v_rs_ws_muscl(is1_muscl%beg:is1_muscl%end, is2_muscl%beg:is2_muscl%end, is3_muscl%beg:is3_muscl%end, 1:sys_size))
63+
@:ALLOCATE(v_rs_ws_muscl(is_muscl%x%beg:is_muscl%x%end, is_muscl%y%beg:is_muscl%y%end, is_muscl%z%beg:is_muscl%z%end, &
64+
& 1:sys_size))
6065

6166
if (n == 0) return
6267

6368
! initializing in y-direction
64-
is2_muscl%beg = -buff_size; is2_muscl%end = n - is2_muscl%beg
65-
is1_muscl%beg = -buff_size; is1_muscl%end = m - is1_muscl%beg
69+
is_muscl%y%beg = -buff_size; is_muscl%y%end = n - is_muscl%y%beg
70+
is_muscl%x%beg = -buff_size; is_muscl%x%end = m - is_muscl%x%beg
6671

6772
if (p == 0) then
68-
is3_muscl%beg = 0
73+
is_muscl%z%beg = 0
6974
else
70-
is3_muscl%beg = -buff_size
75+
is_muscl%z%beg = -buff_size
7176
end if
7277

73-
is3_muscl%end = p - is3_muscl%beg
78+
is_muscl%z%end = p - is_muscl%z%beg
7479

7580
if (p == 0) return
7681

7782
! initializing in z-direction
78-
is2_muscl%beg = -buff_size; is2_muscl%end = n - is2_muscl%beg
79-
is1_muscl%beg = -buff_size; is1_muscl%end = m - is1_muscl%beg
80-
is3_muscl%beg = -buff_size; is3_muscl%end = p - is3_muscl%beg
83+
is_muscl%y%beg = -buff_size; is_muscl%y%end = n - is_muscl%y%beg
84+
is_muscl%x%beg = -buff_size; is_muscl%x%end = m - is_muscl%x%beg
85+
is_muscl%z%beg = -buff_size; is_muscl%z%end = p - is_muscl%z%beg
8186

8287
end subroutine s_initialize_muscl_module
8388

8489
!> Perform MUSCL reconstruction of left and right cell-boundary values from cell-averaged variables
85-
subroutine s_muscl(v_vf, vL_rs_vf_x, vR_rs_vf_x, muscl_dir, is1_muscl_d, &
90+
subroutine s_muscl(v_vf, vL_rs_vf_x, vR_rs_vf_x, muscl_dir, is1_d, &
8691

87-
& is2_muscl_d, is3_muscl_d)
92+
& is2_d, is3_d)
8893

89-
type(scalar_field), dimension(1:), intent(in) :: v_vf
94+
type(scalar_field), dimension(1:), intent(in) :: v_vf
9095
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) :: vL_rs_vf_x, vR_rs_vf_x
91-
integer, intent(in) :: muscl_dir
92-
type(int_bounds_info), intent(in) :: is1_muscl_d, is2_muscl_d, is3_muscl_d
93-
integer :: j, k, l, i
94-
real(wp) :: slopeL, slopeR, slope
96+
integer, intent(in) :: muscl_dir
97+
type(int_bounds_info), intent(in) :: is1_d, is2_d, is3_d
98+
integer :: j, k, l, i
99+
real(wp) :: slopeL, slopeR, slope
95100

96-
is1_muscl = is1_muscl_d
97-
is2_muscl = is2_muscl_d
98-
is3_muscl = is3_muscl_d
101+
is_muscl%x = is1_d
102+
is_muscl%y = is2_d
103+
is_muscl%z = is3_d
99104

100-
$:GPU_UPDATE(device='[is1_muscl, is2_muscl, is3_muscl]')
105+
$:GPU_UPDATE(device='[is_muscl]')
101106

102107
if (muscl_order == 1) then
103108
if (muscl_dir == 1) then
104109
$:GPU_PARALLEL_LOOP(collapse=4)
105110
do i = 1, ubound(v_vf, 1)
106-
do l = is3_muscl%beg, is3_muscl%end
107-
do k = is2_muscl%beg, is2_muscl%end
108-
do j = is1_muscl%beg, is1_muscl%end
111+
do l = is_muscl%z%beg, is_muscl%z%end
112+
do k = is_muscl%y%beg, is_muscl%y%end
113+
do j = is_muscl%x%beg, is_muscl%x%end
109114
vL_rs_vf_x(j, k, l, i) = v_vf(i)%sf(j, k, l)
110115
vR_rs_vf_x(j, k, l, i) = v_vf(i)%sf(j, k, l)
111116
end do
@@ -116,9 +121,9 @@ contains
116121
else if (muscl_dir == 2) then
117122
$:GPU_PARALLEL_LOOP(collapse=4)
118123
do i = 1, ubound(v_vf, 1)
119-
do l = is3_muscl%beg, is3_muscl%end
120-
do j = is1_muscl%beg, is1_muscl%end
121-
do k = is2_muscl%beg, is2_muscl%end
124+
do l = is_muscl%z%beg, is_muscl%z%end
125+
do j = is_muscl%x%beg, is_muscl%x%end
126+
do k = is_muscl%y%beg, is_muscl%y%end
122127
vL_rs_vf_x(k, j, l, i) = v_vf(i)%sf(k, j, l)
123128
vR_rs_vf_x(k, j, l, i) = v_vf(i)%sf(k, j, l)
124129
end do
@@ -129,9 +134,9 @@ contains
129134
else if (muscl_dir == 3) then
130135
$:GPU_PARALLEL_LOOP(collapse=4)
131136
do i = 1, ubound(v_vf, 1)
132-
do j = is1_muscl%beg, is1_muscl%end
133-
do k = is2_muscl%beg, is2_muscl%end
134-
do l = is3_muscl%beg, is3_muscl%end
137+
do j = is_muscl%x%beg, is_muscl%x%end
138+
do k = is_muscl%y%beg, is_muscl%y%end
139+
do l = is_muscl%z%beg, is_muscl%z%end
135140
vL_rs_vf_x(l, k, j, i) = v_vf(i)%sf(l, k, j)
136141
vR_rs_vf_x(l, k, j, i) = v_vf(i)%sf(l, k, j)
137142
end do
@@ -162,9 +167,9 @@ contains
162167
if (muscl_order == 2) then
163168
! MUSCL Reconstruction
164169
#:for MUSCL_DIR, XYZ, STENCIL_VAR, COORDS, X_BND, Y_BND, Z_BND in &
165-
[(1, 'x', 'j', '{STENCIL_IDX}, k, l', 'is1_muscl', 'is2_muscl', 'is3_muscl'), &
166-
(2, 'y', 'k', 'j, {STENCIL_IDX}, l', 'is2_muscl', 'is1_muscl', 'is3_muscl'), &
167-
(3, 'z', 'l', 'j, k, {STENCIL_IDX}', 'is3_muscl', 'is2_muscl', 'is1_muscl')]
170+
[(1, 'x', 'j', '{STENCIL_IDX}, k, l', 'is_muscl%x', 'is_muscl%y', 'is_muscl%z'), &
171+
(2, 'y', 'k', 'j, {STENCIL_IDX}, l', 'is_muscl%y', 'is_muscl%x', 'is_muscl%z'), &
172+
(3, 'z', 'l', 'j, k, {STENCIL_IDX}', 'is_muscl%z', 'is_muscl%y', 'is_muscl%x')]
168173
#:set SV = STENCIL_VAR
169174
#:set SF = lambda offs: COORDS.format(STENCIL_IDX = SV + offs)
170175
if (muscl_dir == ${MUSCL_DIR}$) then
@@ -223,7 +228,7 @@ contains
223228
call nvtxStartRange("WENO-INTCOMP")
224229
#:for MUSCL_DIR, XYZ in [(1, 'x'), (2, 'y'), (3, 'z')]
225230
if (muscl_dir == ${MUSCL_DIR}$) then
226-
call s_thinc_compression(v_rs_ws_muscl, vL_rs_vf_x, vR_rs_vf_x, muscl_dir, is1_muscl, is2_muscl, is3_muscl)
231+
call s_thinc_compression(v_rs_ws_muscl, vL_rs_vf_x, vR_rs_vf_x, muscl_dir, is_muscl%x, is_muscl%y, is_muscl%z)
227232
end if
228233
#:endfor
229234
call nvtxEndRange()

0 commit comments

Comments
 (0)