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 10957 for NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/BDY/bdydyn3d.F90 – NEMO

Ignore:
Timestamp:
2019-05-10T12:26:38+02:00 (5 years ago)
Author:
davestorkey
Message:

2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps : BDY module (including some removal of redundant code). Tested in AMM12 and ORCA1 (not full SETTE test at this stage).

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/BDY/bdydyn3d.F90

    r10529 r10957  
    3333CONTAINS 
    3434 
    35    SUBROUTINE bdy_dyn3d( kt ) 
     35   SUBROUTINE bdy_dyn3d( kt, Kbb, puu, pvv, Kaa ) 
    3636      !!---------------------------------------------------------------------- 
    3737      !!                  ***  SUBROUTINE bdy_dyn3d  *** 
     
    4040      !! 
    4141      !!---------------------------------------------------------------------- 
    42       INTEGER, INTENT(in) ::   kt   ! Main time step counter 
     42      INTEGER                             , INTENT( in    ) ::   kt        ! Main time step counter 
     43      INTEGER                             , INTENT( in    ) ::   Kbb, Kaa  ! Time level indices 
     44      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT( inout ) ::   puu, pvv  ! Ocean velocities (to be updated at open boundaries) 
    4345      ! 
    4446      INTEGER ::   ib_bdy   ! loop index 
     
    4951         SELECT CASE( cn_dyn3d(ib_bdy) ) 
    5052         CASE('none')        ;   CYCLE 
    51          CASE('frs' )        ;   CALL bdy_dyn3d_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy ) 
    52          CASE('specified')   ;   CALL bdy_dyn3d_spe( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy ) 
    53          CASE('zero')        ;   CALL bdy_dyn3d_zro( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy ) 
    54          CASE('orlanski' )   ;   CALL bdy_dyn3d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.false. ) 
    55          CASE('orlanski_npo');   CALL bdy_dyn3d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.true. ) 
    56          CASE('zerograd')    ;   CALL bdy_dyn3d_zgrad( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy ) 
    57          CASE('neumann')     ;   CALL bdy_dyn3d_nmn( idx_bdy(ib_bdy), ib_bdy ) 
     53         CASE('frs' )        ;   CALL bdy_dyn3d_frs(           puu, pvv, Kaa, idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy ) 
     54         CASE('specified')   ;   CALL bdy_dyn3d_spe(           puu, pvv, Kaa, idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy ) 
     55         CASE('zero')        ;   CALL bdy_dyn3d_zro(           puu, pvv, Kaa, idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy ) 
     56         CASE('orlanski' )   ;   CALL bdy_dyn3d_orlanski( Kbb, puu, pvv, Kaa, idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.false. ) 
     57         CASE('orlanski_npo');   CALL bdy_dyn3d_orlanski( Kbb, puu, pvv, Kaa, idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.true. ) 
     58         CASE('zerograd')    ;   CALL bdy_dyn3d_zgrad(         puu, pvv, Kaa, idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy ) 
     59         CASE('neumann')     ;   CALL bdy_dyn3d_nmn(           puu, pvv, Kaa, idx_bdy(ib_bdy), ib_bdy ) 
    5860         CASE DEFAULT        ;   CALL ctl_stop( 'bdy_dyn3d : unrecognised option for open boundaries for baroclinic velocities' ) 
    5961         END SELECT 
     
    6365 
    6466 
    65    SUBROUTINE bdy_dyn3d_spe( idx, dta, kt , ib_bdy ) 
     67   SUBROUTINE bdy_dyn3d_spe( puu, pvv, Kaa, idx, dta, ib_bdy ) 
    6668      !!---------------------------------------------------------------------- 
    6769      !!                  ***  SUBROUTINE bdy_dyn3d_spe  *** 
     
    7173      !! 
    7274      !!---------------------------------------------------------------------- 
    73       INTEGER        , INTENT(in) ::   kt      ! time step index 
    74       TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices 
    75       TYPE(OBC_DATA) , INTENT(in) ::   dta     ! OBC external data 
    76       INTEGER        , INTENT(in) ::   ib_bdy  ! BDY set index 
     75      INTEGER                             , INTENT( in    ) ::   Kaa       ! Time level index 
     76      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT( inout ) ::   puu, pvv  ! Ocean velocities (to be updated at open boundaries) 
     77      TYPE(OBC_INDEX)                     , INTENT( in    ) ::   idx       ! OBC indices 
     78      TYPE(OBC_DATA)                      , INTENT( in    ) ::   dta       ! OBC external data 
     79      INTEGER                             , INTENT( in    ) ::   ib_bdy    ! BDY set index 
    7780      ! 
    7881      INTEGER  ::   jb, jk         ! dummy loop indices 
     
    8689            ii   = idx%nbi(jb,igrd) 
    8790            ij   = idx%nbj(jb,igrd) 
    88             ua(ii,ij,jk) = dta%u3d(jb,jk) * umask(ii,ij,jk) 
     91            puu(ii,ij,jk,Kaa) = dta%u3d(jb,jk) * umask(ii,ij,jk) 
    8992         END DO 
    9093      END DO 
     
    9598            ii   = idx%nbi(jb,igrd) 
    9699            ij   = idx%nbj(jb,igrd) 
    97             va(ii,ij,jk) = dta%v3d(jb,jk) * vmask(ii,ij,jk) 
    98          END DO 
    99       END DO 
    100       CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )   ! Boundary points should be updated   
    101       CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1., ib_bdy )    
    102       ! 
    103       IF( kt == nit000 )   CLOSE( unit = 102 ) 
     100            pvv(ii,ij,jk,Kaa) = dta%v3d(jb,jk) * vmask(ii,ij,jk) 
     101         END DO 
     102      END DO 
     103      CALL lbc_bdy_lnk( 'bdydyn3d', puu(:,:,:,Kaa), 'U', -1., ib_bdy )   ! Boundary points should be updated   
     104      CALL lbc_bdy_lnk( 'bdydyn3d', pvv(:,:,:,Kaa), 'V', -1., ib_bdy )    
    104105      ! 
    105106   END SUBROUTINE bdy_dyn3d_spe 
    106107 
    107108 
    108    SUBROUTINE bdy_dyn3d_zgrad( idx, dta, kt , ib_bdy ) 
     109   SUBROUTINE bdy_dyn3d_zgrad( puu, pvv, Kaa, idx, dta, ib_bdy ) 
    109110      !!---------------------------------------------------------------------- 
    110111      !!                  ***  SUBROUTINE bdy_dyn3d_zgrad  *** 
     
    113114      !! 
    114115      !!---------------------------------------------------------------------- 
    115       INTEGER                     ::   kt 
    116       TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices 
    117       TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data 
    118       INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index 
     116      INTEGER                             , INTENT( in    ) ::   Kaa       ! Time level index 
     117      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT( inout ) ::   puu, pvv  ! Ocean velocities (to be updated at open boundaries) 
     118      TYPE(OBC_INDEX)                     , INTENT( in    ) ::   idx       ! OBC indices 
     119      TYPE(OBC_DATA)                      , INTENT( in    ) ::   dta       ! OBC external data 
     120      INTEGER                             , INTENT( in    ) ::   ib_bdy    ! BDY set index 
    119121      !! 
    120122      INTEGER  ::   jb, jk         ! dummy loop indices 
     
    130132            ij   = idx%nbj(jb,igrd) 
    131133            fu   = ABS( ABS (NINT( idx%flagu(jb,igrd) ) ) - 1 ) 
    132             ua(ii,ij,jk) = ua(ii,ij,jk) * REAL( 1 - fu) + ( ua(ii,ij+fu,jk) * umask(ii,ij+fu,jk) & 
    133                         &+ ua(ii,ij-fu,jk) * umask(ii,ij-fu,jk) ) * umask(ii,ij,jk) * REAL( fu ) 
     134            puu(ii,ij,jk,Kaa) = puu(ii,ij,jk,Kaa) * REAL( 1 - fu) + ( puu(ii,ij+fu,jk,Kaa) * umask(ii,ij+fu,jk) & 
     135                        &+ puu(ii,ij-fu,jk,Kaa) * umask(ii,ij-fu,jk) ) * umask(ii,ij,jk) * REAL( fu ) 
    134136         END DO 
    135137      END DO 
     
    141143            ij   = idx%nbj(jb,igrd) 
    142144            fv   = ABS( ABS (NINT( idx%flagv(jb,igrd) ) ) - 1 ) 
    143             va(ii,ij,jk) = va(ii,ij,jk) * REAL( 1 - fv ) + ( va(ii+fv,ij,jk) * vmask(ii+fv,ij,jk) & 
    144                         &+ va(ii-fv,ij,jk) * vmask(ii-fv,ij,jk) ) * vmask(ii,ij,jk) * REAL( fv ) 
    145          END DO 
    146       END DO 
    147       CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )   ! Boundary points should be updated   
    148       CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1., ib_bdy )    
    149       ! 
    150       IF( kt == nit000 )   CLOSE( unit = 102 ) 
     145            pvv(ii,ij,jk,Kaa) = pvv(ii,ij,jk,Kaa) * REAL( 1 - fv ) + ( pvv(ii+fv,ij,jk,Kaa) * vmask(ii+fv,ij,jk) & 
     146                        &+ pvv(ii-fv,ij,jk,Kaa) * vmask(ii-fv,ij,jk) ) * vmask(ii,ij,jk) * REAL( fv ) 
     147         END DO 
     148      END DO 
     149      CALL lbc_bdy_lnk( 'bdydyn3d', puu(:,:,:,Kaa), 'U', -1., ib_bdy )   ! Boundary points should be updated   
     150      CALL lbc_bdy_lnk( 'bdydyn3d', pvv(:,:,:,Kaa), 'V', -1., ib_bdy )    
    151151      ! 
    152152   END SUBROUTINE bdy_dyn3d_zgrad 
    153153 
    154154 
    155    SUBROUTINE bdy_dyn3d_zro( idx, dta, kt, ib_bdy ) 
     155   SUBROUTINE bdy_dyn3d_zro( puu, pvv, Kaa, idx, dta, ib_bdy ) 
    156156      !!---------------------------------------------------------------------- 
    157157      !!                  ***  SUBROUTINE bdy_dyn3d_zro  *** 
     
    160160      !! 
    161161      !!---------------------------------------------------------------------- 
    162       INTEGER        , INTENT(in) ::   kt      ! time step index 
    163       TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices 
    164       TYPE(OBC_DATA) , INTENT(in) ::   dta     ! OBC external data 
    165       INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index 
     162      INTEGER                             , INTENT( in    ) ::   Kaa       ! Time level index 
     163      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT( inout ) ::   puu, pvv  ! Ocean velocities (to be updated at open boundaries) 
     164      TYPE(OBC_INDEX)                     , INTENT( in    ) ::   idx       ! OBC indices 
     165      TYPE(OBC_DATA)                      , INTENT( in    ) ::   dta       ! OBC external data 
     166      INTEGER                             , INTENT( in    ) ::   ib_bdy    ! BDY set index 
    166167      ! 
    167168      INTEGER  ::   ib, ik         ! dummy loop indices 
     
    175176         ij = idx%nbj(ib,igrd) 
    176177         DO ik = 1, jpkm1 
    177             ua(ii,ij,ik) = 0._wp 
     178            puu(ii,ij,ik,Kaa) = 0._wp 
    178179         END DO 
    179180      END DO 
     
    184185         ij = idx%nbj(ib,igrd) 
    185186         DO ik = 1, jpkm1 
    186             va(ii,ij,ik) = 0._wp 
    187          END DO 
    188       END DO 
    189       ! 
    190       CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )   ;   CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1.,ib_bdy )   ! Boundary points should be updated 
    191       ! 
    192       IF( kt == nit000 )   CLOSE( unit = 102 ) 
     187            pvv(ii,ij,ik,Kaa) = 0._wp 
     188         END DO 
     189      END DO 
     190      ! 
     191      CALL lbc_bdy_lnk( 'bdydyn3d', puu(:,:,:,Kaa), 'U', -1., ib_bdy )   ;   CALL lbc_bdy_lnk( 'bdydyn3d', pvv(:,:,:,Kaa), 'V', -1.,ib_bdy )   ! Boundary points should be updated 
    193192      ! 
    194193   END SUBROUTINE bdy_dyn3d_zro 
    195194 
    196195 
    197    SUBROUTINE bdy_dyn3d_frs( idx, dta, kt, ib_bdy ) 
     196   SUBROUTINE bdy_dyn3d_frs( puu, pvv, Kaa, idx, dta, ib_bdy ) 
    198197      !!---------------------------------------------------------------------- 
    199198      !!                  ***  SUBROUTINE bdy_dyn3d_frs  *** 
     
    206205      !!               topography. Tellus, 365-382. 
    207206      !!---------------------------------------------------------------------- 
    208       INTEGER        , INTENT(in) ::   kt      ! time step index 
    209       TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices 
    210       TYPE(OBC_DATA) , INTENT(in) ::   dta     ! OBC external data 
    211       INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index 
     207      INTEGER                             , INTENT( in    ) ::   Kaa       ! Time level index 
     208      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT( inout ) ::   puu, pvv  ! Ocean velocities (to be updated at open boundaries) 
     209      TYPE(OBC_INDEX)                     , INTENT( in    ) ::   idx       ! OBC indices 
     210      TYPE(OBC_DATA)                      , INTENT( in    ) ::   dta       ! OBC external data 
     211      INTEGER                             , INTENT( in    ) ::   ib_bdy    ! BDY set index 
    212212      ! 
    213213      INTEGER  ::   jb, jk         ! dummy loop indices 
     
    222222            ij   = idx%nbj(jb,igrd) 
    223223            zwgt = idx%nbw(jb,igrd) 
    224             ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta%u3d(jb,jk) - ua(ii,ij,jk) ) ) * umask(ii,ij,jk) 
     224            puu(ii,ij,jk,Kaa) = ( puu(ii,ij,jk,Kaa) + zwgt * ( dta%u3d(jb,jk) - puu(ii,ij,jk,Kaa) ) ) * umask(ii,ij,jk) 
    225225         END DO 
    226226      END DO 
     
    232232            ij   = idx%nbj(jb,igrd) 
    233233            zwgt = idx%nbw(jb,igrd) 
    234             va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta%v3d(jb,jk) - va(ii,ij,jk) ) ) * vmask(ii,ij,jk) 
     234            pvv(ii,ij,jk,Kaa) = ( pvv(ii,ij,jk,Kaa) + zwgt * ( dta%v3d(jb,jk) - pvv(ii,ij,jk,Kaa) ) ) * vmask(ii,ij,jk) 
    235235         END DO 
    236236      END DO  
    237       CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )    ! Boundary points should be updated 
    238       CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1., ib_bdy )    
    239       ! 
    240       IF( kt == nit000 )   CLOSE( unit = 102 ) 
     237      CALL lbc_bdy_lnk( 'bdydyn3d', puu(:,:,:,Kaa), 'U', -1., ib_bdy )    ! Boundary points should be updated 
     238      CALL lbc_bdy_lnk( 'bdydyn3d', pvv(:,:,:,Kaa), 'V', -1., ib_bdy )    
    241239      ! 
    242240   END SUBROUTINE bdy_dyn3d_frs 
    243241 
    244242 
    245    SUBROUTINE bdy_dyn3d_orlanski( idx, dta, ib_bdy, ll_npo ) 
     243   SUBROUTINE bdy_dyn3d_orlanski( Kbb, puu, pvv, Kaa, idx, dta, ib_bdy, ll_npo ) 
    246244      !!---------------------------------------------------------------------- 
    247245      !!                 ***  SUBROUTINE bdy_dyn3d_orlanski  *** 
     
    253251      !! References:  Marchesiello, McWilliams and Shchepetkin, Ocean Modelling vol. 3 (2001)     
    254252      !!---------------------------------------------------------------------- 
    255       TYPE(OBC_INDEX),              INTENT(in) ::   idx  ! OBC indices 
    256       TYPE(OBC_DATA),               INTENT(in) ::   dta  ! OBC external data 
    257       INTEGER,                      INTENT(in) ::   ib_bdy  ! BDY set index 
    258       LOGICAL,                      INTENT(in) ::   ll_npo  ! switch for NPO version 
     253      INTEGER                             , INTENT( in    ) ::   Kbb, Kaa  ! Time level indices 
     254      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT( inout ) ::   puu, pvv  ! Ocean velocities (to be updated at open boundaries) 
     255      TYPE(OBC_INDEX)                     , INTENT( in    ) ::   idx  ! OBC indices 
     256      TYPE(OBC_DATA)                      , INTENT( in    ) ::   dta  ! OBC external data 
     257      INTEGER                             , INTENT( in    ) ::   ib_bdy  ! BDY set index 
     258      LOGICAL                             , INTENT( in    ) ::   ll_npo  ! switch for NPO version 
    259259 
    260260      INTEGER  ::   jb, igrd                               ! dummy loop indices 
    261261      !!---------------------------------------------------------------------- 
    262262      ! 
    263       !! Note that at this stage the ub and ua arrays contain the baroclinic velocities.  
     263      !! Note that at this stage the puu(:,:,:,Kbb) and puu(:,:,:,Kaa) arrays contain the baroclinic velocities.  
    264264      ! 
    265265      igrd = 2      ! Orlanski bc on u-velocity;  
    266266      !             
    267       CALL bdy_orlanski_3d( idx, igrd, ub, ua, dta%u3d, ll_npo ) 
     267      CALL bdy_orlanski_3d( idx, igrd, puu(:,:,:,Kbb), puu(:,:,:,Kaa), dta%u3d, ll_npo ) 
    268268 
    269269      igrd = 3      ! Orlanski bc on v-velocity 
    270270      !   
    271       CALL bdy_orlanski_3d( idx, igrd, vb, va, dta%v3d, ll_npo ) 
    272       ! 
    273       CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )    ! Boundary points should be updated 
    274       CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1., ib_bdy )    
     271      CALL bdy_orlanski_3d( idx, igrd, pvv(:,:,:,Kbb), pvv(:,:,:,Kaa), dta%v3d, ll_npo ) 
     272      ! 
     273      CALL lbc_bdy_lnk( 'bdydyn3d', puu(:,:,:,Kaa), 'U', -1., ib_bdy )    ! Boundary points should be updated 
     274      CALL lbc_bdy_lnk( 'bdydyn3d', pvv(:,:,:,Kaa), 'V', -1., ib_bdy )    
    275275      ! 
    276276   END SUBROUTINE bdy_dyn3d_orlanski 
    277277 
    278278 
    279    SUBROUTINE bdy_dyn3d_dmp( kt ) 
     279   SUBROUTINE bdy_dyn3d_dmp( kt, Kbb, puu, pvv, Krhs ) 
    280280      !!---------------------------------------------------------------------- 
    281281      !!                  ***  SUBROUTINE bdy_dyn3d_dmp  *** 
     
    284284      !! 
    285285      !!---------------------------------------------------------------------- 
    286       INTEGER, INTENT(in) ::   kt   ! time step index 
     286      INTEGER                             , INTENT( in    ) ::   kt        ! time step 
     287      INTEGER                             , INTENT( in    ) ::   Kbb, Krhs ! Time level indices 
     288      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT( inout ) ::   puu, pvv  ! Ocean velocities and trends (to be updated at open boundaries) 
    287289      ! 
    288290      INTEGER  ::   jb, jk         ! dummy loop indices 
     
    302304               zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd) 
    303305               DO jk = 1, jpkm1 
    304                   ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%u3d(jb,jk) - & 
    305                                    ub(ii,ij,jk) + ub_b(ii,ij)) ) * umask(ii,ij,jk) 
     306                  puu(ii,ij,jk,Krhs) = ( puu(ii,ij,jk,Krhs) + zwgt * ( dta_bdy(ib_bdy)%u3d(jb,jk) - & 
     307                                   puu(ii,ij,jk,Kbb) + uu_b(ii,ij,Kbb)) ) * umask(ii,ij,jk) 
    306308               END DO 
    307309            END DO 
     
    313315               zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd) 
    314316               DO jk = 1, jpkm1 
    315                   va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%v3d(jb,jk) -  & 
    316                                    vb(ii,ij,jk) + vb_b(ii,ij)) ) * vmask(ii,ij,jk) 
     317                  pvv(ii,ij,jk,Krhs) = ( pvv(ii,ij,jk,Krhs) + zwgt * ( dta_bdy(ib_bdy)%v3d(jb,jk) -  & 
     318                                   pvv(ii,ij,jk,Kbb) + vv_b(ii,ij,Kbb)) ) * vmask(ii,ij,jk) 
    317319               END DO 
    318320            END DO 
     
    320322      END DO 
    321323      ! 
    322       CALL lbc_lnk_multi( 'bdydyn3d', ua, 'U', -1.,  va, 'V', -1. )   ! Boundary points should be updated 
     324      CALL lbc_lnk_multi( 'bdydyn3d', puu(:,:,:,Krhs), 'U', -1.,  pvv(:,:,:,Krhs), 'V', -1. )   ! Boundary points should be updated 
    323325      ! 
    324326      IF( ln_timing )   CALL timing_stop('bdy_dyn3d_dmp') 
     
    327329 
    328330 
    329    SUBROUTINE bdy_dyn3d_nmn( idx, ib_bdy ) 
     331   SUBROUTINE bdy_dyn3d_nmn( puu, pvv, Kaa, idx, ib_bdy ) 
    330332      !!---------------------------------------------------------------------- 
    331333      !!                 ***  SUBROUTINE bdy_dyn3d_nmn  *** 
     
    336338      !! 
    337339      !!---------------------------------------------------------------------- 
    338       TYPE(OBC_INDEX),              INTENT(in) ::   idx  ! OBC indices 
    339       INTEGER,                      INTENT(in) ::   ib_bdy  ! BDY set index 
     340      INTEGER                             , INTENT( in    ) ::   Kaa       ! Time level index 
     341      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT( inout ) ::   puu, pvv  ! Ocean velocities (to be updated at open boundaries) 
     342      TYPE(OBC_INDEX)                     , INTENT( in    ) ::   idx       ! OBC indices 
     343      INTEGER                             , INTENT( in    ) ::   ib_bdy    ! BDY set index 
    340344 
    341345      INTEGER  ::   jb, igrd                               ! dummy loop indices 
    342346      !!---------------------------------------------------------------------- 
    343347      ! 
    344       !! Note that at this stage the ub and ua arrays contain the baroclinic velocities.  
     348      !! Note that at this stage the puu(:,:,:,Kbb) and puu(:,:,:,Kaa) arrays contain the baroclinic velocities.  
    345349      ! 
    346350      igrd = 2      ! Neumann bc on u-velocity;  
    347351      !             
    348       CALL bdy_nmn( idx, igrd, ua ) 
     352      CALL bdy_nmn( idx, igrd, puu(:,:,:,Kaa) ) 
    349353 
    350354      igrd = 3      ! Neumann bc on v-velocity 
    351355      !   
    352       CALL bdy_nmn( idx, igrd, va ) 
    353       ! 
    354       CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )    ! Boundary points should be updated 
    355       CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1., ib_bdy ) 
     356      CALL bdy_nmn( idx, igrd, pvv(:,:,:,Kaa) ) 
     357      ! 
     358      CALL lbc_bdy_lnk( 'bdydyn3d', puu(:,:,:,Kaa), 'U', -1., ib_bdy )    ! Boundary points should be updated 
     359      CALL lbc_bdy_lnk( 'bdydyn3d', pvv(:,:,:,Kaa), 'V', -1., ib_bdy ) 
    356360      ! 
    357361   END SUBROUTINE bdy_dyn3d_nmn 
Note: See TracChangeset for help on using the changeset viewer.