diff --git a/src/core_atmosphere/dynamics/mpas_atm_time_integration.F b/src/core_atmosphere/dynamics/mpas_atm_time_integration.F index 319ca96244..657cc98ff4 100644 --- a/src/core_atmosphere/dynamics/mpas_atm_time_integration.F +++ b/src/core_atmosphere/dynamics/mpas_atm_time_integration.F @@ -233,6 +233,7 @@ subroutine mpas_atm_dynamics_init(domain) real (kind=RKIND), dimension(:,:,:), pointer :: zb3_cell real (kind=RKIND), dimension(:), pointer :: fzm real (kind=RKIND), dimension(:), pointer :: fzp + real (kind=RKIND), dimension(:), pointer :: meshScalingRegionalCell #endif @@ -356,6 +357,9 @@ subroutine mpas_atm_dynamics_init(domain) call mpas_pool_get_array(mesh, 'fzp', fzp) !$acc enter data copyin(fzp) + call mpas_pool_get_array(mesh, 'meshScalingRegionalCell', meshScalingRegionalCell) + !$acc enter data copyin(meshScalingRegionalCell) + #endif end subroutine mpas_atm_dynamics_init @@ -425,6 +429,7 @@ subroutine mpas_atm_dynamics_finalize(domain) real (kind=RKIND), dimension(:,:,:), pointer :: zb3_cell real (kind=RKIND), dimension(:), pointer :: fzm real (kind=RKIND), dimension(:), pointer :: fzp + real (kind=RKIND), dimension(:), pointer :: meshScalingRegionalCell #endif @@ -547,6 +552,10 @@ subroutine mpas_atm_dynamics_finalize(domain) call mpas_pool_get_array(mesh, 'fzp', fzp) !$acc exit data delete(fzp) + + call mpas_pool_get_array(mesh, 'meshScalingRegionalCell', meshScalingRegionalCell) + !$acc exit data delete(meshScalingRegionalCell) + #endif end subroutine mpas_atm_dynamics_finalize @@ -6780,17 +6789,29 @@ subroutine atm_bdy_adjust_scalars_work( scalars_new, scalars_driving, dt, dt_rk, integer :: iCell, iEdge, iScalar, i, k, cell1, cell2 !--- + MPAS_ACC_TIMER_START('atm_bdy_adjust_scalars [ACC_data_xfer]') + !$acc enter data create(scalars_tmp) & + !$acc copyin(scalars_driving, scalars_new) + MPAS_ACC_TIMER_STOP('atm_bdy_adjust_scalars [ACC_data_xfer]') + !$acc parallel default(present) + !$acc loop gang worker do iCell = cellSolveStart, cellSolveEnd ! threaded over cells if ( (bdyMaskCell(iCell) > 1) .and. (bdyMaskCell(iCell) <= nRelaxZone) ) then ! relaxation zone laplacian_filter_coef = dt_rk*(real(bdyMaskCell(iCell)) - 1.)/real(nRelaxZone)/(10.*dt*meshScalingRegionalCell(iCell)) rayleigh_damping_coef = laplacian_filter_coef/5.0 - scalars_tmp(1:num_scalars,1:nVertLevels,iCell) = scalars_new(1:num_scalars,1:nVertLevels,iCell) + !$acc loop vector collapse(2) + do k=1,nVertLevels + do iScalar=1,num_scalars + scalars_tmp(iScalar,k,iCell) = scalars_new(iScalar,k,iCell) + end do + end do ! first, we compute the 2nd-order laplacian filter ! + !$acc loop seq do i=1,nEdgesOnCell(iCell) iEdge = edgesOnCell(i,iCell) ! edge_sign = r_areaCell*edgesOnCell_sign(i,iCell) * dvEdge(iEdge) * invDcEdge(iEdge) * laplacian_filter_coef @@ -6799,6 +6820,7 @@ subroutine atm_bdy_adjust_scalars_work( scalars_new, scalars_driving, dt, dt_rk, cell1 = cellsOnEdge(1,iEdge) cell2 = cellsOnEdge(2,iEdge) !DIR$ IVDEP + !$acc loop vector collapse(2) do k=1,nVertLevels do iScalar = 1, num_scalars filter_flux = edge_sign*( (scalars_new(iScalar,k,cell2)-scalars_driving(iScalar,k,cell2)) & @@ -6811,6 +6833,7 @@ subroutine atm_bdy_adjust_scalars_work( scalars_new, scalars_driving, dt, dt_rk, ! second, we compute the Rayleigh damping component ! !DIR$ IVDEP + !$acc loop vector collapse(2) do k=1,nVertLevels do iScalar = 1, num_scalars scalars_tmp(iScalar,k,iCell) =scalars_tmp(iScalar,k,iCell) - rayleigh_damping_coef * (scalars_new(iScalar,k,iCell)-scalars_driving(iScalar,k,iCell)) @@ -6818,10 +6841,11 @@ subroutine atm_bdy_adjust_scalars_work( scalars_new, scalars_driving, dt, dt_rk, end do else if ( bdyMaskCell(iCell) > nRelaxZone) then ! specified zone - + ! update the specified-zone values ! !DIR$ IVDEP + !$acc loop vector collapse(2) do k=1,nVertLevels do iScalar = 1, num_scalars scalars_tmp(iScalar,k,iCell) = scalars_driving(iScalar,k,iCell) @@ -6831,12 +6855,16 @@ subroutine atm_bdy_adjust_scalars_work( scalars_new, scalars_driving, dt, dt_rk, end if end do ! updates now in temp storage - + !$acc end parallel + !$OMP BARRIER + !$acc parallel default(present) + !$acc loop gang worker do iCell = cellSolveStart, cellSolveEnd ! threaded over cells if (bdyMaskCell(iCell) > 1) then ! update values !DIR$ IVDEP + !$acc loop vector collapse(2) do k=1,nVertLevels do iScalar = 1, num_scalars scalars_new(iScalar,k,iCell) = scalars_tmp(iScalar,k,iCell) @@ -6844,6 +6872,12 @@ subroutine atm_bdy_adjust_scalars_work( scalars_new, scalars_driving, dt, dt_rk, end do end if end do + !$acc end parallel + + MPAS_ACC_TIMER_START('atm_bdy_adjust_scalars [ACC_data_xfer]') + !$acc exit data delete(scalars_tmp, scalars_driving) & + !$acc copyout(scalars_new) + MPAS_ACC_TIMER_STOP('atm_bdy_adjust_scalars [ACC_data_xfer]') end subroutine atm_bdy_adjust_scalars_work