New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 13882 for NEMO/branches/2020/dev_r13508_HPC-09_loop_fusion/src/OCE/TRA/traadv_fct_lf.F90 – NEMO

Ignore:
Timestamp:
2020-11-26T10:52:00+01:00 (3 years ago)
Author:
francesca
Message:

loop fusion v2 - mus and fct advection schemes - ticket #2367

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/dev_r13508_HPC-09_loop_fusion/src/OCE/TRA/traadv_fct_lf.F90

    r13881 r13882  
    4747#  include "domzgr_substitute.h90" 
    4848 
     49#define tracer_flux_i(out,zfp,zfm,ji,jj,jk) \ 
     50        zfp = pU(ji,jj,jk) + ABS( pU(ji,jj,jk) ) ; \ 
     51        zfm = pU(ji,jj,jk) - ABS( pU(ji,jj,jk) ) ; \ 
     52        out = 0.5 * ( zfp * pt(ji,jj,jk,jn,Kbb) + zfm * pt(ji+1,jj,jk,jn,Kbb) ) 
     53 
     54#define tracer_flux_j(out,zfp,zfm,ji,jj,jk) \ 
     55        zfp = pV(ji,jj,jk) + ABS( pV(ji,jj,jk) ) ; \ 
     56        zfm = pV(ji,jj,jk) - ABS( pV(ji,jj,jk) ) ; \ 
     57        out = 0.5 * ( zfp * pt(ji,jj,jk,jn,Kbb) + zfm * pt(ji,jj+1,jk,jn,Kbb) ) 
     58 
    4959#define search_in_neighbour(out,OP,vec,ji,jj,jk) \ 
    50 out = OP(vec(ji,jj,jk),vec(ji-1,jj,jk),vec(ji+1,jj,jk),vec(ji,jj-1,jk),vec(ji,jj+1,jk),vec(ji,jj,MAX(jk-1,1)),vec(ji,jj,jk+1)) 
     60        out = OP(vec(ji,jj,jk),vec(ji-1,jj,jk),vec(ji+1,jj,jk), \ 
     61        vec(ji,jj-1,jk),vec(ji,jj+1,jk),vec(ji,jj,MAX(jk-1,1)),vec(ji,jj,jk+1)) 
    5162 
    5263#define pos_part_of_flux(ji,jj,jk,out) \ 
    53 out = MAX(0.,paa_in(ji-1,jj,jk)) - MIN(0.,paa_in(ji,jj,jk)) \ 
    54 + MAX(0.,pbb_in(ji,jj-1,jk)) - MIN(0.,pbb_in(ji,jj,jk)) \ 
    55 + MAX(0.,pcc_in(ji,jj,jk+1)) - MIN(0.,pcc_in(ji,jj,jk))  
     64        out = MAX(0.,paa_in(ji-1,jj,jk)) - MIN(0.,paa_in(ji,jj,jk)) \ 
     65        + MAX(0.,pbb_in(ji,jj-1,jk)) - MIN(0.,pbb_in(ji,jj,jk)) \ 
     66        + MAX(0.,pcc_in(ji,jj,jk+1)) - MIN(0.,pcc_in(ji,jj,jk))  
    5667 
    5768#define neg_part_of_flux(ji,jj,jk,out) \ 
    58 out = MAX( 0.,paa_in(ji,jj,jk) ) - MIN( 0., paa_in(ji-1,jj,jk)) \  
    59 + MAX( 0.,pbb_in(ji,jj,jk) ) - MIN( 0., pbb_in(ji,jj-1,jk)) \ 
    60 + MAX( 0.,pcc_in(ji,jj,jk) ) - MIN( 0., pcc_in(ji,jj,jk+1)) 
     69        out = MAX( 0.,paa_in(ji,jj,jk) ) - MIN( 0., paa_in(ji-1,jj,jk)) \  
     70        + MAX( 0.,pbb_in(ji,jj,jk) ) - MIN( 0., pbb_in(ji,jj-1,jk)) \ 
     71        + MAX( 0.,pcc_in(ji,jj,jk) ) - MIN( 0., pcc_in(ji,jj,jk+1)) 
    6172 
    6273#define beta_terms(bt,betup,betdo,up,pos,do,neg,ji,jj,jk) \ 
    63 bt = e1e2t(ji,jj) * e3t(ji,jj,jk,Kmm) / p2dt ; \ 
    64 betup = ( up            - paft(ji,jj,jk) ) / ( pos + zrtrn ) * bt ; \ 
    65 betdo = ( paft(ji,jj,jk) - do            ) / ( neg + zrtrn ) * bt 
     74        bt = e1e2t(ji,jj) * e3t(ji,jj,jk,Kmm) / p2dt ; \ 
     75        betup = ( up            - paft(ji,jj,jk) ) / ( pos + zrtrn ) * bt ; \ 
     76        betdo = ( paft(ji,jj,jk) - do            ) / ( neg + zrtrn ) * bt 
    6677 
    6778#define monotonic_flux(a,b,c,betup_p1,betdo_p1,vec,vec_in,jk) \ 
    68 a = MIN( 1._wp, zbetdo, betup_p1 ) ; \ 
    69 b = MIN( 1._wp, zbetup, betdo_p1 ) ; \ 
    70 c = ( 0.5_wp  + SIGN( 0.5_wp , vec_in(ji,jj,jk) ) ) ; \ 
    71 vec(ji,jj,jk) = vec_in(ji,jj,jk) * ( c * a + ( 1._wp - c) * b ) 
    72  
    73 #define monotonic_flux_k(a,b,c,betup_p1,betdo_p1,vec,vec_in,jk) \ 
    74 a = MIN( 1._wp, betdo_p1, zbetup ) ; \ 
    75 b = MIN( 1._wp, betup_p1, zbetdo ) ; \ 
    76 c = ( 0.5 + SIGN( 0.5_wp , vec_in(ji,jj,jk+1) ) ) ; \ 
    77 vec(ji,jj,jk+1) = vec_in(ji,jj,jk+1) * ( c * a + ( 1._wp - c) * b ) 
     79        a = MIN( 1._wp, zbetdo(ji,jj), betup_p1 ) ; \ 
     80        b = MIN( 1._wp, zbetup(ji,jj), betdo_p1 ) ; \ 
     81        c = ( 0.5_wp  + SIGN( 0.5_wp , vec_in(ji,jj,jk) ) ) ; \ 
     82        vec(ji,jj,jk) = vec_in(ji,jj,jk) * ( c * a + ( 1._wp - c) * b ) 
     83 
     84#define monotonic_flux_k(a,b,c,betup,betdo,vec,vec_in,jk) \ 
     85        a = MIN( 1._wp, betdo, zbetup_ptr(ji,jj) ) ; \ 
     86        b = MIN( 1._wp, betup, zbetdo_ptr(ji,jj) ) ; \ 
     87        c = ( 0.5 + SIGN( 0.5_wp , vec_in(ji,jj,jk) ) ) ; \ 
     88        vec(ji,jj,jk) = vec_in(ji,jj,jk) * ( c * a + ( 1._wp - c) * b ) 
    7889 
    7990   !!---------------------------------------------------------------------- 
     
    114125      INTEGER  ::   ji, jj, jk, jn                           ! dummy loop indices   
    115126      REAL(wp) ::   ztra                                     ! local scalar 
    116       REAL(wp) ::   zwx, zwx_im1, zfp_ui, zfp_ui_m1, zfp_vj, zfp_vj_m1, zfp_wk, zC2t_u, zC4t_u   !   -      - 
    117       REAL(wp) ::   zwy, zwy_jm1, zfm_ui, zfm_ui_m1, zfm_vj, zfm_vj_m1, zfm_wk, zC2t_v, zC4t_v   !   -      - 
    118       REAL(wp) ::   ztu_im1, ztu_ip1, ztv_jm1, ztv_jp1, zltu_ip1, zltv_jp1, zltu, zltv 
     127      REAL(wp) ::   zwx_im1, zfp_ui, zfp_ui_m1, zfp_vj, zfp_vj_m1, zfp_wk, zC2t_u, zC4t_u   !   -      - 
     128      REAL(wp) ::   zwy_jm1, zfm_ui, zfm_ui_m1, zfm_vj, zfm_vj_m1, zfm_wk, zC2t_v, zC4t_v   !   -      - 
     129      REAL(wp) ::   ztu_im1, ztu_ip1, ztv_jm1, ztv_jp1   
    119130      REAL(wp), DIMENSION(jpi,jpj,jpk)        ::   zwi, zwx_3d, zwy_3d, zwz, ztw, zltu_3d, zltv_3d, ztu, ztv 
    120131      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ztrdx, ztrdy, ztrdz, zptry 
     
    129140      ENDIF 
    130141      !! -- init to 0 
    131       !! -- init to 0 
    132       zwi(:,:,:) = 0._wp 
    133142      zwx_3d(:,:,:) = 0._wp 
    134143      zwy_3d(:,:,:) = 0._wp 
     
    173182         ! 
    174183         !        !==  upstream advection with initial mass fluxes & intermediate update  ==! 
    175          !                    !* upstream tracer flux in the i and j direction  
    176184         !                               !* upstream tracer flux in the k direction *! 
    177185         DO_3D( 1, 1, 1, 1, 2, jpkm1 )      ! Interior value ( multiplied by wmask) 
     
    192200         ENDIF 
    193201         !                
     202         !                    !* upstream tracer flux in the i and j direction  
    194203         DO jk = 1, jpkm1 
    195204            DO jj = 1, jpj-1 
    196                zfp_ui = pU(1,jj,jk) + ABS( pU(1,jj,jk) ) 
    197                zfm_ui = pU(1,jj,jk) - ABS( pU(1,jj,jk) ) 
    198                zwx = 0.5 * ( zfp_ui * pt(1,jj,jk,jn,Kbb) + zfm_ui * pt(2,jj  ,jk,jn,Kbb) ) 
    199                zwx_3d(1,jj,jk) = zwx 
    200  
    201                zfp_vj = pV(1,jj,jk) + ABS( pV(1,jj,jk) ) 
    202                zfm_vj = pV(1,jj,jk) - ABS( pV(1,jj,jk) ) 
    203                zwy = 0.5 * ( zfp_vj * pt(1,jj,jk,jn,Kbb) + zfm_vj * pt(1  ,jj+1,jk,jn,Kbb) ) 
    204                zwy_3d(1,jj,jk) = zwy 
     205               tracer_flux_i(zwx_3d(1,jj,jk),zfp_ui,zfm_ui,1,jj,jk) 
     206               tracer_flux_j(zwy_3d(1,jj,jk),zfp_vj,zfm_vj,1,jj,jk) 
    205207            END DO 
    206208            DO ji = 1, jpi-1 
    207                zfp_ui = pU(ji,1,jk) + ABS( pU(ji,1,jk) ) 
    208                zfm_ui = pU(ji,1,jk) - ABS( pU(ji,1,jk) ) 
    209                zwx = 0.5 * ( zfp_ui * pt(ji,1,jk,jn,Kbb) + zfm_ui * pt(ji+1,1  ,jk,jn,Kbb) ) 
    210                zwx_3d(ji,1,jk) = zwx 
    211  
    212                zfp_vj = pV(ji,1,jk) + ABS( pV(ji,1,jk) ) 
    213                zfm_vj = pV(ji,1,jk) - ABS( pV(ji,1,jk) ) 
    214                zwy = 0.5 * ( zfp_vj * pt(ji,1,jk,jn,Kbb) + zfm_vj * pt(ji  ,2,jk,jn,Kbb) ) 
    215                zwy_3d(ji,1,jk) = zwy 
     209               tracer_flux_i(zwx_3d(ji,1,jk),zfp_ui,zfm_ui,ji,1,jk) 
     210               tracer_flux_j(zwy_3d(ji,1,jk),zfp_vj,zfm_vj,ji,1,jk) 
    216211            END DO 
    217212            DO_2D( 1, 1, 1, 1 ) 
    218                zfp_ui = pU(ji,jj,jk) + ABS( pU(ji,jj,jk) ) 
    219                zfm_ui = pU(ji,jj,jk) - ABS( pU(ji,jj,jk) ) 
    220                zwx = 0.5 * ( zfp_ui * pt(ji,jj,jk,jn,Kbb) + zfm_ui * pt(ji+1,jj  ,jk,jn,Kbb) ) 
    221                zwx_3d(ji,jj,jk) = zwx 
    222  
    223                zfp_ui_m1 = pU(ji-1,jj,jk) + ABS( pU(ji-1,jj,jk) ) 
    224                zfm_ui_m1 = pU(ji-1,jj,jk) - ABS( pU(ji-1,jj,jk) ) 
    225                zwx_im1 = 0.5 * ( zfp_ui_m1 * pt(ji-1,jj,jk,jn,Kbb) + zfm_ui_m1 * pt(ji,jj  ,jk,jn,Kbb) ) 
    226  
    227                zfp_vj = pV(ji,jj,jk) + ABS( pV(ji,jj,jk) ) 
    228                zfm_vj = pV(ji,jj,jk) - ABS( pV(ji,jj,jk) ) 
    229                zwy = 0.5 * ( zfp_vj * pt(ji,jj,jk,jn,Kbb) + zfm_vj * pt(ji  ,jj+1,jk,jn,Kbb) ) 
    230                zwy_3d(ji,jj,jk) = zwy 
    231  
    232                zfp_vj_m1 = pV(ji,jj-1,jk) + ABS( pV(ji,jj-1,jk) ) 
    233                zfm_vj_m1 = pV(ji,jj-1,jk) - ABS( pV(ji,jj-1,jk) ) 
    234                zwy_jm1 = 0.5 * ( zfp_vj_m1 * pt(ji,jj-1,jk,jn,Kbb) + zfm_vj_m1 * pt(ji,jj,jk,jn,Kbb) ) 
    235                !                               ! total intermediate advective trends 
    236                ztra = - (  zwx - zwx_im1 + zwy - zwy_jm1 + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj) 
     213               tracer_flux_i(zwx_3d(ji,jj,jk),zfp_ui,zfm_ui,ji,jj,jk) 
     214               tracer_flux_i(zwx_im1,zfp_ui_m1,zfm_ui_m1,ji-1,jj,jk) 
     215               tracer_flux_j(zwy_3d(ji,jj,jk),zfp_vj,zfm_vj,ji,jj,jk) 
     216               tracer_flux_j(zwy_jm1,zfp_vj_m1,zfm_vj_m1,ji,jj-1,jk) 
     217               ztra = - ( zwx_3d(ji,jj,jk) - zwx_im1 + zwy_3d(ji,jj,jk) - zwy_jm1 + zwz(ji,jj,jk) - zwz(ji,jj,jk+1) ) * r1_e1e2t(ji,jj) 
    237218               !                               ! update and guess with monotonic sheme 
    238219               pt(ji,jj,jk,jn,Krhs) =                   pt(ji,jj,jk,jn,Krhs) +       ztra   & 
     
    449430      ! 
    450431      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
    451       INTEGER  ::   ikm1         ! local integer 
    452432      REAL(dp) ::   zpos, zneg, zbt, za, zb, zc, zbig, zrtrn    ! local scalars 
    453433      REAL(dp) ::   zau, zbu, zcu, zav, zbv, zcv, zup, zdo            !   -      - 
    454       REAL(dp) ::   zbetup, zbetdo 
    455434      REAL(dp) ::   zbt_ip1, zpos_ip1, zneg_ip1, zup_ip1, zdo_ip1, zbetup_ip1, zbetdo_ip1  
    456435      REAL(dp) ::   zbt_jp1, zpos_jp1, zneg_jp1, zup_jp1, zdo_jp1, zbetup_jp1, zbetdo_jp1 
    457       REAL(dp) ::   zbt_kp1, zpos_kp1, zneg_kp1, zup_kp1, zdo_kp1, zbetup_kp1, zbetdo_kp1 
     436      REAL(dp), TARGET, DIMENSION(jpi,jpj) :: zbetup_buf, zbetdo_buf, zbetup_ptr_buf, zbetdo_ptr_buf 
     437      REAL(dp), POINTER, DIMENSION(:,:) :: tmp, zbetup, zbetdo, zbetup_ptr, zbetdo_ptr 
    458438      !!---------------------------------------------------------------------- 
    459439      ! 
     
    472452         &        paft * tmask + zbig * ( 1._wp - tmask )  ) 
    473453 
    474       DO_3D( 1, 0, 1, 0, 1, jpk-2 ) 
     454      zbetup => zbetup_buf 
     455      zbetdo => zbetdo_buf 
     456      zbetup_ptr => zbetup_ptr_buf 
     457      zbetdo_ptr => zbetdo_ptr_buf 
     458 
     459      DO_2D( 1, 0, 1, 0 ) 
    475460 
    476461        ! search maximum in neighbourhood 
    477         search_in_neighbour(zup,MAX,zbup,ji,jj,jk) 
    478         search_in_neighbour(zup_ip1,MAX,zbup,ji+1,jj,jk) 
    479         search_in_neighbour(zup_jp1,MAX,zbup,ji,jj+1,jk) 
    480         search_in_neighbour(zup_kp1,MAX,zbup,ji,jj,jk+1) 
     462        search_in_neighbour(zup,MAX,zbup,ji,jj,1) 
     463        search_in_neighbour(zup_ip1,MAX,zbup,ji+1,jj,1) 
     464        search_in_neighbour(zup_jp1,MAX,zbup,ji,jj+1,1) 
    481465 
    482466        ! search minimum in neighbourhood 
    483         search_in_neighbour(zdo,MIN,zbdo,ji,jj,jk) 
    484         search_in_neighbour(zdo_ip1,MIN,zbdo,ji+1,jj,jk) 
    485         search_in_neighbour(zdo_jp1,MIN,zbdo,ji,jj+1,jk) 
    486         search_in_neighbour(zdo_kp1,MIN,zbdo,ji,jj,jk+1) 
     467        search_in_neighbour(zdo,MIN,zbdo,ji,jj,1) 
     468        search_in_neighbour(zdo_ip1,MIN,zbdo,ji+1,jj,1) 
     469        search_in_neighbour(zdo_jp1,MIN,zbdo,ji,jj+1,1) 
    487470 
    488471        ! positive part of the flux 
    489         pos_part_of_flux(ji,jj,jk,zpos) 
    490         pos_part_of_flux(ji+1,jj,jk,zpos_ip1) 
    491         pos_part_of_flux(ji,jj+1,jk,zpos_jp1) 
    492         pos_part_of_flux(ji,jj,jk+1,zpos_kp1) 
     472        pos_part_of_flux(ji,jj,1,zpos) 
     473        pos_part_of_flux(ji+1,jj,1,zpos_ip1) 
     474        pos_part_of_flux(ji,jj+1,1,zpos_jp1) 
    493475 
    494476        ! negative part of the flux 
    495         neg_part_of_flux(ji,jj,jk,zneg) 
    496         neg_part_of_flux(ji+1,jj,jk,zneg_ip1) 
    497         neg_part_of_flux(ji,jj+1,jk,zneg_jp1) 
    498         neg_part_of_flux(ji,jj,jk+1,zneg_kp1) 
     477        neg_part_of_flux(ji,jj,1,zneg) 
     478        neg_part_of_flux(ji+1,jj,1,zneg_ip1) 
     479        neg_part_of_flux(ji,jj+1,1,zneg_jp1) 
    499480 
    500481        ! up & down beta terms 
    501         beta_terms(zbt,zbetup,zbetdo,zup,zpos,zdo,zneg,ji,jj,jk) 
    502         beta_terms(zbt_ip1,zbetup_ip1,zbetdo_ip1,zup_ip1,zpos_ip1,zdo_ip1,zneg_ip1,ji+1,jj,jk) 
    503         beta_terms(zbt_jp1,zbetup_jp1,zbetdo_jp1,zup_jp1,zpos_jp1,zdo_jp1,zneg_jp1,ji,jj+1,jk) 
    504         beta_terms(zbt_kp1,zbetup_kp1,zbetdo_kp1,zup_kp1,zpos_kp1,zdo_kp1,zneg_kp1,ji,jj,jk+1) 
     482        beta_terms(zbt,zbetup(ji,jj),zbetdo(ji,jj),zup,zpos,zdo,zneg,ji,jj,1) 
     483        beta_terms(zbt_ip1,zbetup_ip1,zbetdo_ip1,zup_ip1,zpos_ip1,zdo_ip1,zneg_ip1,ji+1,jj,1) 
     484        beta_terms(zbt_jp1,zbetup_jp1,zbetdo_jp1,zup_jp1,zpos_jp1,zdo_jp1,zneg_jp1,ji,jj+1,1) 
    505485 
    506486        ! 3. monotonic flux in the i & j (paa & pbb) 
    507487        ! ---------------------------------------- 
    508         monotonic_flux(zau,zbu,zcu,zbetup_ip1,zbetdo_ip1,paa,paa_in,jk) 
    509         monotonic_flux(zav,zbv,zcv,zbetup_jp1,zbetdo_jp1,pbb,pbb_in,jk) 
    510  
     488        monotonic_flux(zau,zbu,zcu,zbetup_ip1,zbetdo_ip1,paa,paa_in,1) 
     489        monotonic_flux(zav,zbv,zcv,zbetup_jp1,zbetdo_jp1,pbb,pbb_in,1) 
     490 
     491      END_2D 
     492      tmp => zbetup_ptr 
     493      zbetup_ptr => zbetup 
     494      zbetup => tmp 
     495 
     496      tmp => zbetdo_ptr 
     497      zbetdo_ptr => zbetdo 
     498      zbetdo => tmp 
     499 
     500      DO jk = 2, jpk-1 
     501         DO_2D( 1, 0, 1, 0 ) 
     502 
     503            ! search maximum in neighbourhood 
     504            search_in_neighbour(zup,MAX,zbup,ji,jj,jk) 
     505            search_in_neighbour(zup_ip1,MAX,zbup,ji+1,jj,jk) 
     506            search_in_neighbour(zup_jp1,MAX,zbup,ji,jj+1,jk) 
     507 
     508            ! search minimum in neighbourhood 
     509            search_in_neighbour(zdo,MIN,zbdo,ji,jj,jk) 
     510            search_in_neighbour(zdo_ip1,MIN,zbdo,ji+1,jj,jk) 
     511            search_in_neighbour(zdo_jp1,MIN,zbdo,ji,jj+1,jk) 
     512 
     513            ! positive part of the flux 
     514            pos_part_of_flux(ji,jj,jk,zpos) 
     515            pos_part_of_flux(ji+1,jj,jk,zpos_ip1) 
     516            pos_part_of_flux(ji,jj+1,jk,zpos_jp1) 
     517 
     518            ! negative part of the flux 
     519            neg_part_of_flux(ji,jj,jk,zneg) 
     520            neg_part_of_flux(ji+1,jj,jk,zneg_ip1) 
     521            neg_part_of_flux(ji,jj+1,jk,zneg_jp1) 
     522 
     523            ! up & down beta terms 
     524            beta_terms(zbt,zbetup(ji,jj),zbetdo(ji,jj),zup,zpos,zdo,zneg,ji,jj,jk) 
     525            beta_terms(zbt_ip1,zbetup_ip1,zbetdo_ip1,zup_ip1,zpos_ip1,zdo_ip1,zneg_ip1,ji+1,jj,jk) 
     526            beta_terms(zbt_jp1,zbetup_jp1,zbetdo_jp1,zup_jp1,zpos_jp1,zdo_jp1,zneg_jp1,ji,jj+1,jk) 
     527 
     528            ! 3. monotonic flux in the i & j (paa & pbb) 
     529            ! ---------------------------------------- 
     530            monotonic_flux(zau,zbu,zcu,zbetup_ip1,zbetdo_ip1,paa,paa_in,jk) 
     531            monotonic_flux(zav,zbv,zcv,zbetup_jp1,zbetdo_jp1,pbb,pbb_in,jk) 
     532 
     533            ! monotonic flux in the k direction, i.e. pcc 
     534            ! ------------------------------------------- 
     535            monotonic_flux_k(za,zb,zc,zbetup(ji,jj),zbetdo(ji,jj),pcc,pcc_in,jk) 
     536         END_2D 
     537         tmp => zbetup_ptr 
     538         zbetup_ptr => zbetup 
     539         zbetup => tmp 
     540 
     541         tmp => zbetdo_ptr 
     542         zbetdo_ptr => zbetdo 
     543         zbetdo => tmp 
     544      END DO 
     545      ! 
     546      DO_2D( 1, 0, 1, 0 ) 
    511547        ! monotonic flux in the k direction, i.e. pcc 
    512548        ! ------------------------------------------- 
    513         monotonic_flux_k(za,zb,zc,zbetup_kp1,zbetdo_kp1,pcc,pcc_in,jk) 
    514       END_3D 
    515       ! 
    516       DO_2D( 1, 0, 1, 0 ) 
    517  
    518         ! search maximum in neighbourhood 
    519         search_in_neighbour(zup,MAX,zbup,ji,jj,jpk-1) 
    520         search_in_neighbour(zup_ip1,MAX,zbup,ji+1,jj,jpk-1) 
    521         search_in_neighbour(zup_jp1,MAX,zbup,ji,jj+1,jpk-1) 
    522  
    523         ! search minimum in neighbourhood 
    524         search_in_neighbour(zdo,MIN,zbdo,ji,jj,jk) 
    525         search_in_neighbour(zdo_ip1,MIN,zbdo,ji+1,jj,jpk-1) 
    526         search_in_neighbour(zdo_jp1,MIN,zbdo,ji,jj+1,jpk-1) 
    527  
    528         ! positive part of the flux 
    529         pos_part_of_flux(ji,jj,jpk-1,zpos) 
    530         pos_part_of_flux(ji+1,jj,jpk-1,zpos_ip1) 
    531         pos_part_of_flux(ji,jj+1,jpk-1,zpos_jp1) 
    532  
    533         ! negative part of the flux 
    534         neg_part_of_flux(ji,jj,jpk-1,zneg) 
    535         neg_part_of_flux(ji+1,jj,jpk-1,zneg_ip1) 
    536         neg_part_of_flux(ji,jj+1,jpk-1,zneg_jp1) 
    537  
    538         ! up & down beta terms 
    539         beta_terms(zbt,zbetup,zbetdo,zup,zpos,zdo,zneg,ji,jj,jpk-1) 
    540         beta_terms(zbt_ip1,zbetup_ip1,zbetdo_ip1,zup_ip1,zpos_ip1,zdo_ip1,zneg_ip1,ji+1,jj,jpk-1) 
    541         beta_terms(zbt_jp1,zbetup_jp1,zbetdo_jp1,zup_jp1,zpos_jp1,zdo_jp1,zneg_jp1,ji,jj+1,jpk-1) 
    542  
    543         zbetup_kp1 = 0._dp  
    544         zbetdo_kp1 = 0._dp 
    545  
    546         ! 3. monotonic flux in the i & j direction (paa & pbb) 
    547         ! ---------------------------------------- 
    548         monotonic_flux(zau,zbu,zcu,zbetup_ip1,zbetdo_ip1,paa,paa_in,jpk-1) 
    549         monotonic_flux(zav,zbv,zcv,zbetup_jp1,zbetdo_jp1,pbb,pbb_in,jpk-1) 
    550  
    551         ! monotonic flux in the k direction, i.e. pcc 
    552         ! ------------------------------------------- 
    553         monotonic_flux_k(za,zb,zc,zbetup_kp1,zbetdo_kp1,pcc,pcc_in,jpk-1) 
     549        monotonic_flux_k(za,zb,zc,0._dp,0._dp,pcc,pcc_in,jpk) 
    554550      END_2D 
    555551   END SUBROUTINE nonosc_lf 
Note: See TracChangeset for help on using the changeset viewer.