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 – 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).

Location:
NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/BDY
Files:
5 edited

Legend:

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

    r10647 r10957  
    6565CONTAINS 
    6666 
    67       SUBROUTINE bdy_dta( kt, jit, time_offset ) 
     67      SUBROUTINE bdy_dta( kt, Kmm, time_offset ) 
    6868      !!---------------------------------------------------------------------- 
    6969      !!                   ***  SUBROUTINE bdy_dta  *** 
     
    7575      !!---------------------------------------------------------------------- 
    7676      INTEGER, INTENT(in)           ::   kt           ! ocean time-step index  
    77       INTEGER, INTENT(in), OPTIONAL ::   jit          ! subcycle time-step index (for timesplitting option) 
     77      INTEGER, INTENT(in)           ::   Kmm          ! ocean time level index 
    7878      INTEGER, INTENT(in), OPTIONAL ::   time_offset  ! time offset in units of timesteps. NB. if jit 
    7979      !                                               ! is present then units = subcycle timesteps. 
     
    9494      ! Initialise data arrays once for all from initial conditions where required 
    9595      !--------------------------------------------------------------------------- 
    96       IF( kt == nit000 .AND. .NOT.PRESENT(jit) ) THEN 
     96      IF( kt == nit000 ) THEN 
    9797 
    9898         ! Calculate depth-mean currents 
     
    112112                     ii = idx_bdy(jbdy)%nbi(ib,igrd) 
    113113                     ij = idx_bdy(jbdy)%nbj(ib,igrd) 
    114                      dta_bdy(jbdy)%ssh(ib) = sshn(ii,ij) * tmask(ii,ij,1)          
     114                     dta_bdy(jbdy)%ssh(ib) = ssh(ii,ij,Kmm) * tmask(ii,ij,1)          
    115115                  END DO  
    116116               ENDIF 
     
    120120                     ii = idx_bdy(jbdy)%nbi(ib,igrd) 
    121121                     ij = idx_bdy(jbdy)%nbj(ib,igrd) 
    122                      dta_bdy(jbdy)%u2d(ib) = un_b(ii,ij) * umask(ii,ij,1)          
     122                     dta_bdy(jbdy)%u2d(ib) = uu_b(ii,ij,Kmm) * umask(ii,ij,1)          
    123123                  END DO  
    124124               ENDIF 
     
    128128                     ii = idx_bdy(jbdy)%nbi(ib,igrd) 
    129129                     ij = idx_bdy(jbdy)%nbj(ib,igrd) 
    130                      dta_bdy(jbdy)%v2d(ib) = vn_b(ii,ij) * vmask(ii,ij,1)          
     130                     dta_bdy(jbdy)%v2d(ib) = vv_b(ii,ij,Kmm) * vmask(ii,ij,1)          
    131131                  END DO  
    132132               ENDIF 
     
    141141                        ii = idx_bdy(jbdy)%nbi(ib,igrd) 
    142142                        ij = idx_bdy(jbdy)%nbj(ib,igrd) 
    143                         dta_bdy(jbdy)%u3d(ib,ik) =  ( un(ii,ij,ik) - un_b(ii,ij) ) * umask(ii,ij,ik)          
     143                        dta_bdy(jbdy)%u3d(ib,ik) =  ( uu(ii,ij,ik,Kmm) - uu_b(ii,ij,Kmm) ) * umask(ii,ij,ik)          
    144144                     END DO 
    145145                  END DO  
     
    151151                        ii = idx_bdy(jbdy)%nbi(ib,igrd) 
    152152                        ij = idx_bdy(jbdy)%nbj(ib,igrd) 
    153                         dta_bdy(jbdy)%v3d(ib,ik) =  ( vn(ii,ij,ik) - vn_b(ii,ij) ) * vmask(ii,ij,ik)          
     153                        dta_bdy(jbdy)%v3d(ib,ik) =  ( vv(ii,ij,ik,Kmm) - vv_b(ii,ij,Kmm) ) * vmask(ii,ij,ik)          
    154154                        END DO 
    155155                  END DO  
     
    165165                        ii = idx_bdy(jbdy)%nbi(ib,igrd) 
    166166                        ij = idx_bdy(jbdy)%nbj(ib,igrd) 
    167                         dta_bdy(jbdy)%tem(ib,ik) = tsn(ii,ij,ik,jp_tem) * tmask(ii,ij,ik)          
     167                        dta_bdy(jbdy)%tem(ib,ik) = ts(ii,ij,ik,jp_tem,Kmm) * tmask(ii,ij,ik)          
    168168                     END DO 
    169169                  END DO  
     
    175175                        ii = idx_bdy(jbdy)%nbi(ib,igrd) 
    176176                        ij = idx_bdy(jbdy)%nbj(ib,igrd) 
    177                         dta_bdy(jbdy)%sal(ib,ik) = tsn(ii,ij,ik,jp_sal) * tmask(ii,ij,ik)          
     177                        dta_bdy(jbdy)%sal(ib,ik) = ts(ii,ij,ik,jp_sal,Kmm) * tmask(ii,ij,ik)          
    178178                     END DO 
    179179                  END DO  
     
    227227         dta => dta_bdy(jbdy) 
    228228         IF( nn_dta(jbdy) == 1 ) THEN ! skip this bit if no external data required 
    229        
    230             IF( PRESENT(jit) ) THEN 
    231                ! Update barotropic boundary conditions only 
    232                ! jit is optional argument for fld_read and bdytide_update 
    233                IF( cn_dyn2d(jbdy) /= 'none' ) THEN 
    234                   IF( nn_dyn2d_dta(jbdy) == 2 ) THEN ! tidal harmonic forcing ONLY: initialise arrays 
    235                      IF( dta%ll_ssh ) dta%ssh(:) = 0._wp 
    236                      IF( dta%ll_u2d ) dta%u2d(:) = 0._wp 
    237                      IF( dta%ll_u3d ) dta%v2d(:) = 0._wp 
    238                   ENDIF 
    239                   IF (cn_tra(jbdy) /= 'runoff') THEN 
    240                      IF( nn_dyn2d_dta(jbdy) == 1 .OR. nn_dyn2d_dta(jbdy) == 3 ) THEN 
    241  
    242                         jend = jstart + dta%nread(2) - 1 
    243                         IF( ln_full_vel_array(jbdy) ) THEN 
    244                            CALL fld_read( kt=kt, kn_fsbc=1, sd=bf(jstart:jend), map=nbmap_ptr(jstart:jend),  & 
    245                                      & kit=jit, kt_offset=time_offset , jpk_bdy=nb_jpk_bdy,   & 
    246                                      & fvl=ln_full_vel_array(jbdy)  ) 
    247                         ELSE 
    248                            CALL fld_read( kt=kt, kn_fsbc=1, sd=bf(jstart:jend), map=nbmap_ptr(jstart:jend),  & 
    249                                      & kit=jit, kt_offset=time_offset  ) 
    250                         ENDIF 
    251  
    252                         ! If full velocities in boundary data then extract barotropic velocities from 3D fields 
    253                         IF( ln_full_vel_array(jbdy) .AND.                                             & 
    254                           &    ( nn_dyn2d_dta(jbdy) == 1 .OR. nn_dyn2d_dta(jbdy) == 3 .OR.  & 
    255                           &      nn_dyn3d_dta(jbdy) == 1 ) )THEN 
    256  
    257                            igrd = 2                      ! zonal velocity 
    258                            dta%u2d(:) = 0._wp 
    259                            DO ib = 1, idx_bdy(jbdy)%nblen(igrd) 
    260                               ii   = idx_bdy(jbdy)%nbi(ib,igrd) 
    261                               ij   = idx_bdy(jbdy)%nbj(ib,igrd) 
    262                               DO ik = 1, jpkm1 
    263                                  dta%u2d(ib) = dta%u2d(ib) & 
    264                        &                          + e3u_n(ii,ij,ik) * umask(ii,ij,ik) * dta%u3d(ib,ik) 
    265                               END DO 
    266                               dta%u2d(ib) =  dta%u2d(ib) * r1_hu_n(ii,ij) 
    267                            END DO 
    268                            igrd = 3                      ! meridional velocity 
    269                            dta%v2d(:) = 0._wp 
    270                            DO ib = 1, idx_bdy(jbdy)%nblen(igrd) 
    271                               ii   = idx_bdy(jbdy)%nbi(ib,igrd) 
    272                               ij   = idx_bdy(jbdy)%nbj(ib,igrd) 
    273                               DO ik = 1, jpkm1 
    274                                  dta%v2d(ib) = dta%v2d(ib) & 
    275                        &                       + e3v_n(ii,ij,ik) * vmask(ii,ij,ik) * dta%v3d(ib,ik) 
    276                               END DO 
    277                               dta%v2d(ib) =  dta%v2d(ib) * r1_hv_n(ii,ij) 
    278                            END DO 
    279                         ENDIF                     
    280                      ENDIF 
    281                      IF( nn_dyn2d_dta(jbdy) .ge. 2 ) THEN ! update tidal harmonic forcing 
    282                         CALL bdytide_update( kt=kt, idx=idx_bdy(jbdy), dta=dta, td=tides(jbdy),   &  
    283                           &                 jit=jit, time_offset=time_offset ) 
    284                      ENDIF 
    285                   ENDIF 
    286                ENDIF 
     229            IF (cn_tra(jbdy) == 'runoff') then      ! runoff condition 
     230               jend = nb_bdy_fld(jbdy) 
     231               CALL fld_read( kt=kt, kn_fsbc=1, sd=bf(jstart:jend),  & 
     232                    & map=nbmap_ptr(jstart:jend), kt_offset=time_offset ) 
     233               ! 
     234               igrd = 2                      ! zonal velocity 
     235               DO ib = 1, idx_bdy(jbdy)%nblen(igrd) 
     236                  ii   = idx_bdy(jbdy)%nbi(ib,igrd) 
     237                  ij   = idx_bdy(jbdy)%nbj(ib,igrd) 
     238                  dta%u2d(ib) = dta%u2d(ib) / ( e2u(ii,ij) * hu_0(ii,ij) ) 
     239               END DO 
     240               ! 
     241               igrd = 3                      ! meridional velocity 
     242               DO ib = 1, idx_bdy(jbdy)%nblen(igrd) 
     243                  ii   = idx_bdy(jbdy)%nbi(ib,igrd) 
     244                  ij   = idx_bdy(jbdy)%nbj(ib,igrd) 
     245                  dta%v2d(ib) = dta%v2d(ib) / ( e1v(ii,ij) * hv_0(ii,ij) ) 
     246               END DO 
    287247            ELSE 
    288                IF (cn_tra(jbdy) == 'runoff') then      ! runoff condition 
    289                   jend = nb_bdy_fld(jbdy) 
    290                   CALL fld_read( kt=kt, kn_fsbc=1, sd=bf(jstart:jend),  & 
    291                                & map=nbmap_ptr(jstart:jend), kt_offset=time_offset ) 
    292                   ! 
     248               IF( nn_dyn2d_dta(jbdy) == 2 ) THEN ! tidal harmonic forcing ONLY: initialise arrays 
     249                  IF( dta%ll_ssh ) dta%ssh(:) = 0._wp 
     250                  IF( dta%ll_u2d ) dta%u2d(:) = 0._wp 
     251                  IF( dta%ll_v2d ) dta%v2d(:) = 0._wp 
     252               ENDIF 
     253               IF( dta%nread(1) .gt. 0 ) THEN ! update external data 
     254                  jend = jstart + dta%nread(1) - 1 
     255                  CALL fld_read( kt=kt, kn_fsbc=1, sd=bf(jstart:jend), & 
     256                       & map=nbmap_ptr(jstart:jend), kt_offset=time_offset, jpk_bdy=nb_jpk_bdy,   & 
     257                       & fvl=ln_full_vel_array(jbdy) ) 
     258               ENDIF 
     259               ! If full velocities in boundary data then split into barotropic and baroclinic data 
     260               IF( ln_full_vel_array(jbdy) .and.                                             & 
     261                    & ( nn_dyn2d_dta(jbdy) == 1 .OR. nn_dyn2d_dta(jbdy) == 3 .OR. & 
     262                    &   nn_dyn3d_dta(jbdy) == 1 ) ) THEN 
    293263                  igrd = 2                      ! zonal velocity 
     264                  dta%u2d(:) = 0._wp 
    294265                  DO ib = 1, idx_bdy(jbdy)%nblen(igrd) 
    295266                     ii   = idx_bdy(jbdy)%nbi(ib,igrd) 
    296267                     ij   = idx_bdy(jbdy)%nbj(ib,igrd) 
    297                      dta%u2d(ib) = dta%u2d(ib) / ( e2u(ii,ij) * hu_0(ii,ij) ) 
     268                     DO ik = 1, jpkm1 
     269                        dta%u2d(ib) = dta%u2d(ib) & 
     270                             &                       + e3u(ii,ij,ik,Kmm) * umask(ii,ij,ik) * dta%u3d(ib,ik) 
     271                     END DO 
     272                     dta%u2d(ib) =  dta%u2d(ib) * r1_hu_n(ii,ij) 
     273                     DO ik = 1, jpkm1 
     274                        dta%u3d(ib,ik) = dta%u3d(ib,ik) - dta%u2d(ib) 
     275                     END DO 
    298276                  END DO 
    299                   ! 
    300277                  igrd = 3                      ! meridional velocity 
     278                  dta%v2d(:) = 0._wp 
    301279                  DO ib = 1, idx_bdy(jbdy)%nblen(igrd) 
    302280                     ii   = idx_bdy(jbdy)%nbi(ib,igrd) 
    303281                     ij   = idx_bdy(jbdy)%nbj(ib,igrd) 
    304                      dta%v2d(ib) = dta%v2d(ib) / ( e1v(ii,ij) * hv_0(ii,ij) ) 
     282                     DO ik = 1, jpkm1 
     283                        dta%v2d(ib) = dta%v2d(ib) & 
     284                             &                       + e3v(ii,ij,ik,Kmm) * vmask(ii,ij,ik) * dta%v3d(ib,ik) 
     285                     END DO 
     286                     dta%v2d(ib) =  dta%v2d(ib) * r1_hv_n(ii,ij) 
     287                     DO ik = 1, jpkm1 
     288                        dta%v3d(ib,ik) = dta%v3d(ib,ik) - dta%v2d(ib) 
     289                     END DO 
    305290                  END DO 
    306                ELSE 
    307                   IF( nn_dyn2d_dta(jbdy) == 2 ) THEN ! tidal harmonic forcing ONLY: initialise arrays 
    308                      IF( dta%ll_ssh ) dta%ssh(:) = 0._wp 
    309                      IF( dta%ll_u2d ) dta%u2d(:) = 0._wp 
    310                      IF( dta%ll_v2d ) dta%v2d(:) = 0._wp 
    311                   ENDIF 
    312                   IF( dta%nread(1) .gt. 0 ) THEN ! update external data 
    313                      jend = jstart + dta%nread(1) - 1 
    314                      CALL fld_read( kt=kt, kn_fsbc=1, sd=bf(jstart:jend), & 
    315                                   & map=nbmap_ptr(jstart:jend), kt_offset=time_offset, jpk_bdy=nb_jpk_bdy,   & 
    316                                   & fvl=ln_full_vel_array(jbdy) ) 
    317                   ENDIF 
    318                   ! If full velocities in boundary data then split into barotropic and baroclinic data 
    319                   IF( ln_full_vel_array(jbdy) .and.                                             & 
    320                     & ( nn_dyn2d_dta(jbdy) == 1 .OR. nn_dyn2d_dta(jbdy) == 3 .OR. & 
    321                     &   nn_dyn3d_dta(jbdy) == 1 ) ) THEN 
    322                      igrd = 2                      ! zonal velocity 
    323                      dta%u2d(:) = 0._wp 
    324                      DO ib = 1, idx_bdy(jbdy)%nblen(igrd) 
    325                         ii   = idx_bdy(jbdy)%nbi(ib,igrd) 
    326                         ij   = idx_bdy(jbdy)%nbj(ib,igrd) 
    327                         DO ik = 1, jpkm1 
    328                            dta%u2d(ib) = dta%u2d(ib) & 
    329                  &                       + e3u_n(ii,ij,ik) * umask(ii,ij,ik) * dta%u3d(ib,ik) 
    330                         END DO 
    331                         dta%u2d(ib) =  dta%u2d(ib) * r1_hu_n(ii,ij) 
    332                         DO ik = 1, jpkm1 
    333                            dta%u3d(ib,ik) = dta%u3d(ib,ik) - dta%u2d(ib) 
    334                         END DO 
    335                      END DO 
    336                      igrd = 3                      ! meridional velocity 
    337                      dta%v2d(:) = 0._wp 
    338                      DO ib = 1, idx_bdy(jbdy)%nblen(igrd) 
    339                         ii   = idx_bdy(jbdy)%nbi(ib,igrd) 
    340                         ij   = idx_bdy(jbdy)%nbj(ib,igrd) 
    341                         DO ik = 1, jpkm1 
    342                            dta%v2d(ib) = dta%v2d(ib) & 
    343                  &                       + e3v_n(ii,ij,ik) * vmask(ii,ij,ik) * dta%v3d(ib,ik) 
    344                         END DO 
    345                         dta%v2d(ib) =  dta%v2d(ib) * r1_hv_n(ii,ij) 
    346                         DO ik = 1, jpkm1 
    347                            dta%v3d(ib,ik) = dta%v3d(ib,ik) - dta%v2d(ib) 
    348                         END DO 
    349                      END DO 
    350                   ENDIF 
    351  
    352                ENDIF 
    353 #if defined key_si3 
    354                ! convert N-cat fields (input) into jpl-cat (output) 
    355                IF( cn_ice(jbdy) /= 'none' .AND. nn_ice_dta(jbdy) == 1 ) THEN 
    356                   jfld_hti = jfld_htit(jbdy) 
    357                   jfld_hts = jfld_htst(jbdy) 
    358                   jfld_ai  = jfld_ait(jbdy) 
    359                   IF    ( jpl /= 1 .AND. nice_cat == 1 ) THEN                       ! case input cat = 1 
    360                      CALL ice_var_itd ( bf(jfld_hti)%fnow(:,1,1), bf(jfld_hts)%fnow(:,1,1), bf(jfld_ai)%fnow(:,1,1), & 
    361                         &               dta_bdy(jbdy)%h_i     , dta_bdy(jbdy)%h_s     , dta_bdy(jbdy)%a_i    ) 
    362                   ELSEIF( jpl /= 1 .AND. nice_cat /= 1 .AND. nice_cat /= jpl ) THEN ! case input cat /=1 and /=jpl 
    363                      CALL ice_var_itd2( bf(jfld_hti)%fnow(:,1,:), bf(jfld_hts)%fnow(:,1,:), bf(jfld_ai)%fnow(:,1,:), & 
    364                         &               dta_bdy(jbdy)%h_i     , dta_bdy(jbdy)%h_s     , dta_bdy(jbdy)%a_i    ) 
    365                   ENDIF 
    366                ENDIF 
    367 #endif 
    368             ENDIF 
     291               ENDIF 
     292 
     293            ENDIF 
     294#if defined key_si3 
     295            ! convert N-cat fields (input) into jpl-cat (output) 
     296            IF( cn_ice(jbdy) /= 'none' .AND. nn_ice_dta(jbdy) == 1 ) THEN 
     297               jfld_hti = jfld_htit(jbdy) 
     298               jfld_hts = jfld_htst(jbdy) 
     299               jfld_ai  = jfld_ait(jbdy) 
     300               IF    ( jpl /= 1 .AND. nice_cat == 1 ) THEN                       ! case input cat = 1 
     301                  CALL ice_var_itd ( bf(jfld_hti)%fnow(:,1,1), bf(jfld_hts)%fnow(:,1,1), bf(jfld_ai)%fnow(:,1,1), & 
     302                       &               dta_bdy(jbdy)%h_i     , dta_bdy(jbdy)%h_s     , dta_bdy(jbdy)%a_i    ) 
     303               ELSEIF( jpl /= 1 .AND. nice_cat /= 1 .AND. nice_cat /= jpl ) THEN ! case input cat /=1 and /=jpl 
     304                  CALL ice_var_itd2( bf(jfld_hti)%fnow(:,1,:), bf(jfld_hts)%fnow(:,1,:), bf(jfld_ai)%fnow(:,1,:), & 
     305                       &               dta_bdy(jbdy)%h_i     , dta_bdy(jbdy)%h_s     , dta_bdy(jbdy)%a_i    ) 
     306               ENDIF 
     307            ENDIF 
     308#endif 
    369309            jstart = jstart + dta%nread(1) 
    370310         ENDIF    ! nn_dta(jbdy) = 1 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/BDY/bdydyn.F90

    r10068 r10957  
    3737CONTAINS 
    3838 
    39    SUBROUTINE bdy_dyn( kt, dyn3d_only ) 
     39   SUBROUTINE bdy_dyn( kt, Kbb, puu, pvv, Kaa, dyn3d_only ) 
    4040      !!---------------------------------------------------------------------- 
    4141      !!                  ***  SUBROUTINE bdy_dyn  *** 
     
    4444      !! 
    4545      !!---------------------------------------------------------------------- 
    46       INTEGER, INTENT(in)           ::   kt           ! Main time step counter 
    47       LOGICAL, INTENT(in), OPTIONAL ::   dyn3d_only   ! T => only update baroclinic velocities 
     46      INTEGER                             , INTENT(in)    ::   kt           ! Main time step counter 
     47      INTEGER                             , INTENT(in)    ::   Kbb, Kaa     ! Ocean time level indices 
     48      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) ::   puu, pvv     ! Ocean velocities (to be updated at open boundaries) 
     49      LOGICAL, OPTIONAL                   , INTENT(in)    ::   dyn3d_only   ! T => only update baroclinic velocities 
    4850      ! 
    4951      INTEGER ::   jk, ii, ij, ib_bdy, ib, igrd     ! Loop counter 
    5052      LOGICAL ::   ll_dyn2d, ll_dyn3d, ll_orlanski 
    51       REAL(wp), DIMENSION(jpi,jpj) :: pua2d, pva2d     ! after barotropic velocities 
     53      REAL(wp), DIMENSION(jpi,jpj) :: zua2d, zva2d     ! after barotropic velocities 
    5254      !!---------------------------------------------------------------------- 
    5355      ! 
     
    7072 
    7173      !                          ! "After" velocities:  
    72       pua2d(:,:) = 0._wp 
    73       pva2d(:,:) = 0._wp      
     74      zua2d(:,:) = 0._wp 
     75      zva2d(:,:) = 0._wp      
    7476      DO jk = 1, jpkm1 
    75          pua2d(:,:) = pua2d(:,:) + e3u_a(:,:,jk) * ua(:,:,jk) * umask(:,:,jk) 
    76          pva2d(:,:) = pva2d(:,:) + e3v_a(:,:,jk) * va(:,:,jk) * vmask(:,:,jk) 
     77         zua2d(:,:) = zua2d(:,:) + e3u(:,:,jk,Kaa) * puu(:,:,jk,Kaa) * umask(:,:,jk) 
     78         zva2d(:,:) = zva2d(:,:) + e3v(:,:,jk,Kaa) * pvv(:,:,jk,Kaa) * vmask(:,:,jk) 
    7779      END DO 
    78       pua2d(:,:) = pua2d(:,:) * r1_hu_a(:,:) 
    79       pva2d(:,:) = pva2d(:,:) * r1_hv_a(:,:) 
     80      zua2d(:,:) = zua2d(:,:) * r1_hu_a(:,:) 
     81      zva2d(:,:) = zva2d(:,:) * r1_hv_a(:,:) 
    8082 
    8183      DO jk = 1 , jpkm1 
    82          ua(:,:,jk) = ( ua(:,:,jk) - pua2d(:,:) ) * umask(:,:,jk) 
    83          va(:,:,jk) = ( va(:,:,jk) - pva2d(:,:) ) * vmask(:,:,jk) 
     84         puu(:,:,jk,Kaa) = ( puu(:,:,jk,Kaa) - zua2d(:,:) ) * umask(:,:,jk) 
     85         pvv(:,:,jk,Kaa) = ( pvv(:,:,jk,Kaa) - zva2d(:,:) ) * vmask(:,:,jk) 
    8486      END DO 
    8587 
     
    8789      IF( ll_orlanski ) THEN     ! "Before" velocities (Orlanski condition only)  
    8890         DO jk = 1 , jpkm1 
    89             ub(:,:,jk) = ( ub(:,:,jk) - ub_b(:,:) ) * umask(:,:,jk) 
    90             vb(:,:,jk) = ( vb(:,:,jk) - vb_b(:,:) ) * vmask(:,:,jk) 
     91            puu(:,:,jk,Kbb) = ( puu(:,:,jk,Kbb) - uu_b(:,:,Kbb) ) * umask(:,:,jk) 
     92            pvv(:,:,jk,Kbb) = ( pvv(:,:,jk,Kbb) - vv_b(:,:,Kbb) ) * vmask(:,:,jk) 
    9193         END DO 
    9294      ENDIF 
     
    9799      !------------------------------------------------------- 
    98100 
    99       IF( ll_dyn2d )   CALL bdy_dyn2d( kt, pua2d, pva2d, ub_b, vb_b, r1_hu_a(:,:), r1_hv_a(:,:), ssha ) 
     101      IF( ll_dyn2d )   CALL bdy_dyn2d( kt, zua2d, zva2d, uu_b(:,:,Kbb), vv_b(:,:,Kbb), r1_hu_a(:,:), r1_hv_a(:,:), ssh(:,:,Kaa) ) 
    100102 
    101       IF( ll_dyn3d )   CALL bdy_dyn3d( kt ) 
     103      IF( ll_dyn3d )   CALL bdy_dyn3d( kt, Kbb, puu, pvv, Kaa ) 
    102104 
    103105      !------------------------------------------------------- 
     
    106108      ! 
    107109      DO jk = 1 , jpkm1 
    108          ua(:,:,jk) = ( ua(:,:,jk) + pua2d(:,:) ) * umask(:,:,jk) 
    109          va(:,:,jk) = ( va(:,:,jk) + pva2d(:,:) ) * vmask(:,:,jk) 
     110         puu(:,:,jk,Kaa) = ( puu(:,:,jk,Kaa) + zua2d(:,:) ) * umask(:,:,jk) 
     111         pvv(:,:,jk,Kaa) = ( pvv(:,:,jk,Kaa) + zva2d(:,:) ) * vmask(:,:,jk) 
    110112      END DO 
    111113      ! 
    112114      IF ( ll_orlanski ) THEN 
    113115         DO jk = 1 , jpkm1 
    114             ub(:,:,jk) = ( ub(:,:,jk) + ub_b(:,:) ) * umask(:,:,jk) 
    115             vb(:,:,jk) = ( vb(:,:,jk) + vb_b(:,:) ) * vmask(:,:,jk) 
     116            puu(:,:,jk,Kbb) = ( puu(:,:,jk,Kbb) + uu_b(:,:,Kbb) ) * umask(:,:,jk) 
     117            pvv(:,:,jk,Kbb) = ( pvv(:,:,jk,Kbb) + vv_b(:,:,Kbb) ) * vmask(:,:,jk) 
    116118         END DO 
    117119      END IF 
  • 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 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/BDY/bdylib.F90

    r10529 r10957  
    3434CONTAINS 
    3535 
    36    SUBROUTINE bdy_frs( idx, pta, dta ) 
     36   SUBROUTINE bdy_frs( idx, phia, dta ) 
    3737      !!---------------------------------------------------------------------- 
    3838      !!                 ***  SUBROUTINE bdy_frs  *** 
     
    4444      TYPE(OBC_INDEX),                     INTENT(in) ::   idx  ! OBC indices 
    4545      REAL(wp), DIMENSION(:,:),            INTENT(in) ::   dta  ! OBC external data 
    46       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pta  ! tracer trend 
     46      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   phia  ! tracer trend 
    4747      !! 
    4848      REAL(wp) ::   zwgt           ! boundary weight 
     
    5757            ij = idx%nbj(ib,igrd) 
    5858            zwgt = idx%nbw(ib,igrd) 
    59             pta(ii,ij,ik) = ( pta(ii,ij,ik) + zwgt * (dta(ib,ik) - pta(ii,ij,ik) ) ) * tmask(ii,ij,ik) 
     59            phia(ii,ij,ik) = ( phia(ii,ij,ik) + zwgt * (dta(ib,ik) - phia(ii,ij,ik) ) ) * tmask(ii,ij,ik) 
    6060         END DO 
    6161      END DO 
     
    6464 
    6565 
    66    SUBROUTINE bdy_spe( idx, pta, dta ) 
     66   SUBROUTINE bdy_spe( idx, phia, dta ) 
    6767      !!---------------------------------------------------------------------- 
    6868      !!                 ***  SUBROUTINE bdy_spe  *** 
     
    7373      TYPE(OBC_INDEX),                     INTENT(in) ::   idx  ! OBC indices 
    7474      REAL(wp), DIMENSION(:,:),            INTENT(in) ::   dta  ! OBC external data 
    75       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pta  ! tracer trend 
     75      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   phia  ! tracer trend 
    7676      !! 
    7777      REAL(wp) ::   zwgt           ! boundary weight 
     
    8585         ij = idx%nbj(ib,igrd) 
    8686         DO ik = 1, jpkm1 
    87             pta(ii,ij,ik) = dta(ib,ik) * tmask(ii,ij,ik) 
     87            phia(ii,ij,ik) = dta(ib,ik) * tmask(ii,ij,ik) 
    8888         END DO 
    8989      END DO 
     
    9292 
    9393 
    94    SUBROUTINE bdy_orl( idx, ptb, pta, dta, ll_npo ) 
     94   SUBROUTINE bdy_orl( idx, phib, phia, dta, ll_npo ) 
    9595      !!---------------------------------------------------------------------- 
    9696      !!                 ***  SUBROUTINE bdy_orl  *** 
     
    102102      TYPE(OBC_INDEX),                     INTENT(in) ::   idx  ! OBC indices 
    103103      REAL(wp), DIMENSION(:,:),            INTENT(in) ::   dta  ! OBC external data 
    104       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   ptb  ! before tracer field 
    105       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pta  ! tracer trend 
     104      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   phib  ! before tracer field 
     105      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   phia  ! tracer trend 
    106106      LOGICAL,                             INTENT(in) ::   ll_npo  ! switch for NPO version 
    107107      !! 
     
    111111      igrd = 1                       ! Everything is at T-points here 
    112112      ! 
    113       CALL bdy_orlanski_3d( idx, igrd, ptb(:,:,:), pta(:,:,:), dta, ll_npo ) 
     113      CALL bdy_orlanski_3d( idx, igrd, phib(:,:,:), phia(:,:,:), dta, ll_npo ) 
    114114      ! 
    115115   END SUBROUTINE bdy_orl 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/BDY/bdytra.F90

    r10529 r10957  
    4040CONTAINS 
    4141 
    42    SUBROUTINE bdy_tra( kt ) 
     42   SUBROUTINE bdy_tra( kt, Kbb, pts, Kaa ) 
    4343      !!---------------------------------------------------------------------- 
    4444      !!                  ***  SUBROUTINE bdy_tra  *** 
     
    4747      !! 
    4848      !!---------------------------------------------------------------------- 
    49       INTEGER, INTENT(in) ::   kt   ! Main time step counter 
     49      INTEGER                                  , INTENT(in)    :: kt        ! Main time step counter 
     50      INTEGER                                  , INTENT(in)    :: Kbb, Kaa  ! time level indices 
     51      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) :: pts       ! tracer fields 
    5052      ! 
    51       INTEGER                        :: ib_bdy, jn, igrd   ! Loop indeces 
     53      INTEGER                        :: ib_bdy, jn, igrd   ! Loop indices 
    5254      TYPE(ztrabdy), DIMENSION(jpts) :: zdta               ! Temporary data structure 
    5355      !!---------------------------------------------------------------------- 
     
    6365            SELECT CASE( TRIM(cn_tra(ib_bdy)) ) 
    6466            CASE('none'        )   ;   CYCLE 
    65             CASE('frs'         )   ;   CALL bdy_frs ( idx_bdy(ib_bdy),                tsa(:,:,:,jn), zdta(jn)%tra ) 
    66             CASE('specified'   )   ;   CALL bdy_spe ( idx_bdy(ib_bdy),                tsa(:,:,:,jn), zdta(jn)%tra ) 
    67             CASE('neumann'     )   ;   CALL bdy_nmn ( idx_bdy(ib_bdy), igrd         , tsa(:,:,:,jn)               ) 
    68             CASE('orlanski'    )   ;   CALL bdy_orl ( idx_bdy(ib_bdy), tsb(:,:,:,jn), tsa(:,:,:,jn), zdta(jn)%tra, ll_npo=.false. ) 
    69             CASE('orlanski_npo')   ;   CALL bdy_orl ( idx_bdy(ib_bdy), tsb(:,:,:,jn), tsa(:,:,:,jn), zdta(jn)%tra, ll_npo=.true. ) 
    70             CASE('runoff'      )   ;   CALL bdy_rnf ( idx_bdy(ib_bdy),                tsa(:,:,:,jn),               jn ) 
     67            CASE('frs'         )   ;   CALL bdy_frs ( idx_bdy(ib_bdy),                pts(:,:,:,jn,Kaa), zdta(jn)%tra ) 
     68            CASE('specified'   )   ;   CALL bdy_spe ( idx_bdy(ib_bdy),                pts(:,:,:,jn,Kaa), zdta(jn)%tra ) 
     69            CASE('neumann'     )   ;   CALL bdy_nmn ( idx_bdy(ib_bdy), igrd         , pts(:,:,:,jn,Kaa)               ) 
     70            CASE('orlanski'    )   ;   CALL bdy_orl ( idx_bdy(ib_bdy), pts(:,:,:,jn,Kbb), pts(:,:,:,jn,Kaa), zdta(jn)%tra, ll_npo=.false. ) 
     71            CASE('orlanski_npo')   ;   CALL bdy_orl ( idx_bdy(ib_bdy), pts(:,:,:,jn,Kbb), pts(:,:,:,jn,Kaa), zdta(jn)%tra, ll_npo=.true. ) 
     72            CASE('runoff'      )   ;   CALL bdy_rnf ( idx_bdy(ib_bdy),                pts(:,:,:,jn,Kaa),               jn ) 
    7173            CASE DEFAULT           ;   CALL ctl_stop( 'bdy_tra : unrecognised option for open boundaries for T and S' ) 
    7274            END SELECT 
    7375            ! Boundary points should be updated 
    74             CALL lbc_bdy_lnk( 'bdytra', tsa(:,:,:,jn), 'T', 1., ib_bdy ) 
     76            CALL lbc_bdy_lnk( 'bdytra', pts(:,:,:,jn,Kaa), 'T', 1., ib_bdy ) 
    7577            !  
    7678         END DO 
     
    8082 
    8183 
    82    SUBROUTINE bdy_rnf( idx, pta, jpa ) 
     84   SUBROUTINE bdy_rnf( idx, pt, jpa ) 
    8385      !!---------------------------------------------------------------------- 
    8486      !!                 ***  SUBROUTINE bdy_rnf  *** 
     
    9092      !!---------------------------------------------------------------------- 
    9193      TYPE(OBC_INDEX),                     INTENT(in) ::   idx  ! OBC indices 
    92       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pta  ! tracer trend 
     94      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pt  ! tracer trend 
    9395      INTEGER,                             INTENT(in) ::   jpa  ! TRA index 
    9496      ! 
     
    105107            ip = bdytmask(ii+1,ij  ) - bdytmask(ii-1,ij  ) 
    106108            jp = bdytmask(ii  ,ij+1) - bdytmask(ii  ,ij-1) 
    107             if (jpa == jp_tem) pta(ii,ij,ik) = pta(ii+ip,ij+jp,ik) * tmask(ii,ij,ik) 
    108             if (jpa == jp_sal) pta(ii,ij,ik) =                 0.1 * tmask(ii,ij,ik) 
     109            if (jpa == jp_tem) pt(ii,ij,ik) = pt(ii+ip,ij+jp,ik) * tmask(ii,ij,ik) 
     110            if (jpa == jp_sal) pt(ii,ij,ik) =                 0.1 * tmask(ii,ij,ik) 
    109111         END DO 
    110112      END DO 
     
    113115 
    114116 
    115    SUBROUTINE bdy_tra_dmp( kt ) 
     117   SUBROUTINE bdy_tra_dmp( kt, Kbb, pts, Krhs ) 
    116118      !!---------------------------------------------------------------------- 
    117119      !!                 ***  SUBROUTINE bdy_tra_dmp  *** 
     
    120122      !!  
    121123      !!---------------------------------------------------------------------- 
    122       INTEGER, INTENT(in) ::   kt   ! 
     124      INTEGER                                  , INTENT(in)    :: kt        ! time step 
     125      INTEGER                                  , INTENT(in)    :: Kbb, Krhs ! time level indices 
     126      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) :: pts       ! active tracers and RHS of tracer equation 
    123127      ! 
    124128      REAL(wp) ::   zwgt           ! boundary weight 
     
    139143               zwgt = idx_bdy(ib_bdy)%nbd(ib,igrd) 
    140144               DO ik = 1, jpkm1 
    141                   zta = zwgt * ( dta_bdy(ib_bdy)%tem(ib,ik) - tsb(ii,ij,ik,jp_tem) ) * tmask(ii,ij,ik) 
    142                   zsa = zwgt * ( dta_bdy(ib_bdy)%sal(ib,ik) - tsb(ii,ij,ik,jp_sal) ) * tmask(ii,ij,ik) 
    143                   tsa(ii,ij,ik,jp_tem) = tsa(ii,ij,ik,jp_tem) + zta 
    144                   tsa(ii,ij,ik,jp_sal) = tsa(ii,ij,ik,jp_sal) + zsa 
     145                  zta = zwgt * ( dta_bdy(ib_bdy)%tem(ib,ik) - pts(ii,ij,ik,jp_tem,Kbb) ) * tmask(ii,ij,ik) 
     146                  zsa = zwgt * ( dta_bdy(ib_bdy)%sal(ib,ik) - pts(ii,ij,ik,jp_sal,Kbb) ) * tmask(ii,ij,ik) 
     147                  pts(ii,ij,ik,jp_tem,Krhs) = pts(ii,ij,ik,jp_tem,Krhs) + zta 
     148                  pts(ii,ij,ik,jp_sal,Krhs) = pts(ii,ij,ik,jp_sal,Krhs) + zsa 
    145149               END DO 
    146150            END DO 
Note: See TracChangeset for help on using the changeset viewer.