diff --git a/src/MNH/sources_neg_correct.f90 b/src/MNH/sources_neg_correct.f90 index 81d49856c651454e06a50734e65120171de54dca..838b896b9570e6a0f0e1aba6b3e7d7ea81bc9ae2 100644 --- a/src/MNH/sources_neg_correct.f90 +++ b/src/MNH/sources_neg_correct.f90 @@ -32,6 +32,7 @@ use modd_param_lima, only: lcold_lima => lcold, lrain_lima => lrain, lspro_lima xctmin_lima => xctmin, xrtmin_lima => xrtmin use mode_budget, only: Budget_store_init, Budget_store_end +use mode_mppdb use mode_msg implicit none @@ -57,6 +58,20 @@ real, dimension(:, :, :), allocatable :: zt, zexn, zlv, zls, zcph, zcor if ( krr == 0 ) return +!$acc data present( ppabst, ptht, prt, prths, prrs, prsvs, prhodj ) + +if ( mppdb_initialized ) then + !Check all IN arrays + call Mppdb_check( ppabst, "Sources_neg_correct beg:ppabst") + call Mppdb_check( ptht, "Sources_neg_correct beg:ptht") + call Mppdb_check( prt, "Sources_neg_correct beg:prt") + if ( Present( prhodj ) ) call Mppdb_check( prhodj, "Sources_neg_correct beg:prhodj") + !Check all INOUT arrays + call Mppdb_check( prths, "Sources_neg_correct beg:prths") + call Mppdb_check( prrs, "Sources_neg_correct beg:prrs") + call Mppdb_check( prsvs, "Sources_neg_correct beg:prsvs") +end if + if ( hbudname /= 'NEADV' .and. hbudname /= 'NECON' .and. hbudname /= 'NEGA' .and. hbudname /= 'NETUR' ) & call Print_msg( NVERB_WARNING, 'GEN', 'Sources_neg_correct', 'budget '//hbudname//' not yet tested' ) @@ -130,15 +145,25 @@ allocate( zt ( Size( prths, 1 ), Size( prths, 2 ), Size( prths, 3 ) ) ) allocate( zexn( Size( prths, 1 ), Size( prths, 2 ), Size( prths, 3 ) ) ) allocate( zlv ( Size( prths, 1 ), Size( prths, 2 ), Size( prths, 3 ) ) ) allocate( zcph( Size( prths, 1 ), Size( prths, 2 ), Size( prths, 3 ) ) ) +if ( hcloud == 'ICE3' .or. hcloud == 'ICE4' .or. hcloud == 'LIMA' ) then + allocate( zls( Size( prths, 1 ), Size( prths, 2 ), Size( prths, 3 ) ) ) +end if +!$acc data create( zt, zexn, zlv, zcph, zls ) + +!$acc kernels zexn(:, :, :) = ( ppabst(:, :, :) / xp00 ) ** (xrd / xcpd ) zt (:, :, :) = ptht(:, :, :) * zexn(:, :, :) zlv (:, :, :) = xlvtt + ( xcpv - xcl ) * ( zt(:, :, :) - xtt ) +!$acc end kernels if ( hcloud == 'ICE3' .or. hcloud == 'ICE4' .or. hcloud == 'LIMA' ) then - allocate( zls( Size( prths, 1 ), Size( prths, 2 ), Size( prths, 3 ) ) ) +!$acc kernels zls(:, :, :) = xlstt + ( xcpv - xci ) * ( zt(:, :, :) - xtt ) +!$acc end kernels end if +!$acc kernels zcph(:, :, :) = xcpd + xcpv * prt(:, :, :, 1) +!$acc end kernels deallocate( zt ) @@ -146,20 +171,25 @@ CLOUD: select case ( hcloud ) case ( 'KESS' ) jrmax = Size( prrs, 4 ) do jr = 2, jrmax +!PW: kernels directive inside do loop on jr because compiler bug... (NVHPC 21.7) +!$acc kernels where ( prrs(:, :, :, jr) < 0. ) prrs(:, :, :, 1) = prrs(:, :, :, 1) + prrs(:, :, :, jr) prths(:, :, :) = prths(:, :, :) - prrs(:, :, :, jr) * zlv(:, :, :) / & ( zcph(:, :, :) * zexn(:, :, :) ) prrs(:, :, :, jr) = 0. end where +!$acc end kernels end do +!$acc kernels where ( prrs(:, :, :, 1) < 0. .and. prrs(:, :, :, 2) > 0. ) prrs(:, :, :, 1) = prrs(:, :, :, 1) + prrs(:, :, :, 2) prths(:, :, :) = prths(:, :, :) - prrs(:, :, :, 2) * zlv(:, :, :) / & ( zcph(:, :, :) * zexn(:, :, :) ) prrs(:, :, :, 2) = 0. end where +!$acc end kernels case( 'ICE3', 'ICE4' ) @@ -169,12 +199,15 @@ CLOUD: select case ( hcloud ) jrmax = Size( prrs, 4 ) end if do jr = 4, jrmax +!$acc kernels +!PW: kernels directive inside do loop on jr because compiler bug... (NVHPC 21.7) where ( prrs(:, :, :, jr) < 0.) prrs(:, :, :, 1) = prrs(:, :, :, 1) + prrs(:, :, :, jr) prths(:, :, :) = prths(:, :, :) - prrs(:, :, :, jr) * zls(:, :, :) / & ( zcph(:, :, :) * zexn(:, :, :) ) prrs(:, :, :, jr) = 0. end where +!$acc end kernels end do ! ! cloud @@ -184,25 +217,31 @@ CLOUD: select case ( hcloud ) jrmax = 3 end if do jr = 2, jrmax +!$acc kernels +!PW: kernels directive inside do loop on jr because compiler bug... (NVHPC 21.7) where ( prrs(:, :, :, jr) < 0.) prrs(:, :, :, 1) = prrs(:, :, :, 1) + prrs(:, :, :, jr) prths(:, :, :) = prths(:, :, :) - prrs(:, :, :, jr) * zlv(:, :, :) / & ( zcph(:, :, :) * zexn(:, :, :) ) prrs(:, :, :, jr) = 0. end where +!$acc end kernels end do ! ! if rc or ri are positive, we can correct negative rv ! cloud +!$acc kernels where ( prrs(:, :, :, 1) < 0. .and. prrs(:, :, :, 2) > 0. ) prrs(:, :, :, 1) = prrs(:, :, :, 1) + prrs(:, :, :, 2) prths(:, :, :) = prths(:, :, :) - prrs(:, :, :, 2) * zlv(:, :, :) / & ( zcph(:, :, :) * zexn(:, :, :) ) prrs(:, :, :, 2) = 0. end where +!$acc end kernels ! ice if ( krr > 3 ) then allocate( zcor( Size( prths, 1 ), Size( prths, 2 ), Size( prths, 3 ) ) ) +!$acc kernels create( zcor ) where ( prrs(:, :, :, 1) < 0. .and. prrs(:, :, :, 4) > 0. ) zcor(:, :, :) = Min( -prrs(:, :, :, 1), prrs(:, :, :, 4) ) prrs(:, :, :, 1) = prrs(:, :, :, 1) + zcor(:, :, :) @@ -210,14 +249,20 @@ CLOUD: select case ( hcloud ) ( zcph(:, :, :) * zexn(:, :, :) ) prrs(:, :, :, 4) = prrs(:, :, :, 4) - zcor(:, :, :) end where +!$acc end kernels + deallocate( zcor ) end if ! ! case( 'C2R2', 'KHKO' ) +!$acc kernels where ( prrs(:, :, :, 2) < 0. .or. prsvs(:, :, :, nsv_c2r2beg + 1) < 0. ) prsvs(:, :, :, nsv_c2r2beg) = 0. end where +!$acc end kernels do jsv = 2, 3 +!$acc kernels +!PW: kernels directive inside do loop on jr because compiler bug... (NVHPC 21.7) where ( prrs(:, :, :, jsv) < 0. .or. prsvs(:, :, :, nsv_c2r2beg - 1 + jsv) < 0. ) prrs(:, :, :, 1) = prrs(:, :, :, 1) + prrs(:, :, :, jsv) prths(:, :, :) = prths(:, :, :) - prrs(:, :, :, jsv) * zlv(:, :, :) / & @@ -225,7 +270,9 @@ CLOUD: select case ( hcloud ) prrs(:, :, :, jsv) = 0. prsvs(:, :, :, nsv_c2r2beg - 1 + jsv) = 0. end where +!$acc end kernels end do +!$acc kernels where ( prrs(:, :, :, 1) < 0. .and. prrs(:, :, :, 2) > 0. ) prrs(:, :, :, 1) = prrs(:, :, :, 1) + prrs(:, :, :, 2) prths(:, :, :) = prths(:, :, :) - prrs(:, :, :, 2) * zlv(:, :, :) / & @@ -233,9 +280,11 @@ CLOUD: select case ( hcloud ) prrs(:, :, :, 2) = 0. prsvs(:, :, :, nsv_c2r2beg + 1) = 0. end where +!$acc end kernels ! ! case( 'LIMA' ) +!$acc kernels ! Correction where rc<0 or Nc<0 if ( lwarm_lima ) then where ( prrs(:, :, :, 2) < xrtmin_lima(2) / ptstep .or. prsvs(:, :, :, nsv_lima_nc) < xctmin_lima(2) / ptstep ) @@ -263,8 +312,10 @@ CLOUD: select case ( hcloud ) prsvs(:, :, :, nsv_lima_nr) = 0. end where end if +!$acc end kernels ! Correction where ri<0 or Ni<0 if ( lcold_lima ) then +!$acc kernels where ( prrs(:, :, :, 4) < xrtmin_lima(4) / ptstep .or. prsvs(:, :, :, nsv_lima_ni) < xctmin_lima(4) / ptstep ) prrs(:, :, :, 1) = prrs(:, :, :, 1) + prrs(:, :, :, 4) prths(:, :, :) = prths(:, :, :) - prrs(:, :, :, 4) * zls(:, :, :) / & @@ -272,18 +323,23 @@ CLOUD: select case ( hcloud ) prrs(:, :, :, 4) = 0. prsvs(:, :, :, nsv_lima_ni) = 0. end where +!$acc end kernels if ( hbudname /= 'NETUR' ) then do jr = 5, Size( prrs, 4 ) +!PW: kernels directive inside do loop on jr because compiler bug... (NVHPC 21.7) +!$acc kernels where ( prrs(:, :, :, jr) < 0. ) prrs(:, :, :, 1) = prrs(:, :, :, 1) + prrs(:, :, :, jr) prths(:, :, :) = prths(:, :, :) - prrs(:, :, :, jr) * zls(:, :, :) / & ( zcph(:, :, :) * zexn(:, :, :) ) prrs(:, :, :, jr) = 0. end where +!$acc end kernels end do end if if(krr > 3) then allocate( zcor( Size( prths, 1 ), Size( prths, 2 ), Size( prths, 3 ) ) ) +!$acc kernels create( zcor ) where ( prrs(:, :, :, 1) < 0. .and. prrs(:, :, :, 4) > 0. ) zcor(:, :, :) = Min( -prrs(:, :, :, 1), prrs(:, :, :, 4) ) prrs(:, :, :, 1) = prrs(:, :, :, 1) + zcor(:, :, :) @@ -291,14 +347,18 @@ CLOUD: select case ( hcloud ) ( zcph(:, :, :) * zexn(:, :, :) ) prrs(:, :, :, 4) = prrs(:, :, :, 4) - zcor(:, :, :) end where +!$acc end kernels deallocate( zcor ) end if end if +!$acc kernels prsvs(:, :, :, nsv_lima_beg : isv_lima_end) = Max( 0.0, prsvs(:, :, :, nsv_lima_beg : isv_lima_end) ) +!$acc end kernels end select CLOUD +!$acc end data if ( hbudname /= 'NECON' .and. hbudname /= 'NEGA' ) then if ( hcloud == 'KESS' .or. hcloud == 'ICE3' .or. hcloud == 'ICE4' .or. & @@ -354,6 +414,15 @@ else !NECON + NEGA end if end if +if ( mppdb_initialized ) then + !Check all INOUT arrays + call Mppdb_check( prths, "Sources_neg_correct end:prths") + call Mppdb_check( prrs, "Sources_neg_correct end:prrs") + call Mppdb_check( prsvs, "Sources_neg_correct end:prsvs") +end if + +!$acc end data + end subroutine Sources_neg_correct end module mode_sources_neg_correct