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 13427 for NEMO/branches/2020/dev_r13327_KERNEL-06_2_techene_e3/src/OCE/DOM – NEMO

Ignore:
Timestamp:
2020-08-21T18:26:57+02:00 (4 years ago)
Author:
techene
Message:

debug in order to pass non linear SETTE test when agrif in not activated see #2385

Location:
NEMO/branches/2020/dev_r13327_KERNEL-06_2_techene_e3/src/OCE/DOM
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/dev_r13327_KERNEL-06_2_techene_e3/src/OCE/DOM/dom_oce.F90

    r13286 r13427  
    287287         &       e3w_0(jpi,jpj,jpk) , e3uw_0(jpi,jpj,jpk) , e3vw_0(jpi,jpj,jpk)                      ,  STAT=ierr(ii) ) 
    288288         ! 
    289 #if ! defined key_qco 
     289#if defined key_qco 
     290      ii = ii+1 
     291      ALLOCATE( r3t  (jpi,jpj,jpt)   , r3u  (jpi,jpj,jpt)    , r3v  (jpi,jpj,jpt)    , r3f  (jpi,jpj) ,      & 
     292         &      r3t_f(jpi,jpj)       , r3u_f(jpi,jpj)        , r3v_f(jpi,jpj)                         ,  STAT=ierr(ii) )              
     293#else 
    290294      ii = ii+1 
    291295      ALLOCATE( e3t(jpi,jpj,jpk,jpt) , e3u (jpi,jpj,jpk,jpt) , e3v (jpi,jpj,jpk,jpt) , e3f(jpi,jpj,jpk) ,      & 
     
    294298         ! 
    295299      ii = ii+1 
    296       ALLOCATE( r3t  (jpi,jpj,jpt)   , r3u  (jpi,jpj,jpt)    , r3v  (jpi,jpj,jpt)    , r3f  (jpi,jpj) ,  & 
    297          &      r3t_f(jpi,jpj)       , r3u_f(jpi,jpj)        , r3v_f(jpi,jpj)                         ,  STAT=ierr(ii) )        
    298          ! 
    299       ii = ii+1 
    300300      ALLOCATE( ht_0(jpi,jpj) ,    hu_0(jpi,jpj)    ,    hv_0(jpi,jpj)     , hf_0(jpi,jpj) ,       & 
    301301         &   r1_ht_0(jpi,jpj) , r1_hu_0(jpi,jpj) ,    r1_hv_0(jpi,jpj),   r1_hf_0(jpi,jpj) ,   STAT=ierr(ii)  ) 
     
    304304      ii = ii+1 
    305305      ALLOCATE( ht  (jpi,jpj) ,    hu  (jpi,jpj,jpt),    hv  (jpi,jpj,jpt)                 ,       & 
    306          &                      r1_hu  (jpi,jpj,jpt), r1_hv  (jpi,jpj,jpt)                 ,   STAT=ierr(ii)  ) 
    307 #else 
    308       ii = ii+1 
    309       ALLOCATE(                    hu  (jpi,jpj,jpt),    hv  (jpi,jpj,jpt)                 ,       & 
    310306         &                      r1_hu  (jpi,jpj,jpt), r1_hv  (jpi,jpj,jpt)                 ,   STAT=ierr(ii)  ) 
    311307#endif 
  • NEMO/branches/2020/dev_r13327_KERNEL-06_2_techene_e3/src/OCE/DOM/domain.F90

    r13286 r13427  
    177177      ! 
    178178      IF( ln_linssh ) THEN       != Fix in time : set to the reference one for all 
    179       ! 
     179         ! 
    180180         DO jt = 1, jpt                         ! depth of t- and w-grid-points 
    181181            gdept(:,:,:,jt) = gdept_0(:,:,:) 
     
    185185         ! 
    186186         DO jt = 1, jpt                         ! vertical scale factors 
    187             e3t(:,:,:,jt) =  e3t_0(:,:,:) 
    188             e3u(:,:,:,jt) =  e3u_0(:,:,:) 
    189             e3v(:,:,:,jt) =  e3v_0(:,:,:) 
    190             e3w(:,:,:,jt) =  e3w_0(:,:,:) 
     187            e3t (:,:,:,jt) =  e3t_0(:,:,:) 
     188            e3u (:,:,:,jt) =  e3u_0(:,:,:) 
     189            e3v (:,:,:,jt) =  e3v_0(:,:,:) 
     190            e3w (:,:,:,jt) =  e3w_0(:,:,:) 
    191191            e3uw(:,:,:,jt) = e3uw_0(:,:,:) 
    192192            e3vw(:,:,:,jt) = e3vw_0(:,:,:) 
    193193         END DO 
    194             e3f(:,:,:)    =  e3f_0(:,:,:) 
     194            e3f (:,:,:)    =  e3f_0(:,:,:) 
    195195         ! 
    196196         DO jt = 1, jpt                         ! water column thickness and its inverse 
    197             hu(:,:,jt)    =    hu_0(:,:) 
    198             hv(:,:,jt)    =    hv_0(:,:) 
     197               hu(:,:,jt) =    hu_0(:,:) 
     198               hv(:,:,jt) =    hv_0(:,:) 
    199199            r1_hu(:,:,jt) = r1_hu_0(:,:) 
    200200            r1_hv(:,:,jt) = r1_hv_0(:,:) 
    201201         END DO 
    202             ht(:,:) =    ht_0(:,:) 
     202               ht   (:,:) =    ht_0(:,:) 
    203203         ! 
    204204      ELSE                       != time varying : initialize before/now/after variables 
  • NEMO/branches/2020/dev_r13327_KERNEL-06_2_techene_e3/src/OCE/DOM/domqco.F90

    r13295 r13427  
    88   !!            3.3  !  2011-10  (M. Leclair) totally rewrote domvvl: vvl option includes z_star and z_tilde coordinates 
    99   !!            3.6  !  2014-11  (P. Mathiot) add ice shelf capability 
    10    !!            4.1  !  2019-08  (A. Coward, D. Storkey) rename dom_vvl_sf_swp -> dom_vvl_sf_update for new timestepping 
    11    !!            4.x  !  2020-02  (G. Madec, S. Techene) pure z* (quasi-eulerian) coordinate 
     10   !!            4.1  !  2019-08  (A. Coward, D. Storkey) add time level indices for prognostic variables 
     11   !!            4.x  !  2020-02  (S. Techene, G. Madec) quasi-eulerian coordinate (z* or s*) from domvvl 
    1212   !!---------------------------------------------------------------------- 
    1313 
    1414   !!---------------------------------------------------------------------- 
    15    !!   dom_qe_init   : define initial vertical scale factors, depths and column thickness 
    16    !!   dom_qe_r3c    : Compute ssh/h_0 ratioat t-, u-, v-, and optionally f-points 
    17    !!       qe_rst_read : read/write restart file 
    18    !!   dom_qe_ctl    : Check the vvl options 
     15   !!   dom_qco_init   : define initial vertical scale factors, depths and column thickness 
     16   !!   dom_qco_zgr    : Set ssh/h_0 ratio at t 
     17   !!   dom_qco_r3c    : Compute ssh/h_0 ratio at t-, u-, v-, and optionally f-points 
     18   !!       qco_rst_read : read/write restart file 
     19   !!       qco_ctl    : Check the vvl options 
    1920   !!---------------------------------------------------------------------- 
    2021   USE oce            ! ocean dynamics and tracers 
     
    7980      !! 
    8081      !!---------------------------------------------------------------------- 
    81       INTEGER, INTENT(in) :: Kbb, Kmm, Kaa 
     82      INTEGER, INTENT(in) ::   Kbb, Kmm, Kaa   ! time level indices 
     83      !!---------------------------------------------------------------------- 
    8284      ! 
    8385      IF(lwp) WRITE(numout,*) 
     
    8587      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~' 
    8688      ! 
    87       CALL dom_qco_ctl     ! choose vertical coordinate (z_star, z_tilde or layer) 
    88       ! 
    89       !                    ! Read or initialize e3t_(b/n), tilde_e3t_(b/n) and hdiv_lf 
    90       CALL qe_rst_read( nit000, Kbb, Kmm ) 
    91       ! 
    92       CALL dom_qco_zgr(Kbb, Kmm, Kaa) ! interpolation scale factor, depth and water column 
    93       ! 
    94       ! IF(lwxios) THEN   ! define variables in restart file when writing with XIOS 
    95       !    CALL iom_set_rstw_var_active('e3t_b') 
    96       !    CALL iom_set_rstw_var_active('e3t_n') 
    97       ! ENDIF 
     89      CALL qco_ctl                            ! choose vertical coordinate (z_star, z_tilde or layer) 
     90      !              ! CAUTION COM A METTRE !!! 
     91!!st      CALL qco_rst_read2( nit000, Kbb, Kmm )  ! Read or initialize ssh_(Kbb/Kmm) and r3 
     92!!st CAUTION if read2 removed change restart.F90 !  
     93      ! 
     94      CALL qco_rst_read( nit000, Kbb, Kmm )   ! Read or initialize ssh_(Kbb/Kmm) 
     95      ! 
     96      CALL dom_qco_zgr( Kbb, Kmm, Kaa )       ! interpolation scale factor, depth and water column 
     97      ! 
     98      IF(lwxios) THEN   ! define variables in restart file when writing with XIOS 
     99         CALL iom_set_rstw_var_active('sshb') 
     100         CALL iom_set_rstw_var_active('sshn') 
     101      ENDIF 
    98102      ! 
    99103   END SUBROUTINE dom_qco_init 
    100104 
    101105 
    102    SUBROUTINE dom_qco_zgr(Kbb, Kmm, Kaa) 
     106   SUBROUTINE dom_qco_zgr( Kbb, Kmm, Kaa ) 
    103107      !!---------------------------------------------------------------------- 
    104108      !!                ***  ROUTINE dom_qco_init  *** 
    105109      !! 
    106       !! ** Purpose :  Initialization of all ssh. to h._0 ratio 
    107       !! 
    108       !! ** Method  :  - interpolate scale factors 
     110      !! ** Purpose :  Initialization of all ssh./h._0 ratio 
     111      !! 
     112      !! ** Method  :  - call domqco using Kbb and Kmm 
    109113      !! 
    110114      !! ** Action  : - r3(t/u/v)_b 
    111115      !!              - r3(t/u/v/f)_n 
    112       !! 
    113       !! Reference  : Leclair, M., and G. Madec, 2011, Ocean Modelling. 
    114       !!---------------------------------------------------------------------- 
    115       INTEGER, INTENT(in) :: Kbb, Kmm, Kaa 
     116      !!---------------------------------------------------------------------- 
     117      INTEGER, INTENT(in) ::   Kbb, Kmm, Kaa   ! time level indices 
    116118      !!---------------------------------------------------------------------- 
    117119      ! 
     
    121123      CALL dom_qco_r3c( ssh(:,:,Kmm), r3t(:,:,Kmm), r3u(:,:,Kmm), r3v(:,:,Kmm), r3f(:,:) ) 
    122124      ! 
     125#if defined key_agrif 
     126      ! We need to define r3[tuv](Kaa) for AGRIF initialisation (should not be a problem for the restartability...) 
     127      r3t(:,:,Kaa) = r3t(:,:,Kmm) 
     128      r3u(:,:,Kaa) = r3u(:,:,Kmm) 
     129      r3v(:,:,Kaa) = r3v(:,:,Kmm) 
     130#endif 
     131      ! 
    123132   END SUBROUTINE dom_qco_zgr 
    124133 
     
    148157      !                                      !==  ratio at u-,v-point  ==! 
    149158      ! 
    150       IF( ln_dynadv_vec ) THEN                     !- Vector Form   (thickness weighted averaging) 
     159!!st      IF( ln_dynadv_vec ) THEN                     !- Vector Form   (thickness weighted averaging) 
     160#if ! defined key_qcoTest_FluxForm 
     161      !                                ! no 'key_qcoTest_FluxForm' : surface weighted ssh average 
    151162         DO_2D( 0, 0, 0, 0 ) 
    152163            pr3u(ji,jj) = 0.5_wp * (  e1e2t(ji  ,jj) * pssh(ji  ,jj)  & 
     
    155166               &                    + e1e2t(ji,jj+1) * pssh(ji,jj+1)  ) * r1_hv_0(ji,jj) * r1_e1e2v(ji,jj) 
    156167         END_2D 
    157       ELSE                                         !- Flux Form   (simple averaging) 
     168!!st      ELSE                                         !- Flux Form   (simple averaging) 
     169#else 
    158170         DO_2D( 0, 0, 0, 0 ) 
    159             pr3u(ji,jj) = 0.5_wp * (  pssh(ji  ,jj) + pssh(ji+1,jj)  ) * r1_hu_0(ji,jj) 
    160             pr3v(ji,jj) = 0.5_wp * (  pssh(ji,jj  ) + pssh(ji,jj+1)  ) * r1_hv_0(ji,jj) 
     171            pr3u(ji,jj) = 0.5_wp * (  pssh(ji,jj) + pssh(ji+1,jj  )  ) * r1_hu_0(ji,jj) 
     172            pr3v(ji,jj) = 0.5_wp * (  pssh(ji,jj) + pssh(ji  ,jj+1)  ) * r1_hv_0(ji,jj) 
    161173         END_2D 
    162       ENDIF 
     174!!st      ENDIF 
     175#endif          
    163176      ! 
    164177      IF( .NOT.PRESENT( pr3f ) ) THEN              !- lbc on ratio at u-, v-points only 
     
    168181      ELSE                                   !==  ratio at f-point  ==! 
    169182         ! 
    170          IF( ln_dynadv_vec )   THEN                !- Vector Form   (thickness weighted averaging) 
     183!!st         IF( ln_dynadv_vec )   THEN                !- Vector Form   (thickness weighted averaging) 
     184#if ! defined key_qcoTest_FluxForm 
     185         !                                ! no 'key_qcoTest_FluxForm' : surface weighted ssh average 
     186 
    171187            DO_2D( 1, 0, 1, 0 )                               ! start from 1 since lbc_lnk('F') doesn't update the 1st row/line 
    172188               pr3f(ji,jj) = 0.25_wp * (  e1e2t(ji  ,jj  ) * pssh(ji  ,jj  )  & 
     
    175191                  &                     + e1e2t(ji+1,jj+1) * pssh(ji+1,jj+1)  ) * r1_hf_0(ji,jj) * r1_e1e2f(ji,jj) 
    176192            END_2D 
    177          ELSE                                      !- Flux Form   (simple averaging) 
     193!!st         ELSE                                      !- Flux Form   (simple averaging) 
     194#else 
    178195            DO_2D( 1, 0, 1, 0 )                               ! start from 1 since lbc_lnk('F') doesn't update the 1st row/line 
    179                pr3f(ji,jj) = 0.25_wp * (  pssh(ji  ,jj  ) + pssh(ji+1,jj  )  & 
    180                   &                     + pssh(ji  ,jj+1) + pssh(ji+1,jj+1)  ) * r1_hf_0(ji,jj) 
     196               pr3f(ji,jj) = 0.25_wp * (  pssh(ji,jj  ) + pssh(ji+1,jj  )  & 
     197                  &                     + pssh(ji,jj+1) + pssh(ji+1,jj+1)  ) * r1_hf_0(ji,jj) 
    181198            END_2D 
    182          ENDIF 
     199!!st         ENDIF 
     200#endif 
    183201         !                                                 ! lbc on ratio at u-,v-,f-points 
    184202         CALL lbc_lnk_multi( 'dom_qco_r3c', pr3u, 'U', 1._wp, pr3v, 'V', 1._wp, pr3f, 'F', 1._wp ) 
     
    189207 
    190208 
    191    SUBROUTINE qe_rst_read( kt, Kbb, Kmm ) 
     209   SUBROUTINE qco_rst_read( kt, Kbb, Kmm ) 
    192210      !!--------------------------------------------------------------------- 
    193       !!                   ***  ROUTINE qe_rst_read  *** 
     211      !!                   ***  ROUTINE qco_rst_read  *** 
    194212      !! 
    195213      !! ** Purpose :   Read ssh in restart file 
     
    199217      !!                it is set to the _0 values. 
    200218      !!---------------------------------------------------------------------- 
    201       INTEGER         , INTENT(in) ::   kt        ! ocean time-step 
    202       INTEGER         , INTENT(in) ::   Kbb, Kmm  ! ocean time level indices 
     219      INTEGER, INTENT(in) ::   kt         ! ocean time-step 
     220      INTEGER, INTENT(in) ::   Kbb, Kmm   ! ocean time level indices 
    203221      ! 
    204222      INTEGER ::   ji, jj, jk 
     
    206224      !!---------------------------------------------------------------------- 
    207225      ! 
    208          IF( ln_rstart ) THEN                   !* Read the restart file 
    209             CALL rst_read_open                  !  open the restart file if necessary 
    210             ! 
    211             id1 = iom_varid( numror, 'sshb', ldstop = .FALSE. ) 
    212             id2 = iom_varid( numror, 'sshn', ldstop = .FALSE. ) 
    213             ! 
    214             !                             ! --------- ! 
    215             !                             ! all cases ! 
    216             !                             ! --------- ! 
    217             ! 
    218             IF( MIN( id1, id2 ) > 0 ) THEN       ! all required arrays exist 
    219                CALL iom_get( numror, jpdom_auto, 'sshb'   , ssh(:,:,Kbb), ldxios = lrxios    ) 
    220                CALL iom_get( numror, jpdom_auto, 'sshn'   , ssh(:,:,Kmm), ldxios = lrxios    ) 
    221                ! needed to restart if land processor not computed 
    222                IF(lwp) write(numout,*) 'qe_rst_read : ssh(:,:,Kbb) and ssh(:,:,Kmm) found in restart files' 
    223                WHERE ( ssmask(:,:) == 0.0_wp )   !!gm/st ==> sm should not be necessary on ssh when it was required on e3 
    224                   ssh(:,:,Kmm) = 0._wp 
    225                   ssh(:,:,Kbb) = 0._wp 
    226                END WHERE 
    227                IF( l_1st_euler ) THEN 
    228                   ssh(:,:,Kbb) = ssh(:,:,Kmm) 
    229                ENDIF 
    230             ELSE IF( id1 > 0 ) THEN 
    231                IF(lwp) write(numout,*) 'qe_rst_read WARNING : ssh(:,:,Kmm) not found in restart files' 
    232                IF(lwp) write(numout,*) 'sshn set equal to sshb.' 
    233                IF(lwp) write(numout,*) 'neuler is forced to 0' 
    234                CALL iom_get( numror, jpdom_auto, 'sshb', ssh(:,:,Kbb), ldxios = lrxios ) 
    235                ssh(:,:,Kmm) = ssh(:,:,Kbb) 
    236                l_1st_euler = .TRUE. 
    237             ELSE IF( id2 > 0 ) THEN 
    238                IF(lwp) write(numout,*) 'qe_rst_read WARNING : ssh(:,:,Kbb) not found in restart files' 
    239                IF(lwp) write(numout,*) 'sshb set equal to sshn.' 
    240                IF(lwp) write(numout,*) 'neuler is forced to 0' 
    241                CALL iom_get( numror, jpdom_auto, 'sshn', ssh(:,:,Kmm), ldxios = lrxios ) 
     226      IF( ln_rstart ) THEN                   !* Read the restart file 
     227         CALL rst_read_open                  !  open the restart file if necessary 
     228         ! 
     229         id1 = iom_varid( numror, 'sshb', ldstop = .FALSE. ) 
     230         id2 = iom_varid( numror, 'sshn', ldstop = .FALSE. ) 
     231         ! 
     232         !                             ! --------- ! 
     233         !                             ! all cases ! 
     234         !                             ! --------- ! 
     235         ! 
     236         IF( MIN( id1, id2 ) > 0 ) THEN       ! all required arrays exist 
     237            CALL iom_get( numror, jpdom_auto, 'sshb'   , ssh(:,:,Kbb), ldxios = lrxios    ) 
     238            CALL iom_get( numror, jpdom_auto, 'sshn'   , ssh(:,:,Kmm), ldxios = lrxios    ) 
     239            ! needed to restart if land processor not computed 
     240            IF(lwp) write(numout,*) 'qco_rst_read : ssh(:,:,Kbb) and ssh(:,:,Kmm) found in restart files' 
     241            !!WHERE ( ssmask(:,:) == 0.0_wp )   !!gm/st ==> sm should not be necessary on ssh while it was required on e3 
     242            !!   ssh(:,:,Kmm) = 0._wp 
     243            !!   ssh(:,:,Kbb) = 0._wp 
     244            !!END WHERE 
     245            IF( l_1st_euler ) THEN 
    242246               ssh(:,:,Kbb) = ssh(:,:,Kmm) 
    243                l_1st_euler = .TRUE. 
    244             ELSE 
    245                IF(lwp) write(numout,*) 'qe_rst_read WARNING : ssh(:,:,Kmm) not found in restart file' 
    246                IF(lwp) write(numout,*) 'ssh_b and ssh_n set to zero' 
    247                IF(lwp) write(numout,*) 'neuler is forced to 0' 
    248                ssh(:,:,:) = 0._wp 
    249                l_1st_euler = .TRUE. 
    250247            ENDIF 
    251             ! 
    252          ELSE                                   !* Initialize at "rest" 
    253             ! 
    254             IF( ll_wd ) THEN   ! MJB ll_wd edits start here - these are essential 
     248         ELSE IF( id1 > 0 ) THEN 
     249            IF(lwp) write(numout,*) 'qco_rst_read WARNING : ssh(:,:,Kmm) not found in restart files' 
     250            IF(lwp) write(numout,*) 'sshn set equal to sshb.' 
     251            IF(lwp) write(numout,*) 'neuler is forced to 0' 
     252            CALL iom_get( numror, jpdom_auto, 'sshb', ssh(:,:,Kbb), ldxios = lrxios ) 
     253            ssh(:,:,Kmm) = ssh(:,:,Kbb) 
     254            l_1st_euler = .TRUE. 
     255         ELSE IF( id2 > 0 ) THEN 
     256            IF(lwp) write(numout,*) 'qco_rst_read WARNING : ssh(:,:,Kbb) not found in restart files' 
     257            IF(lwp) write(numout,*) 'sshb set equal to sshn.' 
     258            IF(lwp) write(numout,*) 'neuler is forced to 0' 
     259            CALL iom_get( numror, jpdom_auto, 'sshn', ssh(:,:,Kmm), ldxios = lrxios ) 
     260            ssh(:,:,Kbb) = ssh(:,:,Kmm) 
     261            l_1st_euler = .TRUE. 
     262         ELSE 
     263            IF(lwp) write(numout,*) 'qco_rst_read WARNING : ssh(:,:,Kmm) not found in restart file' 
     264            IF(lwp) write(numout,*) 'ssh_b and ssh_n set to zero' 
     265            IF(lwp) write(numout,*) 'neuler is forced to 0' 
     266            ssh(:,:,:) = 0._wp 
     267            l_1st_euler = .TRUE. 
     268         ENDIF 
     269         ! 
     270      ELSE                                   !* Initialize at "rest" 
     271         ! 
     272         IF( ll_wd ) THEN   ! MJB ll_wd edits start here - these are essential 
     273            ! 
     274            IF( cn_cfg == 'wad' ) THEN            ! Wetting and drying test case 
     275               CALL usr_def_istate( gdept(:,:,:,Kbb), tmask, ts(:,:,:,:,Kbb), uu(:,:,:,Kbb), vv(:,:,:,Kbb), ssh(:,:,Kbb)  ) 
     276               ts (:,:,:,:,Kmm) = ts (:,:,:,:,Kbb)       ! set now values from to before ones 
     277               ssh(:,:    ,Kmm) = ssh(:,:    ,Kbb) 
     278               uu (:,:,:  ,Kmm) = uu (:,:,:  ,Kbb) 
     279               vv (:,:,:  ,Kmm) = vv (:,:,:  ,Kbb) 
     280            ELSE                                  ! if not test case 
     281               ssh(:,:,Kmm) = -ssh_ref 
     282               ssh(:,:,Kbb) = -ssh_ref 
    255283               ! 
    256                IF( cn_cfg == 'wad' ) THEN            ! Wetting and drying test case 
    257                   CALL usr_def_istate( gdept(:,:,:,Kbb), tmask, ts(:,:,:,:,Kbb), uu(:,:,:,Kbb), vv(:,:,:,Kbb), ssh(:,:,Kbb)  ) 
    258                   ts (:,:,:,:,Kmm) = ts (:,:,:,:,Kbb)       ! set now values from to before ones 
    259                   ssh(:,:    ,Kmm) = ssh(:,:    ,Kbb) 
    260                   uu (:,:,:  ,Kmm) = uu (:,:,:  ,Kbb) 
    261                   vv (:,:,:  ,Kmm) = vv (:,:,:  ,Kbb) 
    262                ELSE                                  ! if not test case 
    263                   ssh(:,:,Kmm) = -ssh_ref 
    264                   ssh(:,:,Kbb) = -ssh_ref 
    265                   ! 
    266                   DO_2D( 1, 1, 1, 1 ) 
    267                      IF( ht_0(ji,jj)-ssh_ref <  rn_wdmin1 ) THEN ! if total depth is less than min depth 
    268                         ssh(ji,jj,Kbb) = rn_wdmin1 - (ht_0(ji,jj) ) 
    269                         ssh(ji,jj,Kmm) = rn_wdmin1 - (ht_0(ji,jj) ) 
    270                      ENDIF 
    271                   END_2D 
    272                ENDIF 
    273  
    274                DO ji = 1, jpi 
    275                   DO jj = 1, jpj 
    276                      IF ( ht_0(ji,jj) .LE. 0.0 .AND. NINT( ssmask(ji,jj) ) .EQ. 1) THEN 
    277                        CALL ctl_stop( 'qe_rst_read: ht_0 must be positive at potentially wet points' ) 
    278                      ENDIF 
    279                   END DO 
     284               DO_2D( 1, 1, 1, 1 ) 
     285                  IF( ht_0(ji,jj)-ssh_ref <  rn_wdmin1 ) THEN ! if total depth is less than min depth 
     286                     ssh(ji,jj,Kbb) = rn_wdmin1 - (ht_0(ji,jj) ) 
     287                     ssh(ji,jj,Kmm) = rn_wdmin1 - (ht_0(ji,jj) ) 
     288                  ENDIF 
     289               END_2D 
     290            ENDIF 
     291            ! 
     292            DO ji = 1, jpi 
     293               DO jj = 1, jpj 
     294                  IF ( ht_0(ji,jj) .LE. 0.0 .AND. NINT( ssmask(ji,jj) ) .EQ. 1) THEN 
     295                    CALL ctl_stop( 'qco_rst_read: ht_0 must be positive at potentially wet points' ) 
     296                  ENDIF 
    280297               END DO 
     298            END DO 
     299            ! 
     300         ELSE 
     301            ! 
     302            ! Just to read set ssh in fact, called latter once vertical grid 
     303            ! is set up: 
     304!           CALL usr_def_istate( gdept_0, tmask, ts(:,:,:,:,Kbb), uu(:,:,:,Kbb), vv(:,:,:,Kbb), ssh(:,:,Kbb)  ) 
     305!           ! 
     306            ssh(:,:,:) = 0._wp 
     307            ! 
     308         ENDIF           ! end of ll_wd edits 
     309         ! 
     310      ENDIF 
     311      ! 
     312   END SUBROUTINE qco_rst_read 
     313 
     314    
     315   SUBROUTINE qco_rst_read2( kt, Kbb, Kmm ) 
     316      !!--------------------------------------------------------------------- 
     317      !!                   ***  ROUTINE qco_rst_read  *** 
     318      !! 
     319      !! ** Purpose :   Read ssh in restart file 
     320      !! 
     321      !! ** Method  :   use of IOM library 
     322      !!                if the restart does not contain ssh, 
     323      !!                it is set to the _0 values. 
     324      !!---------------------------------------------------------------------- 
     325      INTEGER, INTENT(in) ::   kt         ! ocean time-step 
     326      INTEGER, INTENT(in) ::   Kbb, Kmm   ! ocean time level indices 
     327      ! 
     328      INTEGER ::   ji, jj, jk 
     329      INTEGER ::   id1, id2     ! local integers 
     330      !!---------------------------------------------------------------------- 
     331      ! 
     332      IF( ln_rstart ) THEN                   !* Read the restart file 
     333         CALL rst_read_open                  !  open the restart file if necessary 
     334         ! 
     335         id1 = iom_varid( numror, 'sshb', ldstop = .FALSE. ) 
     336         id2 = iom_varid( numror, 'sshn', ldstop = .FALSE. ) 
     337         ! 
     338         !                             ! --------- ! 
     339         !                             ! all cases ! 
     340         !                             ! --------- ! 
     341         ! 
     342         IF( MIN( id1, id2 ) > 0 ) THEN       ! all required arrays exist 
     343            CALL iom_get( numror, jpdom_auto, 'sshb'   , ssh(:,:,Kbb), ldxios = lrxios    ) 
     344            CALL iom_get( numror, jpdom_auto, 'sshn'   , ssh(:,:,Kmm), ldxios = lrxios    ) 
     345            CALL iom_get( numror, jpdom_auto, 'r3tb'   , r3t(:,:,Kbb), ldxios = lrxios    ) 
     346            CALL iom_get( numror, jpdom_auto, 'r3tn'   , r3t(:,:,Kmm), ldxios = lrxios    ) 
     347            CALL iom_get( numror, jpdom_auto, 'r3ub'   , r3u(:,:,Kbb), ldxios = lrxios, cd_type = 'U'    ) 
     348            CALL iom_get( numror, jpdom_auto, 'r3un'   , r3u(:,:,Kmm), ldxios = lrxios, cd_type = 'U'    ) 
     349            CALL iom_get( numror, jpdom_auto, 'r3vb'   , r3v(:,:,Kbb), ldxios = lrxios, cd_type = 'V'    ) 
     350            CALL iom_get( numror, jpdom_auto, 'r3vn'   , r3v(:,:,Kmm), ldxios = lrxios, cd_type = 'V'    ) 
     351            CALL iom_get( numror, jpdom_auto, 'r3f'    , r3f(:,:)    , ldxios = lrxios, cd_type = 'F'    ) 
     352             
     353            ! needed to restart if land processor not computed 
     354            IF(lwp) write(numout,*) 'qco_rst_read : ssh(:,:,Kbb) and ssh(:,:,Kmm) found in restart files' 
     355            !!WHERE ( ssmask(:,:) == 0.0_wp )   !!gm/st ==> sm should not be necessary on ssh while it was required on e3 
     356            !!   ssh(:,:,Kmm) = 0._wp 
     357            !!   ssh(:,:,Kbb) = 0._wp 
     358            !!END WHERE 
     359            IF( l_1st_euler ) THEN 
     360               ssh(:,:,Kbb) = ssh(:,:,Kmm) 
     361            ENDIF 
     362         ELSE IF( id1 > 0 ) THEN 
     363            IF(lwp) write(numout,*) 'qco_rst_read WARNING : ssh(:,:,Kmm) not found in restart files' 
     364            IF(lwp) write(numout,*) 'sshn set equal to sshb.' 
     365            IF(lwp) write(numout,*) 'neuler is forced to 0' 
     366            CALL iom_get( numror, jpdom_auto, 'sshb', ssh(:,:,Kbb), ldxios = lrxios ) 
     367            ssh(:,:,Kmm) = ssh(:,:,Kbb) 
     368            l_1st_euler = .TRUE. 
     369         ELSE IF( id2 > 0 ) THEN 
     370            IF(lwp) write(numout,*) 'qco_rst_read WARNING : ssh(:,:,Kbb) not found in restart files' 
     371            IF(lwp) write(numout,*) 'sshb set equal to sshn.' 
     372            IF(lwp) write(numout,*) 'neuler is forced to 0' 
     373            CALL iom_get( numror, jpdom_auto, 'sshn', ssh(:,:,Kmm), ldxios = lrxios ) 
     374            ssh(:,:,Kbb) = ssh(:,:,Kmm) 
     375            l_1st_euler = .TRUE. 
     376         ELSE 
     377            IF(lwp) write(numout,*) 'qco_rst_read WARNING : ssh(:,:,Kmm) not found in restart file' 
     378            IF(lwp) write(numout,*) 'ssh_b and ssh_n set to zero' 
     379            IF(lwp) write(numout,*) 'neuler is forced to 0' 
     380            ssh(:,:,:) = 0._wp 
     381            l_1st_euler = .TRUE. 
     382         ENDIF 
     383         ! 
     384      ELSE                                   !* Initialize at "rest" 
     385         ! 
     386         IF( ll_wd ) THEN   ! MJB ll_wd edits start here - these are essential 
     387            ! 
     388            IF( cn_cfg == 'wad' ) THEN            ! Wetting and drying test case 
     389               CALL usr_def_istate( gdept(:,:,:,Kbb), tmask, ts(:,:,:,:,Kbb), uu(:,:,:,Kbb), vv(:,:,:,Kbb), ssh(:,:,Kbb)  ) 
     390               ts (:,:,:,:,Kmm) = ts (:,:,:,:,Kbb)       ! set now values from to before ones 
     391               ssh(:,:    ,Kmm) = ssh(:,:    ,Kbb) 
     392               uu (:,:,:  ,Kmm) = uu (:,:,:  ,Kbb) 
     393               vv (:,:,:  ,Kmm) = vv (:,:,:  ,Kbb) 
     394            ELSE                                  ! if not test case 
     395               ssh(:,:,Kmm) = -ssh_ref 
     396               ssh(:,:,Kbb) = -ssh_ref 
    281397               ! 
    282             ELSE 
    283                ! 
    284                ! Just to read set ssh in fact, called latter once vertical grid 
    285                ! is set up: 
    286 !               CALL usr_def_istate( gdept_0, tmask, ts(:,:,:,:,Kbb), uu(:,:,:,Kbb), vv(:,:,:,Kbb), ssh(:,:,Kbb)  ) 
    287 !               ! 
    288                 ssh(:,:,:) = 0._wp 
    289                ! 
    290             ENDIF           ! end of ll_wd edits 
    291             ! 
    292          ENDIF 
    293       ! 
    294    END SUBROUTINE qe_rst_read 
    295  
    296  
    297    SUBROUTINE dom_qco_ctl 
     398               DO_2D( 1, 1, 1, 1 ) 
     399                  IF( ht_0(ji,jj)-ssh_ref <  rn_wdmin1 ) THEN ! if total depth is less than min depth 
     400                     ssh(ji,jj,Kbb) = rn_wdmin1 - (ht_0(ji,jj) ) 
     401                     ssh(ji,jj,Kmm) = rn_wdmin1 - (ht_0(ji,jj) ) 
     402                  ENDIF 
     403               END_2D 
     404            ENDIF 
     405            ! 
     406            DO ji = 1, jpi 
     407               DO jj = 1, jpj 
     408                  IF ( ht_0(ji,jj) .LE. 0.0 .AND. NINT( ssmask(ji,jj) ) .EQ. 1) THEN 
     409                    CALL ctl_stop( 'qco_rst_read: ht_0 must be positive at potentially wet points' ) 
     410                  ENDIF 
     411               END DO 
     412            END DO 
     413            ! 
     414         ELSE 
     415            ! 
     416            ! Just to read set ssh in fact, called latter once vertical grid 
     417            ! is set up: 
     418!           CALL usr_def_istate( gdept_0, tmask, ts(:,:,:,:,Kbb), uu(:,:,:,Kbb), vv(:,:,:,Kbb), ssh(:,:,Kbb)  ) 
     419!           ! 
     420            ssh(:,:,:) = 0._wp 
     421            r3t(:,:,:) = 0._wp 
     422            r3u(:,:,:) = 0._wp 
     423            r3v(:,:,:) = 0._wp 
     424            r3f(:,:  ) = 0._wp 
     425            ! 
     426         ENDIF           ! end of ll_wd edits 
     427         ! 
     428      ENDIF 
     429      ! 
     430   END SUBROUTINE qco_rst_read2 
     431 
     432 
     433   SUBROUTINE qco_ctl 
    298434      !!--------------------------------------------------------------------- 
    299       !!                  ***  ROUTINE dom_qco_ctl  *** 
     435      !!                  ***  ROUTINE qco_ctl  *** 
    300436      !! 
    301437      !! ** Purpose :   Control the consistency between namelist options 
     
    317453      IF(lwp) THEN                    ! Namelist print 
    318454         WRITE(numout,*) 
    319          WRITE(numout,*) 'dom_qco_ctl : choice/control of the variable vertical coordinate' 
    320          WRITE(numout,*) '~~~~~~~~~~~' 
     455         WRITE(numout,*) 'qco_ctl : choice/control of the variable vertical coordinate' 
     456         WRITE(numout,*) '~~~~~~~~' 
    321457         WRITE(numout,*) '   Namelist nam_vvl : chose a vertical coordinate' 
    322458         WRITE(numout,*) '      zstar                      ln_vvl_zstar   = ', ln_vvl_zstar 
     
    362498#endif 
    363499      ! 
    364    END SUBROUTINE dom_qco_ctl 
     500   END SUBROUTINE qco_ctl 
    365501 
    366502   !!====================================================================== 
  • NEMO/branches/2020/dev_r13327_KERNEL-06_2_techene_e3/src/OCE/DOM/domzgr_substitute.h90

    r13237 r13427  
    1919#   define  e3uw(i,j,k,t)  (e3uw_0(i,j,k)*(1._wp+r3u(i,j,t))) 
    2020#   define  e3vw(i,j,k,t)  (e3vw_0(i,j,k)*(1._wp+r3v(i,j,t))) 
    21 #   define  ht(i,j)        (ht_0(i,j)+ssh(i,j,Kmm)) 
     21#   define  ht(i,j)        (ht_0(i,j)*(1._wp+r3t(i,j,Kmm))) 
    2222#   define  hu(i,j,t)      (hu_0(i,j)*(1._wp+r3u(i,j,t))) 
    2323#   define  hv(i,j,t)      (hv_0(i,j)*(1._wp+r3v(i,j,t))) 
     
    2929#endif 
    3030!!---------------------------------------------------------------------- 
     31!!#   define  e3t_f(i,j,k)   (e3t_0(i,j,k)*(1._wp+r3t_f(i,j)*tmask(i,j,k))) 
     32!!#   define  e3u_f(i,j,k)   (e3u_0(i,j,k)*(1._wp+r3u_f(i,j)*umask(i,j,k))) 
     33!!#   define  e3v_f(i,j,k)   (e3v_0(i,j,k)*(1._wp+r3v_f(i,j)*vmask(i,j,k))) 
  • NEMO/branches/2020/dev_r13327_KERNEL-06_2_techene_e3/src/OCE/DOM/istate.F90

    r13295 r13427  
    4242   PRIVATE 
    4343 
    44    PUBLIC   istate_init   ! routine called by step.F90 
     44   PUBLIC   istate_init   ! routine called by nemogcm.F90 
    4545 
    4646   !! * Substitutions 
     
    8080!!gm 
    8181 
    82       rhd  (:,:,:  ) = 0._wp   ;   rhop (:,:,:  ) = 0._wp      ! set one for all to 0 at level jpk 
    83       rn2b (:,:,:  ) = 0._wp   ;   rn2  (:,:,:  ) = 0._wp      ! set one for all to 0 at levels 1 and jpk 
    84       ts  (:,:,:,:,Kaa) = 0._wp                                   ! set one for all to 0 at level jpk 
    85       rab_b(:,:,:,:) = 0._wp   ;   rab_n(:,:,:,:) = 0._wp      ! set one for all to 0 at level jpk 
     82      rhd  (:,:,:      ) = 0._wp   ;   rhop (:,:,:  ) = 0._wp      ! set one for all to 0 at level jpk 
     83      rn2b (:,:,:      ) = 0._wp   ;   rn2  (:,:,:  ) = 0._wp      ! set one for all to 0 at levels 1 and jpk 
     84      ts   (:,:,:,:,Kaa) = 0._wp                                   ! set one for all to 0 at level jpk 
     85      rab_b(:,:,:,:    ) = 0._wp   ;   rab_n(:,:,:,:) = 0._wp      ! set one for all to 0 at level jpk 
    8686#if defined key_agrif 
    8787      uu   (:,:,:  ,Kaa) = 0._wp   ! used in agrif_oce_sponge at initialization 
     
    9696         CALL agrif_istate( Kbb, Kmm, Kaa )   ! Interp from parent 
    9797         ! 
    98          ts  (:,:,:,:,Kmm) = ts (:,:,:,:,Kbb)  
    99          ssh (:,:,Kmm)     = ssh(:,:,Kbb) 
    100          uu   (:,:,:,Kmm)   = uu  (:,:,:,Kbb) 
    101          vv   (:,:,:,Kmm)   = vv  (:,:,:,Kbb) 
     98         ts (:,:,:,:,Kmm) = ts (:,:,:,:,Kbb)  
     99         ssh(:,:,    Kmm) = ssh(:,:    ,Kbb) 
     100         uu (:,:,:  ,Kmm) = uu (:,:,:  ,Kbb) 
     101         vv (:,:,:  ,Kmm) = vv (:,:,:  ,Kbb) 
    102102      ELSE 
    103103#endif 
     
    117117            CALL dta_tsd( nit000, ts(:,:,:,:,Kbb) )       ! read 3D T and S data at nit000 
    118118            ! 
    119             ssh(:,:,Kbb)  = 0._wp               ! set the ocean at rest 
    120             uu  (:,:,:,Kbb) = 0._wp 
    121             vv  (:,:,:,Kbb) = 0._wp   
     119            ssh(:,:  ,Kbb) = 0._wp               ! set the ocean at rest 
     120            uu (:,:,:,Kbb) = 0._wp 
     121            vv (:,:,:,Kbb) = 0._wp   
    122122            ! 
    123123            IF( ll_wd ) THEN 
     
    139139            CALL usr_def_istate( zgdept, tmask, ts(:,:,:,:,Kbb), uu(:,:,:,Kbb), vv(:,:,:,Kbb), ssh(:,:,Kbb)  )          
    140140         ENDIF 
    141          ts  (:,:,:,:,Kmm) = ts (:,:,:,:,Kbb)       ! set now values from to before ones 
    142          ssh (:,:,Kmm)     = ssh(:,:,Kbb)    
    143          uu   (:,:,:,Kmm)   = uu  (:,:,:,Kbb) 
    144          vv   (:,:,:,Kmm)   = vv  (:,:,:,Kbb) 
     141         ts (:,:,:,:,Kmm) = ts (:,:,:,:,Kbb)       ! set now values from to before ones 
     142         ssh(:,:    ,Kmm) = ssh(:,:    ,Kbb)    
     143         uu (:,:,:  ,Kmm) = uu (:,:,:  ,Kbb) 
     144         vv (:,:,:  ,Kmm) = vv (:,:,:  ,Kbb) 
    145145 
    146146!!gm POTENTIAL BUG : 
Note: See TracChangeset for help on using the changeset viewer.