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 6827 for branches/2016/dev_r6409_SIMPLIF_2_usrdef_tools/NEMOGCM/NEMO/OPA_SRC/DOM/domvvl.F90 – NEMO

Ignore:
Timestamp:
2016-08-01T15:37:15+02:00 (8 years ago)
Author:
flavoni
Message:

#1692 - branch SIMPLIF_2_usrdef: commit routines Fortran to create domain_cfg.nc file, to have domain (input) files for new simplification branch. To be moved in TOOLS directory

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2016/dev_r6409_SIMPLIF_2_usrdef_tools/NEMOGCM/NEMO/OPA_SRC/DOM/domvvl.F90

    r6351 r6827  
    2222   USE phycst          ! physical constant 
    2323   USE dom_oce         ! ocean space and time domain 
    24    USE sbc_oce         ! ocean surface boundary condition 
    25    USE wet_dry         ! wetting and drying 
    26    USE restart         ! ocean restart 
    2724   ! 
    2825   USE in_out_manager  ! I/O manager 
     
    3734 
    3835   PUBLIC  dom_vvl_init       ! called by domain.F90 
    39    PUBLIC  dom_vvl_sf_nxt     ! called by step.F90 
    40    PUBLIC  dom_vvl_sf_swp     ! called by step.F90 
    41    PUBLIC  dom_vvl_interpol   ! called by dynnxt.F90 
    4236 
    4337   !                                                      !!* Namelist nam_vvl 
     
    133127      ! 
    134128      !                    ! Read or initialize e3t_(b/n), tilde_e3t_(b/n) and hdiv_lf 
    135       CALL dom_vvl_rst( nit000, 'READ' ) 
    136129      e3t_a(:,:,jpk) = e3t_0(:,:,jpk)  ! last level always inside the sea floor set one for all 
    137130      ! 
     
    246239 
    247240 
    248    SUBROUTINE dom_vvl_sf_nxt( kt, kcall )  
    249       !!---------------------------------------------------------------------- 
    250       !!                ***  ROUTINE dom_vvl_sf_nxt  *** 
    251       !!                    
    252       !! ** Purpose :  - compute the after scale factors used in tra_zdf, dynnxt, 
    253       !!                 tranxt and dynspg routines 
    254       !! 
    255       !! ** Method  :  - z_star case:  Repartition of ssh INCREMENT proportionnaly to the level thickness. 
    256       !!               - z_tilde_case: after scale factor increment =  
    257       !!                                    high frequency part of horizontal divergence 
    258       !!                                  + retsoring towards the background grid 
    259       !!                                  + thickness difusion 
    260       !!                               Then repartition of ssh INCREMENT proportionnaly 
    261       !!                               to the "baroclinic" level thickness. 
    262       !! 
    263       !! ** Action  :  - hdiv_lf    : restoring towards full baroclinic divergence in z_tilde case 
    264       !!               - tilde_e3t_a: after increment of vertical scale factor  
    265       !!                              in z_tilde case 
    266       !!               - e3(t/u/v)_a 
    267       !! 
    268       !! Reference  : Leclair, M., and Madec, G. 2011, Ocean Modelling. 
    269       !!---------------------------------------------------------------------- 
    270       INTEGER, INTENT( in )           ::   kt      ! time step 
    271       INTEGER, INTENT( in ), OPTIONAL ::   kcall   ! optional argument indicating call sequence 
    272       ! 
    273       INTEGER                ::   ji, jj, jk            ! dummy loop indices 
    274       INTEGER , DIMENSION(3) ::   ijk_max, ijk_min      ! temporary integers 
    275       REAL(wp)               ::   z2dt, z_tmin, z_tmax  ! local scalars 
    276       LOGICAL                ::   ll_do_bclinic         ! local logical 
    277       REAL(wp), POINTER, DIMENSION(:,:,:) ::   ze3t 
    278       REAL(wp), POINTER, DIMENSION(:,:  ) ::   zht, z_scale, zwu, zwv, zhdiv 
    279       !!---------------------------------------------------------------------- 
    280       ! 
    281       IF( ln_linssh )   RETURN      ! No calculation in linear free surface 
    282       ! 
    283       IF( nn_timing == 1 )   CALL timing_start('dom_vvl_sf_nxt') 
    284       ! 
    285       CALL wrk_alloc( jpi,jpj,zht,   z_scale, zwu, zwv, zhdiv ) 
    286       CALL wrk_alloc( jpi,jpj,jpk,   ze3t ) 
    287  
    288       IF( kt == nit000 ) THEN 
    289          IF(lwp) WRITE(numout,*) 
    290          IF(lwp) WRITE(numout,*) 'dom_vvl_sf_nxt : compute after scale factors' 
    291          IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~' 
    292       ENDIF 
    293  
    294       ll_do_bclinic = .TRUE. 
    295       IF( PRESENT(kcall) ) THEN 
    296          IF( kcall == 2 .AND. ln_vvl_ztilde )   ll_do_bclinic = .FALSE. 
    297       ENDIF 
    298  
    299       ! ******************************* ! 
    300       ! After acale factors at t-points ! 
    301       ! ******************************* ! 
    302       !                                                ! --------------------------------------------- ! 
    303       !                                                ! z_star coordinate and barotropic z-tilde part ! 
    304       !                                                ! --------------------------------------------- ! 
    305       ! 
    306       z_scale(:,:) = ( ssha(:,:) - sshb(:,:) ) * ssmask(:,:) / ( ht_0(:,:) + sshn(:,:) + 1. - ssmask(:,:) ) 
    307       DO jk = 1, jpkm1 
    308          ! formally this is the same as e3t_a = e3t_0*(1+ssha/ht_0) 
    309          e3t_a(:,:,jk) = e3t_b(:,:,jk) + e3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk) 
    310       END DO 
    311       ! 
    312       IF( ln_vvl_ztilde .OR. ln_vvl_layer .AND. ll_do_bclinic ) THEN   ! z_tilde or layer coordinate ! 
    313          !                                                            ! ------baroclinic part------ ! 
    314          ! I - initialization 
    315          ! ================== 
    316  
    317          ! 1 - barotropic divergence 
    318          ! ------------------------- 
    319          zhdiv(:,:) = 0._wp 
    320          zht(:,:)   = 0._wp 
    321          DO jk = 1, jpkm1 
    322             zhdiv(:,:) = zhdiv(:,:) + e3t_n(:,:,jk) * hdivn(:,:,jk) 
    323             zht  (:,:) = zht  (:,:) + e3t_n(:,:,jk) * tmask(:,:,jk) 
    324          END DO 
    325          zhdiv(:,:) = zhdiv(:,:) / ( zht(:,:) + 1. - tmask_i(:,:) ) 
    326  
    327          ! 2 - Low frequency baroclinic horizontal divergence  (z-tilde case only) 
    328          ! -------------------------------------------------- 
    329          IF( ln_vvl_ztilde ) THEN 
    330             IF( kt > nit000 ) THEN 
    331                DO jk = 1, jpkm1 
    332                   hdiv_lf(:,:,jk) = hdiv_lf(:,:,jk) - rdt * frq_rst_hdv(:,:)   & 
    333                      &          * ( hdiv_lf(:,:,jk) - e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) ) 
    334                END DO 
    335             ENDIF 
    336          ENDIF 
    337  
    338          ! II - after z_tilde increments of vertical scale factors 
    339          ! ======================================================= 
    340          tilde_e3t_a(:,:,:) = 0._wp  ! tilde_e3t_a used to store tendency terms 
    341  
    342          ! 1 - High frequency divergence term 
    343          ! ---------------------------------- 
    344          IF( ln_vvl_ztilde ) THEN     ! z_tilde case 
    345             DO jk = 1, jpkm1 
    346                tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - ( e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) - hdiv_lf(:,:,jk) ) 
    347             END DO 
    348          ELSE                         ! layer case 
    349             DO jk = 1, jpkm1 
    350                tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) -   e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) * tmask(:,:,jk) 
    351             END DO 
    352          ENDIF 
    353  
    354          ! 2 - Restoring term (z-tilde case only) 
    355          ! ------------------ 
    356          IF( ln_vvl_ztilde ) THEN 
    357             DO jk = 1, jpk 
    358                tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - frq_rst_e3t(:,:) * tilde_e3t_b(:,:,jk) 
    359             END DO 
    360          ENDIF 
    361  
    362          ! 3 - Thickness diffusion term 
    363          ! ---------------------------- 
    364          zwu(:,:) = 0._wp 
    365          zwv(:,:) = 0._wp 
    366          DO jk = 1, jpkm1        ! a - first derivative: diffusive fluxes 
    367             DO jj = 1, jpjm1 
    368                DO ji = 1, fs_jpim1   ! vector opt. 
    369                   un_td(ji,jj,jk) = rn_ahe3 * umask(ji,jj,jk) * e2_e1u(ji,jj)           & 
    370                      &            * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji+1,jj  ,jk) ) 
    371                   vn_td(ji,jj,jk) = rn_ahe3 * vmask(ji,jj,jk) * e1_e2v(ji,jj)           &  
    372                      &            * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji  ,jj+1,jk) ) 
    373                   zwu(ji,jj) = zwu(ji,jj) + un_td(ji,jj,jk) 
    374                   zwv(ji,jj) = zwv(ji,jj) + vn_td(ji,jj,jk) 
    375                END DO 
    376             END DO 
    377          END DO 
    378          DO jj = 1, jpj          ! b - correction for last oceanic u-v points 
    379             DO ji = 1, jpi 
    380                un_td(ji,jj,mbku(ji,jj)) = un_td(ji,jj,mbku(ji,jj)) - zwu(ji,jj) 
    381                vn_td(ji,jj,mbkv(ji,jj)) = vn_td(ji,jj,mbkv(ji,jj)) - zwv(ji,jj) 
    382             END DO 
    383          END DO 
    384          DO jk = 1, jpkm1        ! c - second derivative: divergence of diffusive fluxes 
    385             DO jj = 2, jpjm1 
    386                DO ji = fs_2, fs_jpim1   ! vector opt. 
    387                   tilde_e3t_a(ji,jj,jk) = tilde_e3t_a(ji,jj,jk) + (   un_td(ji-1,jj  ,jk) - un_td(ji,jj,jk)    & 
    388                      &                                          +     vn_td(ji  ,jj-1,jk) - vn_td(ji,jj,jk)    & 
    389                      &                                            ) * r1_e1e2t(ji,jj) 
    390                END DO 
    391             END DO 
    392          END DO 
    393          !                       ! d - thickness diffusion transport: boundary conditions 
    394          !                             (stored for tracer advction and continuity equation) 
    395          CALL lbc_lnk( un_td , 'U' , -1._wp) 
    396          CALL lbc_lnk( vn_td , 'V' , -1._wp) 
    397  
    398          ! 4 - Time stepping of baroclinic scale factors 
    399          ! --------------------------------------------- 
    400          ! Leapfrog time stepping 
    401          ! ~~~~~~~~~~~~~~~~~~~~~~ 
    402          IF( neuler == 0 .AND. kt == nit000 ) THEN 
    403             z2dt =  rdt 
    404          ELSE 
    405             z2dt = 2.0_wp * rdt 
    406          ENDIF 
    407          CALL lbc_lnk( tilde_e3t_a(:,:,:), 'T', 1._wp ) 
    408          tilde_e3t_a(:,:,:) = tilde_e3t_b(:,:,:) + z2dt * tmask(:,:,:) * tilde_e3t_a(:,:,:) 
    409  
    410          ! Maximum deformation control 
    411          ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~ 
    412          ze3t(:,:,jpk) = 0._wp 
    413          DO jk = 1, jpkm1 
    414             ze3t(:,:,jk) = tilde_e3t_a(:,:,jk) / e3t_0(:,:,jk) * tmask(:,:,jk) * tmask_i(:,:) 
    415          END DO 
    416          z_tmax = MAXVAL( ze3t(:,:,:) ) 
    417          IF( lk_mpp )   CALL mpp_max( z_tmax )                 ! max over the global domain 
    418          z_tmin = MINVAL( ze3t(:,:,:) ) 
    419          IF( lk_mpp )   CALL mpp_min( z_tmin )                 ! min over the global domain 
    420          ! - ML - test: for the moment, stop simulation for too large e3_t variations 
    421          IF( ( z_tmax >  rn_zdef_max ) .OR. ( z_tmin < - rn_zdef_max ) ) THEN 
    422             IF( lk_mpp ) THEN 
    423                CALL mpp_maxloc( ze3t, tmask, z_tmax, ijk_max(1), ijk_max(2), ijk_max(3) ) 
    424                CALL mpp_minloc( ze3t, tmask, z_tmin, ijk_min(1), ijk_min(2), ijk_min(3) ) 
    425             ELSE 
    426                ijk_max = MAXLOC( ze3t(:,:,:) ) 
    427                ijk_max(1) = ijk_max(1) + nimpp - 1 
    428                ijk_max(2) = ijk_max(2) + njmpp - 1 
    429                ijk_min = MINLOC( ze3t(:,:,:) ) 
    430                ijk_min(1) = ijk_min(1) + nimpp - 1 
    431                ijk_min(2) = ijk_min(2) + njmpp - 1 
    432             ENDIF 
    433             IF (lwp) THEN 
    434                WRITE(numout, *) 'MAX( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmax 
    435                WRITE(numout, *) 'at i, j, k=', ijk_max 
    436                WRITE(numout, *) 'MIN( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmin 
    437                WRITE(numout, *) 'at i, j, k=', ijk_min             
    438                CALL ctl_warn('MAX( ABS( tilde_e3t_a(:,:,:) ) / e3t_0(:,:,:) ) too high') 
    439             ENDIF 
    440          ENDIF 
    441          ! - ML - end test 
    442          ! - ML - Imposing these limits will cause a baroclinicity error which is corrected for below 
    443          tilde_e3t_a(:,:,:) = MIN( tilde_e3t_a(:,:,:),   rn_zdef_max * e3t_0(:,:,:) ) 
    444          tilde_e3t_a(:,:,:) = MAX( tilde_e3t_a(:,:,:), - rn_zdef_max * e3t_0(:,:,:) ) 
    445  
    446          ! 
    447          ! "tilda" change in the after scale factor 
    448          ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 
    449          DO jk = 1, jpkm1 
    450             dtilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - tilde_e3t_b(:,:,jk) 
    451          END DO 
    452          ! III - Barotropic repartition of the sea surface height over the baroclinic profile 
    453          ! ================================================================================== 
    454          ! add ( ssh increment + "baroclinicity error" ) proportionly to e3t(n) 
    455          ! - ML - baroclinicity error should be better treated in the future 
    456          !        i.e. locally and not spread over the water column. 
    457          !        (keep in mind that the idea is to reduce Eulerian velocity as much as possible) 
    458          zht(:,:) = 0. 
    459          DO jk = 1, jpkm1 
    460             zht(:,:)  = zht(:,:) + tilde_e3t_a(:,:,jk) * tmask(:,:,jk) 
    461          END DO 
    462          z_scale(:,:) =  - zht(:,:) / ( ht_0(:,:) + sshn(:,:) + 1. - ssmask(:,:) ) 
    463          DO jk = 1, jpkm1 
    464             dtilde_e3t_a(:,:,jk) = dtilde_e3t_a(:,:,jk) + e3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk) 
    465          END DO 
    466  
    467       ENDIF 
    468  
    469       IF( ln_vvl_ztilde .OR. ln_vvl_layer )  THEN   ! z_tilde or layer coordinate ! 
    470       !                                           ! ---baroclinic part--------- ! 
    471          DO jk = 1, jpkm1 
    472             e3t_a(:,:,jk) = e3t_a(:,:,jk) + dtilde_e3t_a(:,:,jk) * tmask(:,:,jk) 
    473          END DO 
    474       ENDIF 
    475  
    476       IF( ln_vvl_dbg .AND. .NOT. ll_do_bclinic ) THEN   ! - ML - test: control prints for debuging 
    477          ! 
    478          IF( lwp ) WRITE(numout, *) 'kt =', kt 
    479          IF ( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN 
    480             z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( zht(:,:) ) ) 
    481             IF( lk_mpp ) CALL mpp_max( z_tmax )                             ! max over the global domain 
    482             IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(SUM(tilde_e3t_a))) =', z_tmax 
    483          END IF 
    484          ! 
    485          zht(:,:) = 0.0_wp 
    486          DO jk = 1, jpkm1 
    487             zht(:,:) = zht(:,:) + e3t_n(:,:,jk) * tmask(:,:,jk) 
    488          END DO 
    489          z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + sshn(:,:) - zht(:,:) ) ) 
    490          IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain 
    491          IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+sshn-SUM(e3t_n))) =', z_tmax 
    492          ! 
    493          zht(:,:) = 0.0_wp 
    494          DO jk = 1, jpkm1 
    495             zht(:,:) = zht(:,:) + e3t_a(:,:,jk) * tmask(:,:,jk) 
    496          END DO 
    497          z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + ssha(:,:) - zht(:,:) ) ) 
    498          IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain 
    499          IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+ssha-SUM(e3t_a))) =', z_tmax 
    500          ! 
    501          zht(:,:) = 0.0_wp 
    502          DO jk = 1, jpkm1 
    503             zht(:,:) = zht(:,:) + e3t_b(:,:,jk) * tmask(:,:,jk) 
    504          END DO 
    505          z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + sshb(:,:) - zht(:,:) ) ) 
    506          IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain 
    507          IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+sshb-SUM(e3t_b))) =', z_tmax 
    508          ! 
    509          z_tmax = MAXVAL( tmask(:,:,1) *  ABS( sshb(:,:) ) ) 
    510          IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain 
    511          IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(sshb))) =', z_tmax 
    512          ! 
    513          z_tmax = MAXVAL( tmask(:,:,1) *  ABS( sshn(:,:) ) ) 
    514          IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain 
    515          IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(sshn))) =', z_tmax 
    516          ! 
    517          z_tmax = MAXVAL( tmask(:,:,1) *  ABS( ssha(:,:) ) ) 
    518          IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain 
    519          IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ssha))) =', z_tmax 
    520       END IF 
    521  
    522       ! *********************************** ! 
    523       ! After scale factors at u- v- points ! 
    524       ! *********************************** ! 
    525  
    526       CALL dom_vvl_interpol( e3t_a(:,:,:), e3u_a(:,:,:), 'U' ) 
    527       CALL dom_vvl_interpol( e3t_a(:,:,:), e3v_a(:,:,:), 'V' ) 
    528  
    529       ! *********************************** ! 
    530       ! After depths at u- v points         ! 
    531       ! *********************************** ! 
    532  
    533       hu_a(:,:) = e3u_a(:,:,1) * umask(:,:,1) 
    534       hv_a(:,:) = e3v_a(:,:,1) * vmask(:,:,1) 
    535       DO jk = 2, jpkm1 
    536          hu_a(:,:) = hu_a(:,:) + e3u_a(:,:,jk) * umask(:,:,jk) 
    537          hv_a(:,:) = hv_a(:,:) + e3v_a(:,:,jk) * vmask(:,:,jk) 
    538       END DO 
    539       !                                        ! Inverse of the local depth 
    540 !!gm BUG ?  don't understand the use of umask_i here ..... 
    541       r1_hu_a(:,:) = ssumask(:,:) / ( hu_a(:,:) + 1._wp - ssumask(:,:) ) 
    542       r1_hv_a(:,:) = ssvmask(:,:) / ( hv_a(:,:) + 1._wp - ssvmask(:,:) ) 
    543       ! 
    544       CALL wrk_dealloc( jpi,jpj,       zht, z_scale, zwu, zwv, zhdiv ) 
    545       CALL wrk_dealloc( jpi,jpj,jpk,   ze3t ) 
    546       ! 
    547       IF( nn_timing == 1 )  CALL timing_stop('dom_vvl_sf_nxt') 
    548       ! 
    549    END SUBROUTINE dom_vvl_sf_nxt 
    550  
    551  
    552    SUBROUTINE dom_vvl_sf_swp( kt ) 
    553       !!---------------------------------------------------------------------- 
    554       !!                ***  ROUTINE dom_vvl_sf_swp  *** 
    555       !!                    
    556       !! ** Purpose :  compute time filter and swap of scale factors  
    557       !!               compute all depths and related variables for next time step 
    558       !!               write outputs and restart file 
    559       !! 
    560       !! ** Method  :  - swap of e3t with trick for volume/tracer conservation 
    561       !!               - reconstruct scale factor at other grid points (interpolate) 
    562       !!               - recompute depths and water height fields 
    563       !! 
    564       !! ** Action  :  - e3t_(b/n), tilde_e3t_(b/n) and e3(u/v)_n ready for next time step 
    565       !!               - Recompute: 
    566       !!                    e3(u/v)_b        
    567       !!                    e3w_n            
    568       !!                    e3(u/v)w_b       
    569       !!                    e3(u/v)w_n       
    570       !!                    gdept_n, gdepw_n  and gde3w_n 
    571       !!                    h(u/v) and h(u/v)r 
    572       !! 
    573       !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling. 
    574       !!              Leclair, M., and G. Madec, 2011, Ocean Modelling. 
    575       !!---------------------------------------------------------------------- 
    576       INTEGER, INTENT( in ) ::   kt   ! time step 
    577       ! 
    578       INTEGER  ::   ji, jj, jk   ! dummy loop indices 
    579       REAL(wp) ::   zcoef        ! local scalar 
    580       !!---------------------------------------------------------------------- 
    581       ! 
    582       IF( ln_linssh )   RETURN      ! No calculation in linear free surface 
    583       ! 
    584       IF( nn_timing == 1 )  CALL timing_start('dom_vvl_sf_swp') 
    585       ! 
    586       IF( kt == nit000 )   THEN 
    587          IF(lwp) WRITE(numout,*) 
    588          IF(lwp) WRITE(numout,*) 'dom_vvl_sf_swp : - time filter and swap of scale factors' 
    589          IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~   - interpolate scale factors and compute depths for next time step' 
    590       ENDIF 
    591       ! 
    592       ! Time filter and swap of scale factors 
    593       ! ===================================== 
    594       ! - ML - e3(t/u/v)_b are allready computed in dynnxt. 
    595       IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN 
    596          IF( neuler == 0 .AND. kt == nit000 ) THEN 
    597             tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:) 
    598          ELSE 
    599             tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:) &  
    600             &         + atfp * ( tilde_e3t_b(:,:,:) - 2.0_wp * tilde_e3t_n(:,:,:) + tilde_e3t_a(:,:,:) ) 
    601          ENDIF 
    602          tilde_e3t_n(:,:,:) = tilde_e3t_a(:,:,:) 
    603       ENDIF 
    604       gdept_b(:,:,:) = gdept_n(:,:,:) 
    605       gdepw_b(:,:,:) = gdepw_n(:,:,:) 
    606  
    607       e3t_n(:,:,:) = e3t_a(:,:,:) 
    608       e3u_n(:,:,:) = e3u_a(:,:,:) 
    609       e3v_n(:,:,:) = e3v_a(:,:,:) 
    610  
    611       ! Compute all missing vertical scale factor and depths 
    612       ! ==================================================== 
    613       ! Horizontal scale factor interpolations 
    614       ! -------------------------------------- 
    615       ! - ML - e3u_b and e3v_b are allready computed in dynnxt 
    616       ! - JC - hu_b, hv_b, hur_b, hvr_b also 
    617        
    618       CALL dom_vvl_interpol( e3u_n(:,:,:), e3f_n(:,:,:), 'F'  ) 
    619        
    620       ! Vertical scale factor interpolations 
    621       CALL dom_vvl_interpol( e3t_n(:,:,:),  e3w_n(:,:,:), 'W'  ) 
    622       CALL dom_vvl_interpol( e3u_n(:,:,:), e3uw_n(:,:,:), 'UW' ) 
    623       CALL dom_vvl_interpol( e3v_n(:,:,:), e3vw_n(:,:,:), 'VW' ) 
    624       CALL dom_vvl_interpol( e3t_b(:,:,:),  e3w_b(:,:,:), 'W'  ) 
    625       CALL dom_vvl_interpol( e3u_b(:,:,:), e3uw_b(:,:,:), 'UW' ) 
    626       CALL dom_vvl_interpol( e3v_b(:,:,:), e3vw_b(:,:,:), 'VW' ) 
    627  
    628       ! t- and w- points depth (set the isf depth as it is in the initial step) 
    629       gdept_n(:,:,1) = 0.5_wp * e3w_n(:,:,1) 
    630       gdepw_n(:,:,1) = 0.0_wp 
    631       gde3w_n(:,:,1) = gdept_n(:,:,1) - sshn(:,:) 
    632       DO jk = 2, jpk 
    633          DO jj = 1,jpj 
    634             DO ji = 1,jpi 
    635               !    zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))   ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt 
    636                                                                  ! 1 for jk = mikt 
    637                zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk)) 
    638                gdepw_n(ji,jj,jk) = gdepw_n(ji,jj,jk-1) + e3t_n(ji,jj,jk-1) 
    639                gdept_n(ji,jj,jk) =    zcoef  * ( gdepw_n(ji,jj,jk  ) + 0.5 * e3w_n(ji,jj,jk) )  & 
    640                    &             + (1-zcoef) * ( gdept_n(ji,jj,jk-1) +       e3w_n(ji,jj,jk) )  
    641                gde3w_n(ji,jj,jk) = gdept_n(ji,jj,jk) - sshn(ji,jj) 
    642             END DO 
    643          END DO 
    644       END DO 
    645  
    646       ! Local depth and Inverse of the local depth of the water 
    647       ! ------------------------------------------------------- 
    648       hu_n(:,:) = hu_a(:,:)   ;   r1_hu_n(:,:) = r1_hu_a(:,:) 
    649       hv_n(:,:) = hv_a(:,:)   ;   r1_hv_n(:,:) = r1_hv_a(:,:) 
    650       ! 
    651       ht_n(:,:) = e3t_n(:,:,1) * tmask(:,:,1) 
    652       DO jk = 2, jpkm1 
    653          ht_n(:,:) = ht_n(:,:) + e3t_n(:,:,jk) * tmask(:,:,jk) 
    654       END DO 
    655  
    656       ! write restart file 
    657       ! ================== 
    658       IF( lrst_oce )   CALL dom_vvl_rst( kt, 'WRITE' ) 
    659       ! 
    660       IF( nn_timing == 1 )   CALL timing_stop('dom_vvl_sf_swp') 
    661       ! 
    662    END SUBROUTINE dom_vvl_sf_swp 
    663  
    664  
    665241   SUBROUTINE dom_vvl_interpol( pe3_in, pe3_out, pout ) 
    666242      !!--------------------------------------------------------------------- 
     
    684260      IF( nn_timing == 1 )   CALL timing_start('dom_vvl_interpol') 
    685261      ! 
    686       IF(ln_wd) THEN 
    687         zlnwd = 1.0_wp 
    688       ELSE 
    689         zlnwd = 0.0_wp 
    690       END IF 
     262      zlnwd = 0.0_wp 
    691263      ! 
    692264      SELECT CASE ( pout )    !==  type of interpolation  ==! 
     
    772344      ! 
    773345   END SUBROUTINE dom_vvl_interpol 
    774  
    775  
    776    SUBROUTINE dom_vvl_rst( kt, cdrw ) 
    777       !!--------------------------------------------------------------------- 
    778       !!                   ***  ROUTINE dom_vvl_rst  *** 
    779       !!                      
    780       !! ** Purpose :   Read or write VVL file in restart file 
    781       !! 
    782       !! ** Method  :   use of IOM library 
    783       !!                if the restart does not contain vertical scale factors, 
    784       !!                they are set to the _0 values 
    785       !!                if the restart does not contain vertical scale factors increments (z_tilde), 
    786       !!                they are set to 0. 
    787       !!---------------------------------------------------------------------- 
    788       INTEGER         , INTENT(in) ::   kt     ! ocean time-step 
    789       CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag 
    790       ! 
    791       INTEGER ::   ji, jj, jk 
    792       INTEGER ::   id1, id2, id3, id4, id5     ! local integers 
    793       !!---------------------------------------------------------------------- 
    794       ! 
    795       IF( nn_timing == 1 )  CALL timing_start('dom_vvl_rst') 
    796       IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise  
    797          !                                   ! =============== 
    798          IF( ln_rstart ) THEN                   !* Read the restart file 
    799             CALL rst_read_open                  !  open the restart file if necessary 
    800             CALL iom_get( numror, jpdom_autoglo, 'sshn'   , sshn    ) 
    801             ! 
    802             id1 = iom_varid( numror, 'e3t_b', ldstop = .FALSE. ) 
    803             id2 = iom_varid( numror, 'e3t_n', ldstop = .FALSE. ) 
    804             id3 = iom_varid( numror, 'tilde_e3t_b', ldstop = .FALSE. ) 
    805             id4 = iom_varid( numror, 'tilde_e3t_n', ldstop = .FALSE. ) 
    806             id5 = iom_varid( numror, 'hdiv_lf', ldstop = .FALSE. ) 
    807             !                             ! --------- ! 
    808             !                             ! all cases ! 
    809             !                             ! --------- ! 
    810             IF( MIN( id1, id2 ) > 0 ) THEN       ! all required arrays exist 
    811                CALL iom_get( numror, jpdom_autoglo, 'e3t_b', e3t_b(:,:,:) ) 
    812                CALL iom_get( numror, jpdom_autoglo, 'e3t_n', e3t_n(:,:,:) ) 
    813                ! needed to restart if land processor not computed  
    814                IF(lwp) write(numout,*) 'dom_vvl_rst : e3t_b and e3t_n found in restart files' 
    815                WHERE ( tmask(:,:,:) == 0.0_wp )  
    816                   e3t_n(:,:,:) = e3t_0(:,:,:) 
    817                   e3t_b(:,:,:) = e3t_0(:,:,:) 
    818                END WHERE 
    819                IF( neuler == 0 ) THEN 
    820                   e3t_b(:,:,:) = e3t_n(:,:,:) 
    821                ENDIF 
    822             ELSE IF( id1 > 0 ) THEN 
    823                IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_n not found in restart files' 
    824                IF(lwp) write(numout,*) 'e3t_n set equal to e3t_b.' 
    825                IF(lwp) write(numout,*) 'neuler is forced to 0' 
    826                CALL iom_get( numror, jpdom_autoglo, 'e3t_b', e3t_b(:,:,:) ) 
    827                e3t_n(:,:,:) = e3t_b(:,:,:) 
    828                neuler = 0 
    829             ELSE IF( id2 > 0 ) THEN 
    830                IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_b not found in restart files' 
    831                IF(lwp) write(numout,*) 'e3t_b set equal to e3t_n.' 
    832                IF(lwp) write(numout,*) 'neuler is forced to 0' 
    833                CALL iom_get( numror, jpdom_autoglo, 'e3t_n', e3t_n(:,:,:) ) 
    834                e3t_b(:,:,:) = e3t_n(:,:,:) 
    835                neuler = 0 
    836             ELSE 
    837                IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_n not found in restart file' 
    838                IF(lwp) write(numout,*) 'Compute scale factor from sshn' 
    839                IF(lwp) write(numout,*) 'neuler is forced to 0' 
    840                DO jk = 1, jpk 
    841                   e3t_n(:,:,jk) =  e3t_0(:,:,jk) * ( ht_0(:,:) + sshn(:,:) ) & 
    842                       &                          / ( ht_0(:,:) + 1._wp - ssmask(:,:) ) * tmask(:,:,jk)   & 
    843                       &          + e3t_0(:,:,jk)                               * (1._wp -tmask(:,:,jk)) 
    844                END DO 
    845                e3t_b(:,:,:) = e3t_n(:,:,:) 
    846                neuler = 0 
    847             ENDIF 
    848             !                             ! ----------- ! 
    849             IF( ln_vvl_zstar ) THEN       ! z_star case ! 
    850                !                          ! ----------- ! 
    851                IF( MIN( id3, id4 ) > 0 ) THEN 
    852                   CALL ctl_stop( 'dom_vvl_rst: z_star cannot restart from a z_tilde or layer run' ) 
    853                ENDIF 
    854                !                          ! ----------------------- ! 
    855             ELSE                          ! z_tilde and layer cases ! 
    856                !                          ! ----------------------- ! 
    857                IF( MIN( id3, id4 ) > 0 ) THEN  ! all required arrays exist 
    858                   CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_b', tilde_e3t_b(:,:,:) ) 
    859                   CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_n', tilde_e3t_n(:,:,:) ) 
    860                ELSE                            ! one at least array is missing 
    861                   tilde_e3t_b(:,:,:) = 0.0_wp 
    862                   tilde_e3t_n(:,:,:) = 0.0_wp 
    863                ENDIF 
    864                !                          ! ------------ ! 
    865                IF( ln_vvl_ztilde ) THEN   ! z_tilde case ! 
    866                   !                       ! ------------ ! 
    867                   IF( id5 > 0 ) THEN  ! required array exists 
    868                      CALL iom_get( numror, jpdom_autoglo, 'hdiv_lf', hdiv_lf(:,:,:) ) 
    869                   ELSE                ! array is missing 
    870                      hdiv_lf(:,:,:) = 0.0_wp 
    871                   ENDIF 
    872                ENDIF 
    873             ENDIF 
    874             ! 
    875          ELSE                                   !* Initialize at "rest" 
    876             e3t_b(:,:,:) = e3t_0(:,:,:) 
    877             e3t_n(:,:,:) = e3t_0(:,:,:) 
    878             sshn(:,:) = 0.0_wp 
    879  
    880             IF( ln_wd ) THEN 
    881               DO jj = 1, jpj 
    882                 DO ji = 1, jpi 
    883                   IF( e3t_0(ji,jj,1) <= 0.5_wp * rn_wdmin1 ) THEN 
    884                      e3t_b(ji,jj,:) = 0.5_wp * rn_wdmin1  
    885                      e3t_n(ji,jj,:) = 0.5_wp * rn_wdmin1  
    886                      e3t_a(ji,jj,:) = 0.5_wp * rn_wdmin1  
    887                      sshb(ji,jj) = rn_wdmin1 - bathy(ji,jj) 
    888                      sshn(ji,jj) = rn_wdmin1 - bathy(ji,jj) 
    889                      ssha(ji,jj) = rn_wdmin1 - bathy(ji,jj) 
    890                   ENDIF 
    891                 ENDDO 
    892               ENDDO 
    893             END IF 
    894  
    895             IF( ln_vvl_ztilde .OR. ln_vvl_layer) THEN 
    896                tilde_e3t_b(:,:,:) = 0.0_wp 
    897                tilde_e3t_n(:,:,:) = 0.0_wp 
    898                IF( ln_vvl_ztilde ) hdiv_lf(:,:,:) = 0.0_wp 
    899             END IF 
    900          ENDIF 
    901          ! 
    902       ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file 
    903          !                                   ! =================== 
    904          IF(lwp) WRITE(numout,*) '---- dom_vvl_rst ----' 
    905          !                                           ! --------- ! 
    906          !                                           ! all cases ! 
    907          !                                           ! --------- ! 
    908          CALL iom_rstput( kt, nitrst, numrow, 'e3t_b', e3t_b(:,:,:) ) 
    909          CALL iom_rstput( kt, nitrst, numrow, 'e3t_n', e3t_n(:,:,:) ) 
    910          !                                           ! ----------------------- ! 
    911          IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN  ! z_tilde and layer cases ! 
    912             !                                        ! ----------------------- ! 
    913             CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_b', tilde_e3t_b(:,:,:) ) 
    914             CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_n', tilde_e3t_n(:,:,:) ) 
    915          END IF 
    916          !                                           ! -------------!     
    917          IF( ln_vvl_ztilde ) THEN                    ! z_tilde case ! 
    918             !                                        ! ------------ ! 
    919             CALL iom_rstput( kt, nitrst, numrow, 'hdiv_lf', hdiv_lf(:,:,:) ) 
    920          ENDIF 
    921          ! 
    922       ENDIF 
    923       ! 
    924       IF( nn_timing == 1 )  CALL timing_stop('dom_vvl_rst') 
    925       ! 
    926    END SUBROUTINE dom_vvl_rst 
    927346 
    928347 
     
    990409      ! 
    991410      IF( ioptio /= 1 )   CALL ctl_stop( 'Choose ONE vertical coordinate in namelist nam_vvl' ) 
    992       IF( .NOT. ln_vvl_zstar .AND. ln_isf ) CALL ctl_stop( 'Only vvl_zstar has been tested with ice shelf cavity' ) 
    993411      ! 
    994412      IF(lwp) THEN                   ! Print the choice 
Note: See TracChangeset for help on using the changeset viewer.