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 6020 for branches/UKMO/icebergs_restart_single_file/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90 – NEMO

Ignore:
Timestamp:
2015-12-08T12:39:53+01:00 (8 years ago)
Author:
timgraham
Message:

Merged with head of trunk

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/icebergs_restart_single_file/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90

    r6019 r6020  
    1111   !!             3.5  ! 2013-07  (J. Chanut) Switch to Forward-backward time stepping 
    1212   !!             3.6  ! 2013-11  (A. Coward) Update for z-tilde compatibility 
     13   !!             3.7  ! 2015-11  (J. Chanut) free surface simplification 
    1314   !!--------------------------------------------------------------------- 
    14 #if defined key_dynspg_ts   ||   defined key_esopa 
    1515   !!---------------------------------------------------------------------- 
    16    !!   'key_dynspg_ts'         split explicit free surface 
     16   !!                       split explicit free surface 
    1717   !!---------------------------------------------------------------------- 
    1818   !!   dyn_spg_ts  : compute surface pressure gradient trend using a time- 
     
    2323   USE sbc_oce         ! surface boundary condition: ocean 
    2424   USE sbcisf          ! ice shelf variable (fwfisf) 
    25    USE dynspg_oce      ! surface pressure gradient variables 
    2625   USE phycst          ! physical constants 
    2726   USE dynvor          ! vorticity term 
    2827   USE bdy_par         ! for lk_bdy 
    29    USE bdytides        ! open boundary condition data      
     28   USE bdytides        ! open boundary condition data 
    3029   USE bdydyn2d        ! open boundary conditions on barotropic variables 
    3130   USE sbctide         ! tides 
     
    6867   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::  ftsw, ftse   ! (only used with een vorticity scheme) 
    6968 
    70    ! Arrays below are saved to allow testing of the "no time averaging" option 
    71    ! If this option is not retained, these could be replaced by temporary arrays 
    72    REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::  sshbb_e, sshb_e, & ! Instantaneous barotropic arrays 
    73                                                    ubb_e, ub_e,     & 
    74                                                    vbb_e, vb_e 
    75  
    7669   !! * Substitutions 
    7770#  include "domzgr_substitute.h90" 
     
    8881      !!                  ***  routine dyn_spg_ts_alloc  *** 
    8982      !!---------------------------------------------------------------------- 
    90       INTEGER :: ierr(3) 
     83      INTEGER :: ierr(4) 
    9184      !!---------------------------------------------------------------------- 
    9285      ierr(:) = 0 
    9386 
    94       ALLOCATE( sshb_e(jpi,jpj), sshbb_e(jpi,jpj), & 
    95          &      ub_e(jpi,jpj)  , vb_e(jpi,jpj)   , & 
    96          &      ubb_e(jpi,jpj) , vbb_e(jpi,jpj)  , STAT= ierr(1) ) 
     87      ALLOCATE( ssha_e(jpi,jpj),  sshn_e(jpi,jpj), sshb_e(jpi,jpj), sshbb_e(jpi,jpj), & 
     88         &        ua_e(jpi,jpj),    un_e(jpi,jpj),   ub_e(jpi,jpj),   ubb_e(jpi,jpj), & 
     89         &        va_e(jpi,jpj),    vn_e(jpi,jpj),   vb_e(jpi,jpj),   vbb_e(jpi,jpj), & 
     90         &        hu_e(jpi,jpj),   hur_e(jpi,jpj),   hv_e(jpi,jpj),   hvr_e(jpi,jpj), STAT= ierr(1) ) 
    9791 
    9892      ALLOCATE( wgtbtp1(3*nn_baro), wgtbtp2(3*nn_baro), zwz(jpi,jpj), STAT= ierr(2) ) 
    9993 
    100       IF( ln_dynvor_een .or. ln_dynvor_een_old ) ALLOCATE( ftnw(jpi,jpj) , ftne(jpi,jpj) , &  
    101                                                     &      ftsw(jpi,jpj) , ftse(jpi,jpj) , STAT=ierr(3) ) 
     94      IF( ln_dynvor_een ) ALLOCATE( ftnw(jpi,jpj) , ftne(jpi,jpj) , &  
     95         &                          ftsw(jpi,jpj) , ftse(jpi,jpj) , STAT=ierr(3) ) 
     96 
     97      ALLOCATE( ub2_b(jpi,jpj), vb2_b(jpi,jpj), un_adv(jpi,jpj), vn_adv(jpi,jpj), & 
     98#if defined key_agrif 
     99         &      ub2_i_b(jpi,jpj), vb2_i_b(jpi,jpj)                              , & 
     100#endif 
     101         &      STAT= ierr(4)) 
    102102 
    103103      dyn_spg_ts_alloc = MAXVAL(ierr(:)) 
    104104 
    105105      IF( lk_mpp                )   CALL mpp_sum( dyn_spg_ts_alloc ) 
    106       IF( dyn_spg_ts_alloc /= 0 )   CALL ctl_warn('dynspg_oce_alloc: failed to allocate arrays') 
     106      IF( dyn_spg_ts_alloc /= 0 )   CALL ctl_warn('dyn_spg_ts_alloc: failed to allocate arrays') 
    107107      ! 
    108108   END FUNCTION dyn_spg_ts_alloc 
     109 
    109110 
    110111   SUBROUTINE dyn_spg_ts( kt ) 
     
    145146      INTEGER  ::   ikbu, ikbv, noffset      ! local integers 
    146147      REAL(wp) ::   zraur, z1_2dt_b, z2dt_bf    ! local scalars 
    147       REAL(wp) ::   zx1, zy1, zx2, zy2         !   -      - 
    148       REAL(wp) ::   z1_12, z1_8, z1_4, z1_2    !   -      - 
    149       REAL(wp) ::   zu_spg, zv_spg             !   -      - 
    150       REAL(wp) ::   zhura, zhvra               !   -      - 
    151       REAL(wp) ::   za0, za1, za2, za3           !   -      - 
    152       ! 
    153       REAL(wp), POINTER, DIMENSION(:,:) :: zun_e, zvn_e, zsshp2_e 
     148      REAL(wp) ::   zx1, zy1, zx2, zy2          !   -      - 
     149      REAL(wp) ::   z1_12, z1_8, z1_4, z1_2  !   -      - 
     150      REAL(wp) ::   zu_spg, zv_spg              !   -      - 
     151      REAL(wp) ::   zhura, zhvra          !   -      - 
     152      REAL(wp) ::   za0, za1, za2, za3    !   -      - 
     153      ! 
     154      REAL(wp), POINTER, DIMENSION(:,:) :: zsshp2_e 
    154155      REAL(wp), POINTER, DIMENSION(:,:) :: zu_trd, zv_trd, zu_frc, zv_frc, zssh_frc 
    155       REAL(wp), POINTER, DIMENSION(:,:) :: zu_sum, zv_sum, zwx, zwy, zhdiv 
     156      REAL(wp), POINTER, DIMENSION(:,:) :: zwx, zwy, zhdiv 
    156157      REAL(wp), POINTER, DIMENSION(:,:) :: zhup2_e, zhvp2_e, zhust_e, zhvst_e 
    157158      REAL(wp), POINTER, DIMENSION(:,:) :: zsshu_a, zsshv_a 
     
    163164      !                                         !* Allocate temporary arrays 
    164165      CALL wrk_alloc( jpi, jpj, zsshp2_e, zhdiv ) 
    165       CALL wrk_alloc( jpi, jpj, zu_trd, zv_trd, zun_e, zvn_e  ) 
    166       CALL wrk_alloc( jpi, jpj, zwx, zwy, zu_sum, zv_sum, zssh_frc, zu_frc, zv_frc) 
     166      CALL wrk_alloc( jpi, jpj, zu_trd, zv_trd) 
     167      CALL wrk_alloc( jpi, jpj, zwx, zwy, zssh_frc, zu_frc, zv_frc) 
    167168      CALL wrk_alloc( jpi, jpj, zhup2_e, zhvp2_e, zhust_e, zhvst_e) 
    168169      CALL wrk_alloc( jpi, jpj, zsshu_a, zsshv_a                                   ) 
     
    187188      ! 
    188189                                                       ! time offset in steps for bdy data update 
    189       IF (.NOT.ln_bt_fw) THEN ; noffset=-2*nn_baro ; ELSE ;  noffset = 0 ; ENDIF 
     190      IF (.NOT.ln_bt_fw) THEN ; noffset=-nn_baro ; ELSE ;  noffset = 0 ; ENDIF 
    190191      ! 
    191192      IF( kt == nit000 ) THEN                !* initialisation 
     
    219220      ! 
    220221      IF ( kt == nit000 .OR. lk_vvl ) THEN 
    221          IF ( ln_dynvor_een_old ) THEN 
    222             DO jj = 1, jpjm1 
    223                DO ji = 1, jpim1 
    224                   zwz(ji,jj) =   ( ht(ji  ,jj+1) + ht(ji+1,jj+1) +                    & 
    225                         &          ht(ji  ,jj  ) + ht(ji+1,jj  )   ) / 4._wp   
    226                   IF( zwz(ji,jj) /= 0._wp )   zwz(ji,jj) = 1._wp / zwz(ji,jj) 
    227                END DO 
    228             END DO 
     222         IF ( ln_dynvor_een ) THEN              !==  EEN scheme  ==! 
     223            SELECT CASE( nn_een_e3f )              !* ff/e3 at F-point 
     224            CASE ( 0 )                                   ! original formulation  (masked averaging of e3t divided by 4) 
     225               DO jj = 1, jpjm1 
     226                  DO ji = 1, jpim1 
     227                     zwz(ji,jj) =   ( ht(ji  ,jj+1) + ht(ji+1,jj+1) +                    & 
     228                        &             ht(ji  ,jj  ) + ht(ji+1,jj  )   ) / 4._wp   
     229                     IF( zwz(ji,jj) /= 0._wp )   zwz(ji,jj) = ff(ji,jj) / zwz(ji,jj) 
     230                  END DO 
     231               END DO 
     232            CASE ( 1 )                                   ! new formulation  (masked averaging of e3t divided by the sum of mask) 
     233               DO jj = 1, jpjm1 
     234                  DO ji = 1, jpim1 
     235                     zwz(ji,jj) =   ( ht(ji  ,jj+1) + ht(ji+1,jj+1) +                     & 
     236                        &             ht(ji  ,jj  ) + ht(ji+1,jj  )   )                   & 
     237                        &       / ( MAX( 1._wp, tmask(ji  ,jj+1, 1) + tmask(ji+1,jj+1, 1) +    & 
     238                        &                       tmask(ji  ,jj  , 1) + tmask(ji+1,jj  , 1) ) ) 
     239                     IF( zwz(ji,jj) /= 0._wp )   zwz(ji,jj) = ff(ji,jj) / zwz(ji,jj) 
     240                  END DO 
     241               END DO 
     242            END SELECT 
    229243            CALL lbc_lnk( zwz, 'F', 1._wp ) 
    230             zwz(:,:) = ff(:,:) * zwz(:,:) 
    231  
     244            ! 
    232245            ftne(1,:) = 0._wp ; ftnw(1,:) = 0._wp ; ftse(1,:) = 0._wp ; ftsw(1,:) = 0._wp 
    233246            DO jj = 2, jpj 
    234                DO ji = fs_2, jpi   ! vector opt. 
     247               DO ji = 2, jpi 
    235248                  ftne(ji,jj) = zwz(ji-1,jj  ) + zwz(ji  ,jj  ) + zwz(ji  ,jj-1) 
    236249                  ftnw(ji,jj) = zwz(ji-1,jj-1) + zwz(ji-1,jj  ) + zwz(ji  ,jj  ) 
     
    239252               END DO 
    240253            END DO 
    241          ELSE IF ( ln_dynvor_een ) THEN 
    242             DO jj = 1, jpjm1 
    243                DO ji = 1, jpim1 
    244                   zwz(ji,jj) =   ( ht(ji  ,jj+1) + ht(ji+1,jj+1) +                     & 
    245                         &          ht(ji  ,jj  ) + ht(ji+1,jj  )   )                   & 
    246                         &      / ( MAX( 1.0_wp, tmask(ji  ,jj+1, 1) + tmask(ji+1,jj+1, 1) +    & 
    247                         &                       tmask(ji  ,jj  , 1) + tmask(ji+1,jj  , 1) ) ) 
    248                   IF( zwz(ji,jj) /= 0._wp )   zwz(ji,jj) = 1._wp / zwz(ji,jj) 
    249                END DO 
    250             END DO 
    251             CALL lbc_lnk( zwz, 'F', 1._wp ) 
    252             zwz(:,:) = ff(:,:) * zwz(:,:) 
    253  
    254             ftne(1,:) = 0._wp ; ftnw(1,:) = 0._wp ; ftse(1,:) = 0._wp ; ftsw(1,:) = 0._wp 
    255             DO jj = 2, jpj 
    256                DO ji = fs_2, jpi   ! vector opt. 
    257                   ftne(ji,jj) = zwz(ji-1,jj  ) + zwz(ji  ,jj  ) + zwz(ji  ,jj-1) 
    258                   ftnw(ji,jj) = zwz(ji-1,jj-1) + zwz(ji-1,jj  ) + zwz(ji  ,jj  ) 
    259                   ftse(ji,jj) = zwz(ji  ,jj  ) + zwz(ji  ,jj-1) + zwz(ji-1,jj-1) 
    260                   ftsw(ji,jj) = zwz(ji  ,jj-1) + zwz(ji-1,jj-1) + zwz(ji-1,jj  ) 
    261                END DO 
    262             END DO 
    263          ELSE 
     254            ! 
     255         ELSE                                !== all other schemes (ENE, ENS, MIX) 
    264256            zwz(:,:) = 0._wp 
    265             zhf(:,:) = 0. 
     257            zhf(:,:) = 0._wp 
    266258            IF ( .not. ln_sco ) THEN 
     259 
     260!!gm  agree the JC comment  : this should be done in a much clear way 
     261 
    267262! JC: It not clear yet what should be the depth at f-points over land in z-coordinate case 
    268263!     Set it to zero for the time being  
     
    276271 
    277272            DO jj = 1, jpjm1 
    278                zhf(:,jj) = zhf(:,jj)*(1._wp- umask(:,jj,1) * umask(:,jj+1,1)) 
     273               zhf(:,jj) = zhf(:,jj) * (1._wp- umask(:,jj,1) * umask(:,jj+1,1)) 
    279274            END DO 
    280275 
     
    297292      ! If forward start at previous time step, and centered integration,  
    298293      ! then update averaging weights: 
    299       IF ((.NOT.ln_bt_fw).AND.((neuler==0).AND.(kt==nit000+1))) THEN 
     294      IF (.NOT.ln_bt_fw .AND.( neuler==0 .AND. kt==nit000+1 ) ) THEN 
    300295         ll_fw_start=.FALSE. 
    301296         CALL ts_wgt(ln_bt_av, ll_fw_start, icycle, wgtbtp1, wgtbtp2) 
     
    338333         DO jj = 2, jpjm1 
    339334            DO ji = fs_2, fs_jpim1   ! vector opt. 
    340                zy1 = ( zwy(ji,jj-1) + zwy(ji+1,jj-1) ) / e1u(ji,jj) 
    341                zy2 = ( zwy(ji,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj) 
    342                zx1 = ( zwx(ji-1,jj) + zwx(ji-1,jj+1) ) / e2v(ji,jj) 
    343                zx2 = ( zwx(ji  ,jj) + zwx(ji  ,jj+1) ) / e2v(ji,jj) 
     335               zy1 = ( zwy(ji,jj-1) + zwy(ji+1,jj-1) ) * r1_e1u(ji,jj) 
     336               zy2 = ( zwy(ji,jj  ) + zwy(ji+1,jj  ) ) * r1_e1u(ji,jj) 
     337               zx1 = ( zwx(ji-1,jj) + zwx(ji-1,jj+1) ) * r1_e2v(ji,jj) 
     338               zx2 = ( zwx(ji  ,jj) + zwx(ji  ,jj+1) ) * r1_e2v(ji,jj) 
    344339               ! energy conserving formulation for planetary vorticity term 
    345340               zu_trd(ji,jj) = z1_4 * ( zwz(ji  ,jj-1) * zy1 + zwz(ji,jj) * zy2 ) 
     
    352347            DO ji = fs_2, fs_jpim1   ! vector opt. 
    353348               zy1 =   z1_8 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) & 
    354                  &            + zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj) 
     349                 &            + zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) * r1_e1u(ji,jj) 
    355350               zx1 = - z1_8 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) & 
    356                  &            + zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) / e2v(ji,jj) 
     351                 &            + zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) * r1_e2v(ji,jj) 
    357352               zu_trd(ji,jj)  = zy1 * ( zwz(ji  ,jj-1) + zwz(ji,jj) ) 
    358353               zv_trd(ji,jj)  = zx1 * ( zwz(ji-1,jj  ) + zwz(ji,jj) ) 
     
    360355         END DO 
    361356         ! 
    362       ELSEIF ( ln_dynvor_een .or. ln_dynvor_een_old ) THEN  ! enstrophy and energy conserving scheme 
     357      ELSEIF ( ln_dynvor_een ) THEN  ! enstrophy and energy conserving scheme 
    363358         DO jj = 2, jpjm1 
    364359            DO ji = fs_2, fs_jpim1   ! vector opt. 
    365                zu_trd(ji,jj) = + z1_12 / e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) & 
    366                 &                                      + ftnw(ji+1,jj) * zwy(ji+1,jj  ) & 
    367                 &                                      + ftse(ji,jj  ) * zwy(ji  ,jj-1) & 
    368                 &                                      + ftsw(ji+1,jj) * zwy(ji+1,jj-1) ) 
    369                zv_trd(ji,jj) = - z1_12 / e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) & 
    370                 &                                      + ftse(ji,jj+1) * zwx(ji  ,jj+1) & 
    371                 &                                      + ftnw(ji,jj  ) * zwx(ji-1,jj  ) & 
    372                 &                                      + ftne(ji,jj  ) * zwx(ji  ,jj  ) ) 
     360               zu_trd(ji,jj) = + z1_12 * r1_e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) & 
     361                &                                         + ftnw(ji+1,jj) * zwy(ji+1,jj  ) & 
     362                &                                         + ftse(ji,jj  ) * zwy(ji  ,jj-1) & 
     363                &                                         + ftsw(ji+1,jj) * zwy(ji+1,jj-1) ) 
     364               zv_trd(ji,jj) = - z1_12 * r1_e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) & 
     365                &                                         + ftse(ji,jj+1) * zwx(ji  ,jj+1) & 
     366                &                                         + ftnw(ji,jj  ) * zwx(ji-1,jj  ) & 
     367                &                                         + ftne(ji,jj  ) * zwx(ji  ,jj  ) ) 
    373368            END DO 
    374369         END DO 
     
    381376         DO jj = 2, jpjm1  
    382377            DO ji = fs_2, fs_jpim1   ! vector opt. 
    383                zu_trd(ji,jj) = zu_trd(ji,jj) - grav * (  sshn(ji+1,jj  ) - sshn(ji  ,jj  )  ) / e1u(ji,jj) 
    384                zv_trd(ji,jj) = zv_trd(ji,jj) - grav * (  sshn(ji  ,jj+1) - sshn(ji  ,jj  )  ) / e2v(ji,jj) 
     378               zu_trd(ji,jj) = zu_trd(ji,jj) - grav * (  sshn(ji+1,jj  ) - sshn(ji  ,jj  )  ) * r1_e1u(ji,jj) 
     379               zv_trd(ji,jj) = zv_trd(ji,jj) - grav * (  sshn(ji  ,jj+1) - sshn(ji  ,jj  )  ) * r1_e2v(ji,jj) 
    385380            END DO 
    386381         END DO 
     
    431426            DO jj = 2, jpjm1               
    432427               DO ji = fs_2, fs_jpim1   ! vector opt. 
    433                   zu_spg =  grav * (  ssh_ib (ji+1,jj  ) - ssh_ib (ji,jj) ) /e1u(ji,jj) 
    434                   zv_spg =  grav * (  ssh_ib (ji  ,jj+1) - ssh_ib (ji,jj) ) /e2v(ji,jj) 
     428                  zu_spg =  grav * (  ssh_ib (ji+1,jj  ) - ssh_ib (ji,jj) ) * r1_e1u(ji,jj) 
     429                  zv_spg =  grav * (  ssh_ib (ji  ,jj+1) - ssh_ib (ji,jj) ) * r1_e2v(ji,jj) 
    435430                  zu_frc(ji,jj) = zu_frc(ji,jj) + zu_spg 
    436431                  zv_frc(ji,jj) = zv_frc(ji,jj) + zv_spg 
     
    441436               DO ji = fs_2, fs_jpim1   ! vector opt. 
    442437                  zu_spg =  grav * z1_2 * (  ssh_ib (ji+1,jj  ) - ssh_ib (ji,jj)    & 
    443                       &                    + ssh_ibb(ji+1,jj  ) - ssh_ibb(ji,jj)  ) /e1u(ji,jj) 
     438                      &                    + ssh_ibb(ji+1,jj  ) - ssh_ibb(ji,jj)  ) * r1_e1u(ji,jj) 
    444439                  zv_spg =  grav * z1_2 * (  ssh_ib (ji  ,jj+1) - ssh_ib (ji,jj)    & 
    445                       &                    + ssh_ibb(ji  ,jj+1) - ssh_ibb(ji,jj)  ) /e2v(ji,jj) 
     440                      &                    + ssh_ibb(ji  ,jj+1) - ssh_ibb(ji,jj)  ) * r1_e2v(ji,jj) 
    446441                  zu_frc(ji,jj) = zu_frc(ji,jj) + zu_spg 
    447442                  zv_frc(ji,jj) = zv_frc(ji,jj) + zv_spg 
     
    454449      !                                         ! Surface net water flux and rivers 
    455450      IF (ln_bt_fw) THEN 
    456          zssh_frc(:,:) = zraur * ( emp(:,:) - rnf(:,:) + rdivisf * fwfisf(:,:) ) 
     451         zssh_frc(:,:) = zraur * ( emp(:,:) - rnf(:,:) + fwfisf(:,:) ) 
    457452      ELSE 
    458453         zssh_frc(:,:) = zraur * z1_2 * (  emp(:,:) + emp_b(:,:) - rnf(:,:) - rnf_b(:,:)   & 
    459                 &                        + rdivisf * ( fwfisf(:,:) + fwfisf_b(:,:) )       ) 
     454                &                        + fwfisf(:,:) + fwfisf_b(:,:)                     ) 
    460455      ENDIF 
    461456#if defined key_asminc 
     
    465460      ENDIF 
    466461#endif 
    467       !                                   !* Fill boundary data arrays with AGRIF 
    468       !                                   ! ------------------------------------- 
     462      !                                   !* Fill boundary data arrays for AGRIF 
     463      !                                   ! ------------------------------------ 
    469464#if defined key_agrif 
    470465         IF( .NOT.Agrif_Root() ) CALL agrif_dta_ts( kt ) 
     
    490485      IF (ln_bt_fw) THEN                  ! FORWARD integration: start from NOW fields                     
    491486         sshn_e(:,:) = sshn (:,:)             
    492          zun_e (:,:) = un_b (:,:)             
    493          zvn_e (:,:) = vn_b (:,:) 
     487         un_e (:,:) = un_b (:,:)             
     488         vn_e (:,:) = vn_b (:,:) 
    494489         ! 
    495490         hu_e  (:,:) = hu   (:,:)        
     
    499494      ELSE                                ! CENTRED integration: start from BEFORE fields 
    500495         sshn_e(:,:) = sshb (:,:) 
    501          zun_e (:,:) = ub_b (:,:)          
    502          zvn_e (:,:) = vb_b (:,:) 
     496         un_e (:,:) = ub_b (:,:)          
     497         vn_e (:,:) = vb_b (:,:) 
    503498         ! 
    504499         hu_e  (:,:) = hu_b (:,:)        
     
    514509      va_b  (:,:) = 0._wp 
    515510      ssha  (:,:) = 0._wp       ! Sum for after averaged sea level 
    516       zu_sum(:,:) = 0._wp       ! Sum for now transport issued from ts loop 
    517       zv_sum(:,:) = 0._wp 
     511      un_adv(:,:) = 0._wp       ! Sum for now transport issued from ts loop 
     512      vn_adv(:,:) = 0._wp 
    518513      !                                             ! ==================== ! 
    519514      DO jn = 1, icycle                             !  sub-time-step loop  ! 
     
    523518         ! Update only tidal forcing at open boundaries 
    524519#if defined key_tide 
    525          IF ( lk_bdy .AND. lk_tide )      CALL bdy_dta_tides( kt, kit=jn, time_offset=(noffset+1) ) 
    526          IF ( ln_tide_pot .AND. lk_tide ) CALL upd_tide( kt, kit=jn, koffset=noffset ) 
     520         IF ( lk_bdy .AND. lk_tide ) CALL bdy_dta_tides( kt, kit=jn, time_offset=(noffset+1) ) 
     521         IF ( ln_tide_pot .AND. lk_tide ) CALL upd_tide( kt, kit=jn, time_offset=noffset ) 
    527522#endif 
    528523         ! 
     
    539534 
    540535         ! Extrapolate barotropic velocities at step jit+0.5: 
    541          ua_e(:,:) = za1 * zun_e(:,:) + za2 * ub_e(:,:) + za3 * ubb_e(:,:) 
    542          va_e(:,:) = za1 * zvn_e(:,:) + za2 * vb_e(:,:) + za3 * vbb_e(:,:) 
     536         ua_e(:,:) = za1 * un_e(:,:) + za2 * ub_e(:,:) + za3 * ubb_e(:,:) 
     537         va_e(:,:) = za1 * vn_e(:,:) + za2 * vb_e(:,:) + za3 * vbb_e(:,:) 
    543538 
    544539         IF( lk_vvl ) THEN                                !* Update ocean depth (variable volume case only) 
     
    549544            DO jj = 2, jpjm1                                    ! Sea Surface Height at u- & v-points 
    550545               DO ji = 2, fs_jpim1   ! Vector opt. 
    551                   zwx(ji,jj) = z1_2 * umask(ji,jj,1)  * r1_e12u(ji,jj)     & 
    552                      &              * ( e12t(ji  ,jj) * zsshp2_e(ji  ,jj)  & 
    553                      &              +   e12t(ji+1,jj) * zsshp2_e(ji+1,jj) ) 
    554                   zwy(ji,jj) = z1_2 * vmask(ji,jj,1)  * r1_e12v(ji,jj)     & 
    555                      &              * ( e12t(ji,jj  ) * zsshp2_e(ji,jj  )  & 
    556                      &              +   e12t(ji,jj+1) * zsshp2_e(ji,jj+1) ) 
     546                  zwx(ji,jj) = z1_2 * umask(ji,jj,1)  * r1_e1e2u(ji,jj)     & 
     547                     &              * ( e1e2t(ji  ,jj) * zsshp2_e(ji  ,jj)  & 
     548                     &              +   e1e2t(ji+1,jj) * zsshp2_e(ji+1,jj) ) 
     549                  zwy(ji,jj) = z1_2 * vmask(ji,jj,1)  * r1_e1e2v(ji,jj)     & 
     550                     &              * ( e1e2t(ji,jj  ) * zsshp2_e(ji,jj  )  & 
     551                     &              +   e1e2t(ji,jj+1) * zsshp2_e(ji,jj+1) ) 
    557552               END DO 
    558553            END DO 
     
    602597         ! Sum over sub-time-steps to compute advective velocities 
    603598         za2 = wgtbtp2(jn) 
    604          zu_sum  (:,:) = zu_sum  (:,:) + za2 * zwx  (:,:) / e2u  (:,:) 
    605          zv_sum  (:,:) = zv_sum  (:,:) + za2 * zwy  (:,:) / e1v  (:,:) 
     599         un_adv(:,:) = un_adv(:,:) + za2 * zwx(:,:) * r1_e2u(:,:) 
     600         vn_adv(:,:) = vn_adv(:,:) + za2 * zwy(:,:) * r1_e1v(:,:) 
    606601         ! 
    607602         ! Set next sea level: 
     
    609604            DO ji = fs_2, fs_jpim1   ! vector opt. 
    610605               zhdiv(ji,jj) = (   zwx(ji,jj) - zwx(ji-1,jj)   & 
    611                   &             + zwy(ji,jj) - zwy(ji,jj-1)   ) * r1_e12t(ji,jj) 
     606                  &             + zwy(ji,jj) - zwy(ji,jj-1)   ) * r1_e1e2t(ji,jj) 
    612607            END DO 
    613608         END DO 
     
    627622            DO jj = 2, jpjm1 
    628623               DO ji = 2, jpim1      ! NO Vector Opt. 
    629                   zsshu_a(ji,jj) = z1_2 * umask(ji,jj,1)  * r1_e12u(ji,jj)  & 
    630                      &              * ( e12t(ji  ,jj  ) * ssha_e(ji  ,jj  ) & 
    631                      &              +   e12t(ji+1,jj  ) * ssha_e(ji+1,jj  ) ) 
    632                   zsshv_a(ji,jj) = z1_2 * vmask(ji,jj,1)  * r1_e12v(ji,jj)  & 
    633                      &              * ( e12t(ji  ,jj  ) * ssha_e(ji  ,jj  ) & 
    634                      &              +   e12t(ji  ,jj+1) * ssha_e(ji  ,jj+1) ) 
     624                  zsshu_a(ji,jj) = z1_2 * umask(ji,jj,1)  * r1_e1e2u(ji,jj)  & 
     625                     &              * ( e1e2t(ji  ,jj  ) * ssha_e(ji  ,jj  ) & 
     626                     &              +   e1e2t(ji+1,jj  ) * ssha_e(ji+1,jj  ) ) 
     627                  zsshv_a(ji,jj) = z1_2 * vmask(ji,jj,1)  * r1_e1e2v(ji,jj)  & 
     628                     &              * ( e1e2t(ji  ,jj  ) * ssha_e(ji  ,jj  ) & 
     629                     &              +   e1e2t(ji  ,jj+1) * ssha_e(ji  ,jj+1) ) 
    635630               END DO 
    636631            END DO 
     
    666661            DO jj = 2, jpjm1                             
    667662               DO ji = 2, jpim1 
    668                   zx1 = z1_2 * umask(ji  ,jj,1) *  r1_e12u(ji  ,jj)    & 
    669                      &      * ( e12t(ji  ,jj  ) * zsshp2_e(ji  ,jj)    & 
    670                      &      +   e12t(ji+1,jj  ) * zsshp2_e(ji+1,jj  ) ) 
    671                   zy1 = z1_2 * vmask(ji  ,jj,1) *  r1_e12v(ji  ,jj  )  & 
    672                      &       * ( e12t(ji ,jj  ) * zsshp2_e(ji  ,jj  )  & 
    673                      &       +   e12t(ji ,jj+1) * zsshp2_e(ji  ,jj+1) ) 
     663                  zx1 = z1_2 * umask(ji  ,jj,1) *  r1_e1e2u(ji  ,jj)    & 
     664                     &      * ( e1e2t(ji  ,jj  ) * zsshp2_e(ji  ,jj)    & 
     665                     &      +   e1e2t(ji+1,jj  ) * zsshp2_e(ji+1,jj  ) ) 
     666                  zy1 = z1_2 * vmask(ji  ,jj,1) *  r1_e1e2v(ji  ,jj  )  & 
     667                     &       * ( e1e2t(ji ,jj  ) * zsshp2_e(ji  ,jj  )  & 
     668                     &       +   e1e2t(ji ,jj+1) * zsshp2_e(ji  ,jj+1) ) 
    674669                  zhust_e(ji,jj) = hu_0(ji,jj) + zx1  
    675670                  zhvst_e(ji,jj) = hv_0(ji,jj) + zy1 
     
    688683            DO jj = 2, jpjm1 
    689684               DO ji = fs_2, fs_jpim1   ! vector opt. 
    690                   zy1 = ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) ) / e1u(ji,jj) 
    691                   zy2 = ( zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj) 
    692                   zx1 = ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) ) / e2v(ji,jj) 
    693                   zx2 = ( zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) / e2v(ji,jj) 
     685                  zy1 = ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) ) * r1_e1u(ji,jj) 
     686                  zy2 = ( zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) * r1_e1u(ji,jj) 
     687                  zx1 = ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) ) * r1_e2v(ji,jj) 
     688                  zx2 = ( zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) * r1_e2v(ji,jj) 
    694689                  zu_trd(ji,jj) = z1_4 * ( zwz(ji  ,jj-1) * zy1 + zwz(ji,jj) * zy2 ) 
    695690                  zv_trd(ji,jj) =-z1_4 * ( zwz(ji-1,jj  ) * zx1 + zwz(ji,jj) * zx2 ) 
     
    701696               DO ji = fs_2, fs_jpim1   ! vector opt. 
    702697                  zy1 =   z1_8 * ( zwy(ji  ,jj-1) + zwy(ji+1,jj-1) & 
    703                    &             + zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) / e1u(ji,jj) 
     698                   &             + zwy(ji  ,jj  ) + zwy(ji+1,jj  ) ) * r1_e1u(ji,jj) 
    704699                  zx1 = - z1_8 * ( zwx(ji-1,jj  ) + zwx(ji-1,jj+1) & 
    705                    &             + zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) / e2v(ji,jj) 
     700                   &             + zwx(ji  ,jj  ) + zwx(ji  ,jj+1) ) * r1_e2v(ji,jj) 
    706701                  zu_trd(ji,jj)  = zy1 * ( zwz(ji  ,jj-1) + zwz(ji,jj) ) 
    707702                  zv_trd(ji,jj)  = zx1 * ( zwz(ji-1,jj  ) + zwz(ji,jj) ) 
     
    709704            END DO 
    710705            ! 
    711          ELSEIF ( ln_dynvor_een .or. ln_dynvor_een_old ) THEN !==  energy and enstrophy conserving scheme  ==! 
     706         ELSEIF ( ln_dynvor_een ) THEN !==  energy and enstrophy conserving scheme  ==! 
    712707            DO jj = 2, jpjm1 
    713708               DO ji = fs_2, fs_jpim1   ! vector opt. 
    714                   zu_trd(ji,jj) = + z1_12 / e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) & 
    715                      &                                    + ftnw(ji+1,jj) * zwy(ji+1,jj  ) & 
    716                      &                                    + ftse(ji,jj  ) * zwy(ji  ,jj-1) &  
    717                      &                                    + ftsw(ji+1,jj) * zwy(ji+1,jj-1) ) 
    718                   zv_trd(ji,jj) = - z1_12 / e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) &  
    719                      &                                    + ftse(ji,jj+1) * zwx(ji  ,jj+1) & 
    720                      &                                    + ftnw(ji,jj  ) * zwx(ji-1,jj  ) &  
    721                      &                                    + ftne(ji,jj  ) * zwx(ji  ,jj  ) ) 
     709                  zu_trd(ji,jj) = + z1_12 * r1_e1u(ji,jj) * (  ftne(ji,jj  ) * zwy(ji  ,jj  ) & 
     710                     &                                       + ftnw(ji+1,jj) * zwy(ji+1,jj  ) & 
     711                     &                                       + ftse(ji,jj  ) * zwy(ji  ,jj-1) &  
     712                     &                                       + ftsw(ji+1,jj) * zwy(ji+1,jj-1) ) 
     713                  zv_trd(ji,jj) = - z1_12 * r1_e2v(ji,jj) * (  ftsw(ji,jj+1) * zwx(ji-1,jj+1) &  
     714                     &                                       + ftse(ji,jj+1) * zwx(ji  ,jj+1) & 
     715                     &                                       + ftnw(ji,jj  ) * zwx(ji-1,jj  ) &  
     716                     &                                       + ftne(ji,jj  ) * zwx(ji  ,jj  ) ) 
    722717               END DO 
    723718            END DO 
     
    729724            DO jj = 2, jpjm1 
    730725               DO ji = fs_2, fs_jpim1   ! vector opt. 
    731                   zu_spg = grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) / e1u(ji,jj) 
    732                   zv_spg = grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) / e2v(ji,jj) 
     726                  zu_spg = grav * ( pot_astro(ji+1,jj) - pot_astro(ji,jj) ) * r1_e1u(ji,jj) 
     727                  zv_spg = grav * ( pot_astro(ji,jj+1) - pot_astro(ji,jj) ) * r1_e2v(ji,jj) 
    733728                  zu_trd(ji,jj) = zu_trd(ji,jj) + zu_spg 
    734729                  zv_trd(ji,jj) = zv_trd(ji,jj) + zv_spg 
     
    738733         ! 
    739734         ! Add bottom stresses: 
    740          zu_trd(:,:) = zu_trd(:,:) + bfrua(:,:) * zun_e(:,:) * hur_e(:,:) 
    741          zv_trd(:,:) = zv_trd(:,:) + bfrva(:,:) * zvn_e(:,:) * hvr_e(:,:) 
     735         zu_trd(:,:) = zu_trd(:,:) + bfrua(:,:) * un_e(:,:) * hur_e(:,:) 
     736         zv_trd(:,:) = zv_trd(:,:) + bfrva(:,:) * vn_e(:,:) * hvr_e(:,:) 
    742737         ! 
    743738         ! Surface pressure trend: 
     
    745740            DO ji = fs_2, fs_jpim1   ! vector opt. 
    746741               ! Add surface pressure gradient 
    747                zu_spg = - grav * ( zsshp2_e(ji+1,jj) - zsshp2_e(ji,jj) ) / e1u(ji,jj) 
    748                zv_spg = - grav * ( zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj) ) / e2v(ji,jj) 
     742               zu_spg = - grav * ( zsshp2_e(ji+1,jj) - zsshp2_e(ji,jj) ) * r1_e1u(ji,jj) 
     743               zv_spg = - grav * ( zsshp2_e(ji,jj+1) - zsshp2_e(ji,jj) ) * r1_e2v(ji,jj) 
    749744               zwx(ji,jj) = zu_spg 
    750745               zwy(ji,jj) = zv_spg 
     
    756751            DO jj = 2, jpjm1 
    757752               DO ji = fs_2, fs_jpim1   ! vector opt. 
    758                   ua_e(ji,jj) = (                                zun_e(ji,jj)   &  
     753                  ua_e(ji,jj) = (                                 un_e(ji,jj)   &  
    759754                            &     + rdtbt * (                      zwx(ji,jj)   & 
    760755                            &                                 + zu_trd(ji,jj)   & 
     
    762757                            &   ) * umask(ji,jj,1) 
    763758 
    764                   va_e(ji,jj) = (                                zvn_e(ji,jj)   & 
     759                  va_e(ji,jj) = (                                 vn_e(ji,jj)   & 
    765760                            &     + rdtbt * (                      zwy(ji,jj)   & 
    766761                            &                                 + zv_trd(ji,jj)   & 
     
    777772                  zhvra = vmask(ji,jj,1)/(hv_0(ji,jj) + zsshv_a(ji,jj) + 1._wp - vmask(ji,jj,1)) 
    778773 
    779                   ua_e(ji,jj) = (                hu_e(ji,jj)  *  zun_e(ji,jj)   &  
     774                  ua_e(ji,jj) = (                hu_e(ji,jj)  *   un_e(ji,jj)   &  
    780775                            &     + rdtbt * ( zhust_e(ji,jj)  *    zwx(ji,jj)   &  
    781776                            &               + zhup2_e(ji,jj)  * zu_trd(ji,jj)   & 
     
    783778                            &   ) * zhura 
    784779 
    785                   va_e(ji,jj) = (                hv_e(ji,jj)  *  zvn_e(ji,jj)   & 
     780                  va_e(ji,jj) = (                hv_e(ji,jj)  *   vn_e(ji,jj)   & 
    786781                            &     + rdtbt * ( zhvst_e(ji,jj)  *    zwy(ji,jj)   & 
    787782                            &               + zhvp2_e(ji,jj)  * zv_trd(ji,jj)   & 
     
    807802#if defined key_bdy   
    808803                                                           ! open boundaries 
    809          IF( lk_bdy ) CALL bdy_dyn2d( jn, ua_e, va_e, zun_e, zvn_e, hur_e, hvr_e, ssha_e ) 
     804         IF( lk_bdy ) CALL bdy_dyn2d( jn, ua_e, va_e, un_e, vn_e, hur_e, hvr_e, ssha_e ) 
    810805#endif 
    811806#if defined key_agrif                                                            
     
    815810         !                                             !  ---- 
    816811         ubb_e  (:,:) = ub_e  (:,:) 
    817          ub_e   (:,:) = zun_e (:,:) 
    818          zun_e  (:,:) = ua_e  (:,:) 
     812         ub_e   (:,:) = un_e (:,:) 
     813         un_e   (:,:) = ua_e  (:,:) 
    819814         ! 
    820815         vbb_e  (:,:) = vb_e  (:,:) 
    821          vb_e   (:,:) = zvn_e (:,:) 
    822          zvn_e  (:,:) = va_e  (:,:) 
     816         vb_e   (:,:) = vn_e (:,:) 
     817         vn_e   (:,:) = va_e  (:,:) 
    823818         ! 
    824819         sshbb_e(:,:) = sshb_e(:,:) 
     
    845840      ! ----------------------------------------------------------------------------- 
    846841      ! 
    847       ! At this stage ssha holds a time averaged value 
    848       !                                                ! Sea Surface Height at u-,v- and f-points 
    849       IF( lk_vvl ) THEN                                ! (required only in key_vvl case) 
    850          DO jj = 1, jpjm1 
    851             DO ji = 1, jpim1      ! NO Vector Opt. 
    852                zsshu_a(ji,jj) = z1_2 * umask(ji,jj,1)  * r1_e12u(ji,jj) & 
    853                   &              * ( e12t(ji  ,jj) * ssha(ji  ,jj)    & 
    854                   &              +   e12t(ji+1,jj) * ssha(ji+1,jj) ) 
    855                zsshv_a(ji,jj) = z1_2 * vmask(ji,jj,1)  * r1_e12v(ji,jj) & 
    856                   &              * ( e12t(ji,jj  ) * ssha(ji,jj  )    & 
    857                   &              +   e12t(ji,jj+1) * ssha(ji,jj+1) ) 
    858             END DO 
    859          END DO 
    860          CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp, zsshv_a, 'V', 1._wp ) ! Boundary conditions 
    861       ENDIF 
    862       ! 
    863842      ! Set advection velocity correction: 
     843      zwx(:,:) = un_adv(:,:) 
     844      zwy(:,:) = vn_adv(:,:) 
    864845      IF (((kt==nit000).AND.(neuler==0)).OR.(.NOT.ln_bt_fw)) THEN      
    865          un_adv(:,:) = zu_sum(:,:)*hur(:,:) 
    866          vn_adv(:,:) = zv_sum(:,:)*hvr(:,:) 
     846         un_adv(:,:) = zwx(:,:)*hur(:,:) 
     847         vn_adv(:,:) = zwy(:,:)*hvr(:,:) 
    867848      ELSE 
    868          un_adv(:,:) = z1_2 * ( ub2_b(:,:) + zu_sum(:,:)) * hur(:,:) 
    869          vn_adv(:,:) = z1_2 * ( vb2_b(:,:) + zv_sum(:,:)) * hvr(:,:) 
     849         un_adv(:,:) = z1_2 * ( ub2_b(:,:) + zwx(:,:)) * hur(:,:) 
     850         vn_adv(:,:) = z1_2 * ( vb2_b(:,:) + zwy(:,:)) * hvr(:,:) 
    870851      END IF 
    871852 
    872853      IF (ln_bt_fw) THEN ! Save integrated transport for next computation 
    873          ub2_b(:,:) = zu_sum(:,:) 
    874          vb2_b(:,:) = zv_sum(:,:) 
     854         ub2_b(:,:) = zwx(:,:) 
     855         vb2_b(:,:) = zwy(:,:) 
    875856      ENDIF 
    876857      ! 
     
    882863         END DO 
    883864      ELSE 
     865         ! At this stage, ssha has been corrected: compute new depths at velocity points 
     866         DO jj = 1, jpjm1 
     867            DO ji = 1, jpim1      ! NO Vector Opt. 
     868               zsshu_a(ji,jj) = z1_2 * umask(ji,jj,1)  * r1_e1e2u(ji,jj) & 
     869                  &              * ( e1e2t(ji  ,jj) * ssha(ji  ,jj)    & 
     870                  &              +   e1e2t(ji+1,jj) * ssha(ji+1,jj) ) 
     871               zsshv_a(ji,jj) = z1_2 * vmask(ji,jj,1)  * r1_e1e2v(ji,jj) & 
     872                  &              * ( e1e2t(ji,jj  ) * ssha(ji,jj  )    & 
     873                  &              +   e1e2t(ji,jj+1) * ssha(ji,jj+1) ) 
     874            END DO 
     875         END DO 
     876         CALL lbc_lnk_multi( zsshu_a, 'U', 1._wp, zsshv_a, 'V', 1._wp ) ! Boundary conditions 
     877         ! 
    884878         DO jk=1,jpkm1 
    885879            ua(:,:,jk) = ua(:,:,jk) + hur(:,:) * ( ua_b(:,:) - ub_b(:,:) * hu_b(:,:) ) * z1_2dt_b 
     
    900894#if defined key_agrif 
    901895      ! Save time integrated fluxes during child grid integration 
    902       ! (used to update coarse grid transports) 
    903       ! Useless with 2nd order momentum schemes 
     896      ! (used to update coarse grid transports at next time step) 
    904897      ! 
    905898      IF ( (.NOT.Agrif_Root()).AND.(ln_bt_fw) ) THEN 
     
    921914      ! 
    922915      CALL wrk_dealloc( jpi, jpj, zsshp2_e, zhdiv ) 
    923       CALL wrk_dealloc( jpi, jpj, zu_trd, zv_trd, zun_e, zvn_e ) 
    924       CALL wrk_dealloc( jpi, jpj, zwx, zwy, zu_sum, zv_sum, zssh_frc, zu_frc, zv_frc ) 
     916      CALL wrk_dealloc( jpi, jpj, zu_trd, zv_trd ) 
     917      CALL wrk_dealloc( jpi, jpj, zwx, zwy, zssh_frc, zu_frc, zv_frc ) 
    925918      CALL wrk_dealloc( jpi, jpj, zhup2_e, zhvp2_e, zhust_e, zhvst_e ) 
    926919      CALL wrk_dealloc( jpi, jpj, zsshu_a, zsshv_a                                   ) 
     
    10701063      ! 
    10711064      INTEGER  :: ji ,jj 
    1072       INTEGER  ::   ios                 ! Local integer output status for namelist read 
    10731065      REAL(wp) :: zxr2, zyr2, zcmax 
    10741066      REAL(wp), POINTER, DIMENSION(:,:) :: zcu 
    10751067      !! 
    1076       NAMELIST/namsplit/ ln_bt_fw, ln_bt_av, ln_bt_nn_auto, & 
    1077       &                  nn_baro, rn_bt_cmax, nn_bt_flt 
    10781068      !!---------------------------------------------------------------------- 
    10791069      ! 
    1080       REWIND( numnam_ref )              ! Namelist namsplit in reference namelist : time splitting parameters 
    1081       READ  ( numnam_ref, namsplit, IOSTAT = ios, ERR = 901) 
    1082 901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsplit in reference namelist', lwp ) 
    1083  
    1084       REWIND( numnam_cfg )              ! Namelist namsplit in configuration namelist : time splitting parameters 
    1085       READ  ( numnam_cfg, namsplit, IOSTAT = ios, ERR = 902 ) 
    1086 902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsplit in configuration namelist', lwp ) 
    1087       IF(lwm) WRITE ( numond, namsplit ) 
    1088       ! 
    1089       !         ! Max courant number for ext. grav. waves 
     1070      ! Max courant number for ext. grav. waves 
    10901071      ! 
    10911072      CALL wrk_alloc( jpi, jpj, zcu ) 
    10921073      ! 
    1093       IF (lk_vvl) THEN  
    1094          DO jj = 1, jpj 
    1095             DO ji =1, jpi 
    1096                zxr2 = 1./(e1t(ji,jj)*e1t(ji,jj)) 
    1097                zyr2 = 1./(e2t(ji,jj)*e2t(ji,jj)) 
    1098                zcu(ji,jj) = sqrt(grav*ht_0(ji,jj)*(zxr2 + zyr2) ) 
    1099             END DO 
    1100          END DO 
    1101       ELSE 
    1102          DO jj = 1, jpj 
    1103             DO ji =1, jpi 
    1104                zxr2 = 1./(e1t(ji,jj)*e1t(ji,jj)) 
    1105                zyr2 = 1./(e2t(ji,jj)*e2t(ji,jj)) 
    1106                zcu(ji,jj) = sqrt(grav*ht(ji,jj)*(zxr2 + zyr2) ) 
    1107             END DO 
    1108          END DO 
    1109       ENDIF 
    1110  
    1111       zcmax = MAXVAL(zcu(:,:)) 
     1074      DO jj = 1, jpj 
     1075         DO ji =1, jpi 
     1076            zxr2 = r1_e1t(ji,jj) * r1_e1t(ji,jj) 
     1077            zyr2 = r1_e2t(ji,jj) * r1_e2t(ji,jj) 
     1078            zcu(ji,jj) = SQRT( grav * ht_0(ji,jj) * (zxr2 + zyr2) ) 
     1079         END DO 
     1080      END DO 
     1081      ! 
     1082      zcmax = MAXVAL( zcu(:,:) ) 
    11121083      IF( lk_mpp )   CALL mpp_max( zcmax ) 
    11131084 
    11141085      ! Estimate number of iterations to satisfy a max courant number= rn_bt_cmax 
    1115       IF (ln_bt_nn_auto) nn_baro = CEILING( rdt / rn_bt_cmax * zcmax) 
     1086      IF (ln_bt_auto) nn_baro = CEILING( rdt / rn_bt_cmax * zcmax) 
    11161087       
    1117       rdtbt = rdt / FLOAT(nn_baro) 
     1088      rdtbt = rdt / REAL( nn_baro , wp ) 
    11181089      zcmax = zcmax * rdtbt 
    11191090                     ! Print results 
     
    11211092      IF(lwp) WRITE(numout,*) 'dyn_spg_ts : split-explicit free surface' 
    11221093      IF(lwp) WRITE(numout,*) '~~~~~~~~~~' 
    1123       IF( ln_bt_nn_auto ) THEN 
    1124          IF(lwp) WRITE(numout,*) '     ln_ts_nn_auto=.true. Automatically set nn_baro ' 
     1094      IF( ln_bt_auto ) THEN 
     1095         IF(lwp) WRITE(numout,*) '     ln_ts_auto=.true. Automatically set nn_baro ' 
    11251096         IF(lwp) WRITE(numout,*) '     Max. courant number allowed: ', rn_bt_cmax 
    11261097      ELSE 
    1127          IF(lwp) WRITE(numout,*) '     ln_ts_nn_auto=.false.: Use nn_baro in namelist ' 
     1098         IF(lwp) WRITE(numout,*) '     ln_ts_auto=.false.: Use nn_baro in namelist ' 
    11281099      ENDIF 
    11291100 
     
    11701141   END SUBROUTINE dyn_spg_ts_init 
    11711142 
    1172 #else 
    1173    !!--------------------------------------------------------------------------- 
    1174    !!   Default case :   Empty module   No split explicit free surface 
    1175    !!--------------------------------------------------------------------------- 
    1176 CONTAINS 
    1177    INTEGER FUNCTION dyn_spg_ts_alloc()    ! Dummy function 
    1178       dyn_spg_ts_alloc = 0 
    1179    END FUNCTION dyn_spg_ts_alloc 
    1180    SUBROUTINE dyn_spg_ts( kt )            ! Empty routine 
    1181       INTEGER, INTENT(in) :: kt 
    1182       WRITE(*,*) 'dyn_spg_ts: You should not have seen this print! error?', kt 
    1183    END SUBROUTINE dyn_spg_ts 
    1184    SUBROUTINE ts_rst( kt, cdrw )          ! Empty routine 
    1185       INTEGER         , INTENT(in) ::   kt         ! ocean time-step 
    1186       CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag 
    1187       WRITE(*,*) 'ts_rst    : You should not have seen this print! error?', kt, cdrw 
    1188    END SUBROUTINE ts_rst   
    1189    SUBROUTINE dyn_spg_ts_init( kt )       ! Empty routine 
    1190       INTEGER         , INTENT(in) ::   kt         ! ocean time-step 
    1191       WRITE(*,*) 'dyn_spg_ts_init   : You should not have seen this print! error?', kt 
    1192    END SUBROUTINE dyn_spg_ts_init 
    1193 #endif 
    1194     
    11951143   !!====================================================================== 
    11961144END MODULE dynspg_ts 
    1197  
    1198  
    1199  
Note: See TracChangeset for help on using the changeset viewer.