Skip to content

Commit 607d858

Browse files
author
G. Dylan Dickerson
committed
Initial OpenACC port of mpas_vector_reconstruct_2d
Ensures the data needed for the mpas_reconstruct_2d routine has been fetched onto the device (GPU) at the beginning and end of the routine. The time for these transfers are captured in a new timer 'mpas_reconstruct_2d [ACC_data_xfer]'. This is enforced by the default(present) clauses. NOTE: coeffs_reconstruct, nEdgesOnCell, edgesOnCell, latCell, and lonCell are also fetched in mpas_reconstruct_2d. This is because this routine is called before these variables would be uploaded to the device during mpas_atm_dynamics_init as part of atmosphere_core initialization. The copyins will not execute anymore once the model starts timestepping and the OpenACC runtime sees the variables are present on the device.
1 parent 88272f0 commit 607d858

File tree

1 file changed

+32
-0
lines changed

1 file changed

+32
-0
lines changed

src/operators/mpas_vector_reconstruction.F

Lines changed: 32 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -258,10 +258,23 @@ subroutine mpas_reconstruct_2d(meshPool, u, uReconstructX, uReconstructY, uRecon
258258

259259
call mpas_pool_get_config(meshPool, 'on_a_sphere', on_a_sphere)
260260

261+
MPAS_ACC_TIMER_START('mpas_reconstruct_2d [ACC_data_xfer]')
262+
! Only use sections needed, nCells may be all cells or only non-halo cells
263+
!$acc enter data copyin(coeffs_reconstruct(:,:,1:nCells),nEdgesOnCell(1:nCells), &
264+
!$acc edgesOnCell(:,1:nCells),latCell(1:nCells),lonCell(1:nCells))
265+
!$acc enter data copyin(u(:,:))
266+
!$acc enter data create(uReconstructX(:,1:nCells),uReconstructY(:,1:nCells), &
267+
!$acc uReconstructZ(:,1:nCells),uReconstructZonal(:,1:nCells), &
268+
!$acc uReconstructMeridional(:,1:nCells))
269+
MPAS_ACC_TIMER_STOP('mpas_reconstruct_2d [ACC_data_xfer]')
270+
261271
! loop over cell centers
262272
!$omp do schedule(runtime)
273+
!$acc parallel default(present)
274+
!$acc loop gang
263275
do iCell = 1, nCells
264276
! initialize the reconstructed vectors
277+
!$acc loop vector
265278
do k = 1, nVertLevels
266279
uReconstructX(k,iCell) = 0.0
267280
uReconstructY(k,iCell) = 0.0
@@ -270,8 +283,10 @@ subroutine mpas_reconstruct_2d(meshPool, u, uReconstructX, uReconstructY, uRecon
270283

271284
! a more efficient reconstruction where rbf_values*matrix_reconstruct
272285
! has been precomputed in coeffs_reconstruct
286+
!$acc loop seq
273287
do i = 1, nEdgesOnCell(iCell)
274288
iEdge = edgesOnCell(i,iCell)
289+
!$acc loop vector
275290
do k = 1, nVertLevels
276291
uReconstructX(k,iCell) = uReconstructX(k,iCell) &
277292
+ coeffs_reconstruct(1,i,iCell) * u(k,iEdge)
@@ -283,17 +298,21 @@ subroutine mpas_reconstruct_2d(meshPool, u, uReconstructX, uReconstructY, uRecon
283298

284299
enddo
285300
enddo ! iCell
301+
!$acc end parallel
286302
!$omp end do
287303

288304
call mpas_threading_barrier()
289305

290306
if (on_a_sphere) then
291307
!$omp do schedule(runtime)
308+
!$acc parallel default(present)
309+
!$acc loop gang
292310
do iCell = 1, nCells
293311
clat = cos(latCell(iCell))
294312
slat = sin(latCell(iCell))
295313
clon = cos(lonCell(iCell))
296314
slon = sin(lonCell(iCell))
315+
!$acc loop vector
297316
do k = 1, nVertLevels
298317
uReconstructZonal(k,iCell) = -uReconstructX(k,iCell)*slon + &
299318
uReconstructY(k,iCell)*clon
@@ -302,18 +321,31 @@ subroutine mpas_reconstruct_2d(meshPool, u, uReconstructX, uReconstructY, uRecon
302321
+ uReconstructZ(k,iCell)*clat
303322
end do
304323
end do
324+
!$acc end parallel
305325
!$omp end do
306326
else
307327
!$omp do schedule(runtime)
328+
!$acc parallel default(present)
329+
!$acc loop gang vector collapse(2)
308330
do iCell = 1, nCells
309331
do k = 1, nVertLevels
310332
uReconstructZonal (k,iCell) = uReconstructX(k,iCell)
311333
uReconstructMeridional(k,iCell) = uReconstructY(k,iCell)
312334
end do
313335
end do
336+
!$acc end parallel
314337
!$omp end do
315338
end if
316339

340+
MPAS_ACC_TIMER_START('mpas_reconstruct_2d [ACC_data_xfer]')
341+
!$acc exit data delete(coeffs_reconstruct(:,:,1:nCells),nEdgesOnCell(1:nCells), &
342+
!$acc edgesOnCell(:,1:nCells),latCell(1:nCells),lonCell(1:nCells))
343+
!$acc exit data delete(u(:,:))
344+
!$acc exit data copyout(uReconstructX(:,1:nCells),uReconstructY(:,1:nCells), &
345+
!$acc uReconstructZ(:,1:nCells), uReconstructZonal(:,1:nCells), &
346+
!$acc uReconstructMeridional(:,1:nCells))
347+
MPAS_ACC_TIMER_STOP('mpas_reconstruct_2d [ACC_data_xfer]')
348+
317349
end subroutine mpas_reconstruct_2d!}}}
318350

319351

0 commit comments

Comments
 (0)