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 9939 for NEMO/branches/2018/dev_r9838_ENHANCE04_RK3/src/OCE/DOM/domvvl.F90 – NEMO

Ignore:
Timestamp:
2018-07-13T09:28:50+02:00 (6 years ago)
Author:
gm
Message:

#1911 (ENHANCE-04): RK3 branche phased with MLF@9937 branche

Location:
NEMO/branches/2018/dev_r9838_ENHANCE04_RK3
Files:
1 edited
1 copied

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2018/dev_r9838_ENHANCE04_RK3/src/OCE/DOM/domvvl.F90

    r9598 r9939  
    5454   LOGICAL , PUBLIC :: ln_vvl_dbg = .FALSE.                ! debug control prints 
    5555 
    56    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: un_td, vn_td                ! thickness diffusion transport 
    57    REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: hdiv_lf                     ! low frequency part of hz divergence 
    58    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: tilde_e3t_b, tilde_e3t_n    ! baroclinic scale factors 
    59    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: tilde_e3t_a, dtilde_e3t_a   ! baroclinic scale factors 
    60    REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   :: frq_rst_e3t                 ! retoring period for scale factors 
    61    REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   :: frq_rst_hdv                 ! retoring period for low freq. divergence 
     56   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: un_td, vn_td      ! thickness diffusion transport 
     57   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: hdiv_lf           ! low frequency part of hz divergence 
     58   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: te3t_b, te3t_n    ! baroclinic scale factors 
     59   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: te3t_a, dte3t_a   ! baroclinic scale factors 
     60   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   :: frq_rst_e3t       ! retoring period for scale factors 
     61   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   :: frq_rst_hdv       ! retoring period for low freq. divergence 
    6262 
    6363   !! * Substitutions 
     
    7676      IF( ln_vvl_zstar )   dom_vvl_alloc = 0 
    7777      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN 
    78          ALLOCATE( tilde_e3t_b(jpi,jpj,jpk)  , tilde_e3t_n(jpi,jpj,jpk) , tilde_e3t_a(jpi,jpj,jpk) ,   & 
    79             &      dtilde_e3t_a(jpi,jpj,jpk) , un_td  (jpi,jpj,jpk)     , vn_td  (jpi,jpj,jpk)     ,   & 
     78         ALLOCATE( te3t_b(jpi,jpj,jpk)  , te3t_n(jpi,jpj,jpk) , te3t_a(jpi,jpj,jpk) ,   & 
     79            &      dte3t_a(jpi,jpj,jpk) , un_td  (jpi,jpj,jpk)     , vn_td  (jpi,jpj,jpk)     ,   & 
    8080            &      STAT = dom_vvl_alloc        ) 
    8181         IF( lk_mpp             )   CALL mpp_sum ( dom_vvl_alloc ) 
     
    103103      !!               - interpolate scale factors 
    104104      !! 
    105       !! ** Action  : - e3t_(n/b) and tilde_e3t_(n/b) 
     105      !! ** Action  : - e3t_(n/b) and te3t_(n/b) 
    106106      !!              - Regrid: e3(u/v)_n 
    107107      !!                        e3(u/v)_b        
     
    117117      INTEGER ::   ji, jj, jk 
    118118      INTEGER ::   ii0, ii1, ij0, ij1 
    119       REAL(wp)::   zcoef 
     119      REAL(wp)::   zcoef, z1_Dt 
    120120      !!---------------------------------------------------------------------- 
    121121      ! 
     
    129129      IF( dom_vvl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dom_vvl_init : unable to allocate arrays' ) 
    130130      ! 
    131       !                    ! Read or initialize e3t_(b/n), tilde_e3t_(b/n) and hdiv_lf 
     131      !                    ! Read or initialize e3t_(b/n), te3t_(b/n) and hdiv_lf 
    132132      CALL dom_vvl_rst( nit000, 'READ' ) 
    133133      e3t_a(:,:,jpk) = e3t_0(:,:,jpk)  ! last level always inside the sea floor set one for all 
     
    208208         IF( ln_vvl_ztilde_as_zstar ) THEN   ! z-star emulation using z-tile 
    209209            frq_rst_e3t(:,:) = 0._wp               !Ignore namelist settings 
    210             frq_rst_hdv(:,:) = 1._wp / rdt 
     210            frq_rst_hdv(:,:) = 1._wp / rn_Dt 
    211211         ENDIF 
    212212         IF ( ln_vvl_zstar_at_eqtor ) THEN   ! use z-star in vicinity of the Equator 
     213            z1_Dt = 1._wp / rn_Dt 
    213214            DO jj = 1, jpj 
    214215               DO ji = 1, jpi 
     
    216217                  IF( ABS(gphit(ji,jj)) >= 6.) THEN 
    217218                     ! values outside the equatorial band and transition zone (ztilde) 
    218                      frq_rst_e3t(ji,jj) =  2.0_wp * rpi / ( MAX( rn_rst_e3t  , rsmall ) * 86400.e0_wp ) 
    219                      frq_rst_hdv(ji,jj) =  2.0_wp * rpi / ( MAX( rn_lf_cutoff, rsmall ) * 86400.e0_wp ) 
     219                     frq_rst_e3t(ji,jj) =  2._wp * rpi / ( MAX( rn_rst_e3t  , rsmall ) * 86400._wp ) 
     220                     frq_rst_hdv(ji,jj) =  2._wp * rpi / ( MAX( rn_lf_cutoff, rsmall ) * 86400._wp ) 
    220221                  ELSEIF( ABS(gphit(ji,jj)) <= 2.5) THEN    ! Equator strip ==> z-star 
    221222                     ! values inside the equatorial band (ztilde as zstar) 
    222                      frq_rst_e3t(ji,jj) =  0.0_wp 
    223                      frq_rst_hdv(ji,jj) =  1.0_wp / rdt 
     223                     frq_rst_e3t(ji,jj) =  0._wp 
     224                     frq_rst_hdv(ji,jj) =  z1_Dt 
    224225                  ELSE                                      ! transition band (2.5 to 6 degrees N/S) 
    225226                     !                                      ! (linearly transition from z-tilde to z-star) 
    226                      frq_rst_e3t(ji,jj) = 0.0_wp + (frq_rst_e3t(ji,jj)-0.0_wp)*0.5_wp   & 
    227                         &            * (  1.0_wp - COS( rad*(ABS(gphit(ji,jj))-2.5_wp)  & 
    228                         &                                          * 180._wp / 3.5_wp ) ) 
    229                      frq_rst_hdv(ji,jj) = (1.0_wp / rdt)                                & 
    230                         &            + (  frq_rst_hdv(ji,jj)-(1.e0_wp / rdt) )*0.5_wp   & 
    231                         &            * (  1._wp  - COS( rad*(ABS(gphit(ji,jj))-2.5_wp)  & 
    232                         &                                          * 180._wp / 3.5_wp ) ) 
     227                     frq_rst_e3t(ji,jj) = 0._wp + ( frq_rst_e3t(ji,jj) - 0._wp  ) * 0.5_wp                             & 
     228                        &                       * (  1._wp - COS( rad*(ABS(gphit(ji,jj))-2.5_wp) * 180._wp / 3.5_wp )  ) 
     229                     frq_rst_hdv(ji,jj) = z1_Dt + (  frq_rst_hdv(ji,jj) - z1_Dt ) * 0.5_wp                             & 
     230                        &                       * (  1._wp - COS( rad*(ABS(gphit(ji,jj))-2.5_wp) * 180._wp / 3.5_wp )  ) 
    233231                  ENDIF 
    234232               END DO 
     
    237235               ii0 = 103   ;   ii1 = 111        
    238236               ij0 = 128   ;   ij1 = 135   ;    
    239                frq_rst_e3t( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  0.0_wp 
    240                frq_rst_hdv( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  1.e0_wp / rdt 
     237               frq_rst_e3t( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  0._wp 
     238               frq_rst_hdv( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  z1_Dt 
    241239            ENDIF 
    242240         ENDIF 
     
    280278      !! 
    281279      !! ** Action  :  - hdiv_lf    : restoring towards full baroclinic divergence in z_tilde case 
    282       !!               - tilde_e3t_a: after increment of vertical scale factor  
     280      !!               - te3t_a: after increment of vertical scale factor  
    283281      !!                              in z_tilde case 
    284282      !!               - e3(t/u/v)_a 
     
    345343            IF( kt > nit000 ) THEN 
    346344               DO jk = 1, jpkm1 
    347                   hdiv_lf(:,:,jk) = hdiv_lf(:,:,jk) - rdt * frq_rst_hdv(:,:)   & 
     345                  hdiv_lf(:,:,jk) = hdiv_lf(:,:,jk) - rn_Dt * frq_rst_hdv(:,:)   & 
    348346                     &          * ( hdiv_lf(:,:,jk) - e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) ) 
    349347               END DO 
     
    353351         ! II - after z_tilde increments of vertical scale factors 
    354352         ! ======================================================= 
    355          tilde_e3t_a(:,:,:) = 0._wp  ! tilde_e3t_a used to store tendency terms 
     353         te3t_a(:,:,:) = 0._wp  ! te3t_a used to store tendency terms 
    356354 
    357355         ! 1 - High frequency divergence term 
     
    359357         IF( ln_vvl_ztilde ) THEN     ! z_tilde case 
    360358            DO jk = 1, jpkm1 
    361                tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - ( e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) - hdiv_lf(:,:,jk) ) 
     359               te3t_a(:,:,jk) = te3t_a(:,:,jk) - ( e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) - hdiv_lf(:,:,jk) ) 
    362360            END DO 
    363361         ELSE                         ! layer case 
    364362            DO jk = 1, jpkm1 
    365                tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) -   e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) * tmask(:,:,jk) 
     363               te3t_a(:,:,jk) = te3t_a(:,:,jk) -   e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) * tmask(:,:,jk) 
    366364            END DO 
    367365         ENDIF 
     
    371369         IF( ln_vvl_ztilde ) THEN 
    372370            DO jk = 1, jpk 
    373                tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - frq_rst_e3t(:,:) * tilde_e3t_b(:,:,jk) 
     371               te3t_a(:,:,jk) = te3t_a(:,:,jk) - frq_rst_e3t(:,:) * te3t_b(:,:,jk) 
    374372            END DO 
    375373         ENDIF 
     
    383381               DO ji = 1, fs_jpim1   ! vector opt. 
    384382                  un_td(ji,jj,jk) = rn_ahe3 * umask(ji,jj,jk) * e2_e1u(ji,jj)           & 
    385                      &            * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji+1,jj  ,jk) ) 
     383                     &            * ( te3t_b(ji,jj,jk) - te3t_b(ji+1,jj  ,jk) ) 
    386384                  vn_td(ji,jj,jk) = rn_ahe3 * vmask(ji,jj,jk) * e1_e2v(ji,jj)           &  
    387                      &            * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji  ,jj+1,jk) ) 
     385                     &            * ( te3t_b(ji,jj,jk) - te3t_b(ji  ,jj+1,jk) ) 
    388386                  zwu(ji,jj) = zwu(ji,jj) + un_td(ji,jj,jk) 
    389387                  zwv(ji,jj) = zwv(ji,jj) + vn_td(ji,jj,jk) 
     
    400398            DO jj = 2, jpjm1 
    401399               DO ji = fs_2, fs_jpim1   ! vector opt. 
    402                   tilde_e3t_a(ji,jj,jk) = tilde_e3t_a(ji,jj,jk) + (   un_td(ji-1,jj  ,jk) - un_td(ji,jj,jk)    & 
     400                  te3t_a(ji,jj,jk) = te3t_a(ji,jj,jk) + (   un_td(ji-1,jj  ,jk) - un_td(ji,jj,jk)    & 
    403401                     &                                          +     vn_td(ji  ,jj-1,jk) - vn_td(ji,jj,jk)    & 
    404402                     &                                            ) * r1_e1e2t(ji,jj) 
     
    414412         ! Leapfrog time stepping 
    415413         ! ~~~~~~~~~~~~~~~~~~~~~~ 
    416          IF( neuler == 0 .AND. kt == nit000 ) THEN 
    417             z2dt =  rdt 
    418          ELSE 
    419             z2dt = 2.0_wp * rdt 
    420          ENDIF 
    421          CALL lbc_lnk( tilde_e3t_a(:,:,:), 'T', 1._wp ) 
    422          tilde_e3t_a(:,:,:) = tilde_e3t_b(:,:,:) + z2dt * tmask(:,:,:) * tilde_e3t_a(:,:,:) 
     414         CALL lbc_lnk( te3t_a(:,:,:), 'T', 1._wp ) 
     415         te3t_a(:,:,:) = te3t_b(:,:,:) + z2dt * tmask(:,:,:) * te3t_a(:,:,:) 
    423416 
    424417         ! Maximum deformation control 
     
    426419         ze3t(:,:,jpk) = 0._wp 
    427420         DO jk = 1, jpkm1 
    428             ze3t(:,:,jk) = tilde_e3t_a(:,:,jk) / e3t_0(:,:,jk) * tmask(:,:,jk) * tmask_i(:,:) 
     421            ze3t(:,:,jk) = te3t_a(:,:,jk) / e3t_0(:,:,jk) * tmask(:,:,jk) * tmask_i(:,:) 
    429422         END DO 
    430423         z_tmax = MAXVAL( ze3t(:,:,:) ) 
     
    446439            ENDIF 
    447440            IF (lwp) THEN 
    448                WRITE(numout, *) 'MAX( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmax 
     441               WRITE(numout, *) 'MAX( te3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmax 
    449442               WRITE(numout, *) 'at i, j, k=', ijk_max 
    450                WRITE(numout, *) 'MIN( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmin 
     443               WRITE(numout, *) 'MIN( te3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmin 
    451444               WRITE(numout, *) 'at i, j, k=', ijk_min             
    452                CALL ctl_warn('MAX( ABS( tilde_e3t_a(:,:,:) ) / e3t_0(:,:,:) ) too high') 
     445               CALL ctl_warn('MAX( ABS( te3t_a(:,:,:) ) / e3t_0(:,:,:) ) too high') 
    453446            ENDIF 
    454447         ENDIF 
    455448         ! - ML - end test 
    456449         ! - ML - Imposing these limits will cause a baroclinicity error which is corrected for below 
    457          tilde_e3t_a(:,:,:) = MIN( tilde_e3t_a(:,:,:),   rn_zdef_max * e3t_0(:,:,:) ) 
    458          tilde_e3t_a(:,:,:) = MAX( tilde_e3t_a(:,:,:), - rn_zdef_max * e3t_0(:,:,:) ) 
     450         te3t_a(:,:,:) = MIN( te3t_a(:,:,:),   rn_zdef_max * e3t_0(:,:,:) ) 
     451         te3t_a(:,:,:) = MAX( te3t_a(:,:,:), - rn_zdef_max * e3t_0(:,:,:) ) 
    459452 
    460453         ! 
     
    462455         ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 
    463456         DO jk = 1, jpkm1 
    464             dtilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - tilde_e3t_b(:,:,jk) 
     457            dte3t_a(:,:,jk) = te3t_a(:,:,jk) - te3t_b(:,:,jk) 
    465458         END DO 
    466459         ! III - Barotropic repartition of the sea surface height over the baroclinic profile 
     
    470463         !        i.e. locally and not spread over the water column. 
    471464         !        (keep in mind that the idea is to reduce Eulerian velocity as much as possible) 
    472          zht(:,:) = 0. 
     465         zht(:,:) = 0._wp 
    473466         DO jk = 1, jpkm1 
    474             zht(:,:)  = zht(:,:) + tilde_e3t_a(:,:,jk) * tmask(:,:,jk) 
     467            zht(:,:)  = zht(:,:) + te3t_a(:,:,jk) * tmask(:,:,jk) 
    475468         END DO 
    476469         z_scale(:,:) =  - zht(:,:) / ( ht_0(:,:) + sshn(:,:) + 1. - ssmask(:,:) ) 
    477470         DO jk = 1, jpkm1 
    478             dtilde_e3t_a(:,:,jk) = dtilde_e3t_a(:,:,jk) + e3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk) 
     471            dte3t_a(:,:,jk) = dte3t_a(:,:,jk) + e3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk) 
    479472         END DO 
    480473 
     
    484477      !                                           ! ---baroclinic part--------- ! 
    485478         DO jk = 1, jpkm1 
    486             e3t_a(:,:,jk) = e3t_a(:,:,jk) + dtilde_e3t_a(:,:,jk) * tmask(:,:,jk) 
     479            e3t_a(:,:,jk) = e3t_a(:,:,jk) + dte3t_a(:,:,jk) * tmask(:,:,jk) 
    487480         END DO 
    488481      ENDIF 
     
    494487            z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( zht(:,:) ) ) 
    495488            IF( lk_mpp ) CALL mpp_max( z_tmax )                             ! max over the global domain 
    496             IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(SUM(tilde_e3t_a))) =', z_tmax 
     489            IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(SUM(te3t_a))) =', z_tmax 
    497490         END IF 
    498491         ! 
     
    573566      !!               - recompute depths and water height fields 
    574567      !! 
    575       !! ** Action  :  - e3t_(b/n), tilde_e3t_(b/n) and e3(u/v)_n ready for next time step 
     568      !! ** Action  :  - e3t_(b/n), te3t_(b/n) and e3(u/v)_n ready for next time step 
    576569      !!               - Recompute: 
    577570      !!                    e3(u/v)_b        
     
    587580      INTEGER, INTENT( in ) ::   kt   ! time step 
    588581      ! 
    589       INTEGER  ::   ji, jj, jk   ! dummy loop indices 
    590       REAL(wp) ::   zcoef        ! local scalar 
     582      INTEGER  ::   ji, jj, jk    ! dummy loop indices 
     583      REAL(wp) ::   zcoef, ze3f   ! local scalar 
    591584      !!---------------------------------------------------------------------- 
    592585      ! 
     
    605598      ! - ML - e3(t/u/v)_b are allready computed in dynnxt. 
    606599      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN 
    607          IF( neuler == 0 .AND. kt == nit000 ) THEN 
    608             tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:) 
     600         IF( l_1st_euler ) THEN 
     601            te3t_n(:,:,:) = te3t_a(:,:,:) 
    609602         ELSE 
    610             tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:) &  
    611             &         + atfp * ( tilde_e3t_b(:,:,:) - 2.0_wp * tilde_e3t_n(:,:,:) + tilde_e3t_a(:,:,:) ) 
     603            DO jk = 1, jpk 
     604               DO jj = 1, jpj 
     605                  DO ji = 1, jpi 
     606                     ze3f = te3t_n(ji,jj,jk)   & 
     607                        & + rn_atfp * ( te3t_b(ji,jj,jk) - 2.0_wp * te3t_n(ji,jj,jk) + te3t_a(ji,jj,jk) ) 
     608                     te3t_b(ji,jj,jk) = ze3f 
     609                     te3t_n(ji,jj,jk) = te3t_a(ji,jj,jk) 
     610                  END DO 
     611               END DO 
     612            END DO 
    612613         ENDIF 
    613          tilde_e3t_n(:,:,:) = tilde_e3t_a(:,:,:) 
    614614      ENDIF 
    615615      gdept_b(:,:,:) = gdept_n(:,:,:) 
     
    806806            CALL iom_get( numror, jpdom_autoglo, 'sshn'   , sshn, ldxios = lrxios    ) 
    807807            ! 
    808             id1 = iom_varid( numror, 'e3t_b', ldstop = .FALSE. ) 
    809             id2 = iom_varid( numror, 'e3t_n', ldstop = .FALSE. ) 
     808            id1 = iom_varid( numror, 'e3t_b'      , ldstop = .FALSE. ) 
     809            id2 = iom_varid( numror, 'e3t_n'      , ldstop = .FALSE. ) 
    810810            id3 = iom_varid( numror, 'tilde_e3t_b', ldstop = .FALSE. ) 
    811811            id4 = iom_varid( numror, 'tilde_e3t_n', ldstop = .FALSE. ) 
    812             id5 = iom_varid( numror, 'hdiv_lf', ldstop = .FALSE. ) 
     812            id5 = iom_varid( numror, 'hdiv_lf'    , ldstop = .FALSE. ) 
    813813            !                             ! --------- ! 
    814814            !                             ! all cases ! 
     
    823823                  e3t_b(:,:,:) = e3t_0(:,:,:) 
    824824               END WHERE 
    825                IF( neuler == 0 ) THEN 
     825               IF( l_1st_euler ) THEN 
    826826                  e3t_b(:,:,:) = e3t_n(:,:,:) 
    827827               ENDIF 
     
    829829               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_n not found in restart files' 
    830830               IF(lwp) write(numout,*) 'e3t_n set equal to e3t_b.' 
    831                IF(lwp) write(numout,*) 'neuler is forced to 0' 
     831               IF(lwp) write(numout,*) 'l_1st_euler is forced to true' 
    832832               CALL iom_get( numror, jpdom_autoglo, 'e3t_b', e3t_b(:,:,:), ldxios = lrxios ) 
    833833               e3t_n(:,:,:) = e3t_b(:,:,:) 
    834                neuler = 0 
     834               l_1st_euler = .TRUE. 
    835835            ELSE IF( id2 > 0 ) THEN 
    836836               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_b not found in restart files' 
    837837               IF(lwp) write(numout,*) 'e3t_b set equal to e3t_n.' 
    838                IF(lwp) write(numout,*) 'neuler is forced to 0' 
     838               IF(lwp) write(numout,*) 'l_1st_euler is forced to true' 
    839839               CALL iom_get( numror, jpdom_autoglo, 'e3t_n', e3t_n(:,:,:), ldxios = lrxios ) 
    840840               e3t_b(:,:,:) = e3t_n(:,:,:) 
    841                neuler = 0 
     841               l_1st_euler = .TRUE. 
    842842            ELSE 
    843843               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_n not found in restart file' 
    844844               IF(lwp) write(numout,*) 'Compute scale factor from sshn' 
    845                IF(lwp) write(numout,*) 'neuler is forced to 0' 
     845               IF(lwp) write(numout,*) 'l_1st_euler is forced to true' 
    846846               DO jk = 1, jpk 
    847847                  e3t_n(:,:,jk) =  e3t_0(:,:,jk) * ( ht_0(:,:) + sshn(:,:) ) & 
     
    850850               END DO 
    851851               e3t_b(:,:,:) = e3t_n(:,:,:) 
    852                neuler = 0 
     852               l_1st_euler = .TRUE. 
    853853            ENDIF 
    854854            !                             ! ----------- ! 
     
    862862               !                          ! ----------------------- ! 
    863863               IF( MIN( id3, id4 ) > 0 ) THEN  ! all required arrays exist 
    864                   CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_b', tilde_e3t_b(:,:,:), ldxios = lrxios ) 
    865                   CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_n', tilde_e3t_n(:,:,:), ldxios = lrxios ) 
     864                  CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_b', te3t_b(:,:,:), ldxios = lrxios ) 
     865                  CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_n', te3t_n(:,:,:), ldxios = lrxios ) 
    866866               ELSE                            ! one at least array is missing 
    867                   tilde_e3t_b(:,:,:) = 0.0_wp 
    868                   tilde_e3t_n(:,:,:) = 0.0_wp 
     867                  te3t_b(:,:,:) = 0.0_wp 
     868                  te3t_n(:,:,:) = 0.0_wp 
    869869               ENDIF 
    870870               !                          ! ------------ ! 
     
    942942 
    943943            IF( ln_vvl_ztilde .OR. ln_vvl_layer) THEN 
    944                tilde_e3t_b(:,:,:) = 0._wp 
    945                tilde_e3t_n(:,:,:) = 0._wp 
     944               te3t_b(:,:,:) = 0._wp 
     945               te3t_n(:,:,:) = 0._wp 
    946946               IF( ln_vvl_ztilde ) hdiv_lf(:,:,:) = 0._wp 
    947947            END IF 
     
    960960         IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN  ! z_tilde and layer cases ! 
    961961            !                                        ! ----------------------- ! 
    962             CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_b', tilde_e3t_b(:,:,:), ldxios = lwxios) 
    963             CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_n', tilde_e3t_n(:,:,:), ldxios = lwxios) 
     962            CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_b', te3t_b(:,:,:), ldxios = lwxios) 
     963            CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_n', te3t_n(:,:,:), ldxios = lwxios) 
    964964         END IF 
    965965         !                                           ! -------------!     
     
    10161016            WRITE(numout,*) '                         rn_rst_e3t     = 0.e0' 
    10171017            WRITE(numout,*) '            hard-wired : z-tilde cutoff frequency of low-pass filter (days)' 
    1018             WRITE(numout,*) '                         rn_lf_cutoff   = 1.0/rdt' 
     1018            WRITE(numout,*) '                         rn_lf_cutoff   = 1/rn_Dt' 
    10191019         ELSE 
    10201020            WRITE(numout,*) '      z-tilde to zstar restoration timescale (days)        rn_rst_e3t   = ', rn_rst_e3t 
Note: See TracChangeset for help on using the changeset viewer.