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 3862 – NEMO

Changeset 3862


Ignore:
Timestamp:
2013-04-08T13:07:43+02:00 (11 years ago)
Author:
acc
Message:

Branch 2013/dev_r3858_NOC_ZTC #863. Change names of gdept_0, gdepw_0, e3t_0, e3w_0 to gdept_1d, gdepw_1d, e3t_1d, e3w_1d respectively in preparation for the introduction of time invariant 3d arrays using the _0 suffix

Location:
branches/2013/dev_r3858_NOC_ZTC/NEMOGCM
Files:
33 edited

Legend:

Unmodified
Added
Removed
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OFF_SRC/domain.F90

    r3632 r3862  
    246246      !!      vertical scale factors. 
    247247      !! 
    248       !! ** Method  : - reference 1D vertical coordinate (gdep._0, e3._0) 
     248      !! ** Method  : - reference 1D vertical coordinate (gdep._1d, e3._1d) 
    249249      !!              - read/set ocean depth and ocean levels (bathy, mbathy) 
    250250      !!              - vertical coordinate (gdep., e3.) depending on the  
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OFF_SRC/domrea.F90

    r3680 r3862  
    178178            CALL iom_get( inum4, jpdom_data, 'e3w', e3w ) 
    179179 
    180             CALL iom_get( inum4, jpdom_unknown, 'gdept_0', gdept_0 ) ! depth 
    181             CALL iom_get( inum4, jpdom_unknown, 'gdepw_0', gdepw_0 ) 
     180            CALL iom_get( inum4, jpdom_unknown, 'gdept_1d', gdept_1d ) ! depth 
     181            CALL iom_get( inum4, jpdom_unknown, 'gdepw_1d', gdepw_1d ) 
    182182         ENDIF 
    183183 
    184184  
    185185         IF( ln_zps ) THEN                                           ! z-coordinate - partial steps 
    186             CALL iom_get( inum4, jpdom_unknown, 'gdept_0', gdept_0 )    ! reference depth 
    187             CALL iom_get( inum4, jpdom_unknown, 'gdepw_0', gdepw_0 ) 
    188             CALL iom_get( inum4, jpdom_unknown, 'e3t_0'  , e3t_0   )    ! reference scale factors 
    189             CALL iom_get( inum4, jpdom_unknown, 'e3w_0'  , e3w_0   ) 
     186            CALL iom_get( inum4, jpdom_unknown, 'gdept_1d', gdept_1d )  ! reference depth 
     187            CALL iom_get( inum4, jpdom_unknown, 'gdepw_1d', gdepw_1d ) 
     188            CALL iom_get( inum4, jpdom_unknown, 'e3t_1d'  , e3t_1d   )    ! reference scale factors 
     189            CALL iom_get( inum4, jpdom_unknown, 'e3w_1d'  , e3w_1d   ) 
    190190            ! 
    191191            IF( nmsh <= 6 ) THEN                                        ! 3D vertical scale factors 
     
    199199               !                                                        ! deduces the 3D scale factors 
    200200               DO jk = 1, jpk 
    201                   e3t(:,:,jk) = e3t_0(jk)                                     ! set to the ref. factors 
    202                   e3u(:,:,jk) = e3t_0(jk) 
    203                   e3v(:,:,jk) = e3t_0(jk) 
    204                   e3w(:,:,jk) = e3w_0(jk) 
     201                  e3t(:,:,jk) = e3t_1d(jk)                                    ! set to the ref. factors 
     202                  e3u(:,:,jk) = e3t_1d(jk) 
     203                  e3v(:,:,jk) = e3t_1d(jk) 
     204                  e3w(:,:,jk) = e3w_1d(jk) 
    205205               END DO 
    206206               DO jj = 1,jpj                                                  ! adjust the deepest values 
    207207                  DO ji = 1,jpi 
    208208                     ik = mbkt(ji,jj) 
    209                      e3t(ji,jj,ik) = e3tp(ji,jj) * tmask(ji,jj,1) + e3t_0(1) * ( 1._wp - tmask(ji,jj,1) ) 
    210                      e3w(ji,jj,ik) = e3wp(ji,jj) * tmask(ji,jj,1) + e3w_0(1) * ( 1._wp - tmask(ji,jj,1) ) 
     209                     e3t(ji,jj,ik) = e3tp(ji,jj) * tmask(ji,jj,1) + e3t_1d(1) * ( 1._wp - tmask(ji,jj,1) ) 
     210                     e3w(ji,jj,ik) = e3wp(ji,jj) * tmask(ji,jj,1) + e3w_1d(1) * ( 1._wp - tmask(ji,jj,1) ) 
    211211                  END DO 
    212212               END DO 
     
    223223               ! 
    224224               DO jk = 1, jpk                        ! set to z-scale factor if zero (i.e. along closed boundaries) 
    225                   WHERE( e3u(:,:,jk) == 0._wp )   e3u(:,:,jk) = e3t_0(jk) 
    226                   WHERE( e3v(:,:,jk) == 0._wp )   e3v(:,:,jk) = e3t_0(jk) 
     225                  WHERE( e3u(:,:,jk) == 0._wp )   e3u(:,:,jk) = e3t_1d(jk) 
     226                  WHERE( e3v(:,:,jk) == 0._wp )   e3v(:,:,jk) = e3t_1d(jk) 
    227227               END DO 
    228228            END IF 
     
    236236               ! 
    237237               DO jk = 1, jpk                                              ! deduces the 3D depth 
    238                   gdept(:,:,jk) = gdept_0(jk) 
    239                   gdepw(:,:,jk) = gdepw_0(jk) 
     238                  gdept(:,:,jk) = gdept_1d(jk) 
     239                  gdepw(:,:,jk) = gdepw_1d(jk) 
    240240               END DO 
    241241               DO jj = 1, jpj 
     
    254254 
    255255         IF( ln_zco ) THEN           ! Vertical coordinates and scales factors 
    256             CALL iom_get( inum4, jpdom_unknown, 'gdept_0', gdept_0 ) ! depth 
    257             CALL iom_get( inum4, jpdom_unknown, 'gdepw_0', gdepw_0 ) 
    258             CALL iom_get( inum4, jpdom_unknown, 'e3t_0'  , e3t_0   ) 
    259             CALL iom_get( inum4, jpdom_unknown, 'e3w_0'  , e3w_0   ) 
     256            CALL iom_get( inum4, jpdom_unknown, 'gdept_1d', gdept_1d ) ! depth 
     257            CALL iom_get( inum4, jpdom_unknown, 'gdepw_1d', gdepw_1d ) 
     258            CALL iom_get( inum4, jpdom_unknown, 'e3t_1d'  , e3t_1d   ) 
     259            CALL iom_get( inum4, jpdom_unknown, 'e3w_1d'  , e3w_1d   ) 
    260260            DO jk = 1, jpk 
    261                e3t  (:,:,jk) = e3t_0(jk)                                     ! set to the ref. factors 
    262                e3u  (:,:,jk) = e3t_0(jk) 
    263                e3v  (:,:,jk) = e3t_0(jk) 
    264                e3w  (:,:,jk) = e3w_0(jk) 
    265                gdept(:,:,jk) = gdept_0(jk) 
    266                gdepw(:,:,jk) = gdepw_0(jk) 
     261               e3t  (:,:,jk) = e3t_1d(jk)                              ! set to the ref. factors 
     262               e3u  (:,:,jk) = e3t_1d(jk) 
     263               e3v  (:,:,jk) = e3t_1d(jk) 
     264               e3w  (:,:,jk) = e3w_1d(jk) 
     265               gdept(:,:,jk) = gdept_1d(jk) 
     266               gdepw(:,:,jk) = gdepw_1d(jk) 
    267267            END DO 
    268268         ENDIF 
     
    270270!!gm BUG in s-coordinate this does not work! 
    271271      ! deepest/shallowest W level Above/Below ~10m 
    272       zrefdep = 10._wp - ( 0.1_wp * MINVAL(e3w_0) )                  ! ref. depth with tolerance (10% of minimum layer thickness) 
    273       nlb10 = MINLOC( gdepw_0, mask = gdepw_0 > zrefdep, dim = 1 )  ! shallowest W level Below ~10m 
     272      zrefdep = 10._wp - ( 0.1_wp * MINVAL(e3w_1d) )                 ! ref. depth with tolerance (10% of minimum layer thickness) 
     273      nlb10 = MINLOC( gdepw_1d, mask = gdepw_1d > zrefdep, dim = 1 ) ! shallowest W level Below ~10m 
    274274      nla10 = nlb10 - 1                                              ! deepest    W level Above ~10m 
    275275!!gm end bug 
     
    312312         WRITE(numout,*) '              Reference z-coordinate depth and scale factors:' 
    313313         WRITE(numout, "(9x,' level   gdept    gdepw     e3t      e3w  ')" ) 
    314          WRITE(numout, "(10x, i4, 4f9.2)" ) ( jk, gdept_0(jk), gdepw_0(jk), e3t_0(jk), e3w_0(jk), jk = 1, jpk ) 
     314         WRITE(numout, "(10x, i4, 4f9.2)" ) ( jk, gdept_1d(jk), gdepw_1d(jk), e3t_1d(jk), e3w_1d(jk), jk = 1, jpk ) 
    315315      ENDIF 
    316316 
    317317      DO jk = 1, jpk 
    318          IF( e3w_0  (jk) <= 0._wp .OR. e3t_0  (jk) <= 0._wp )   CALL ctl_stop( ' e3w_0 or e3t_0 =< 0 ' ) 
    319          IF( gdepw_0(jk) <  0._wp .OR. gdept_0(jk) <  0._wp )   CALL ctl_stop( ' gdepw_0 or gdept_0 < 0 ' ) 
     318         IF( e3w_1d  (jk) <= 0._wp .OR. e3t_1d  (jk) <= 0._wp )   CALL ctl_stop( ' e3w_1d or e3t_1d =< 0 ' ) 
     319         IF( gdepw_1d(jk) <  0._wp .OR. gdept_1d(jk) <  0._wp )   CALL ctl_stop( ' gdepw_1d or gdept_1d < 0 ' ) 
    320320      END DO 
    321321      !                                     ! ============================ 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90

    r3294 r3862  
    212212               ik = mbkt(ji,jj) 
    213213               IF( ik > 1 ) THEN 
    214                   zztmp = ( gdept_0(ik) - fsdept_0(ji,jj,ik) ) / ( gdept_0(ik) - gdept_0(ik-1) ) 
     214                  zztmp = ( gdept_1d(ik) - fsdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) ) 
    215215                  sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1) 
    216216               ENDIF 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DIA/diadimg.F90

    r3294 r3862  
    112112 
    113113    CASE ( 'T') 
    114        z4dep(:)=gdept_0(:) 
     114       z4dep(:)=gdept_1d(:) 
    115115 
    116116    CASE ( 'W' ) 
    117        z4dep(:)=gdepw_0(:) 
     117       z4dep(:)=gdepw_1d(:) 
    118118 
    119119    CASE ( '2' ) 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DIA/diahth.F90

    r3764 r3862  
    304304      ! ----------------------------- ! 
    305305 
    306       ! find ilevel with (ilevel+1) the deepest W-level above 300m (we assume we can use e3t_0 to do this search...) 
     306      ! find ilevel with (ilevel+1) the deepest W-level above 300m (we assume we can use e3t_1d to do this search...) 
    307307      ilevel   = 0 
    308308      zthick_0 = 0._wp 
    309309      DO jk = 1, jpkm1                       
    310          zthick_0 = zthick_0 + e3t_0(jk) 
     310         zthick_0 = zthick_0 + e3t_1d(jk) 
    311311         IF( zthick_0 < 300. )   ilevel = jk 
    312312      END DO 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DIA/diaptr.F90

    r3764 r3862  
    662662            CALL histbeg(clhstnam, 1, zfoo, jpj, zphi,   & 
    663663               1, 1, 1, jpj, niter, zjulian, zdt*nn_fptr, nhoridz, numptr, domain_id=nidom_ptr) 
    664             ! Vertical grids : gdept_0, gdepw_0 
     664            ! Vertical grids : gdept_1d, gdepw_1d 
    665665            CALL histvert( numptr, "deptht", "Vertical T levels",   & 
    666                &                   "m", jpk, gdept_0, ndepidzt, "down" ) 
     666               &                   "m", jpk, gdept_1d, ndepidzt, "down" ) 
    667667            CALL histvert( numptr, "depthw", "Vertical W levels",   & 
    668                &                   "m", jpk, gdepw_0, ndepidzw, "down" ) 
     668               &                   "m", jpk, gdepw_1d, ndepidzw, "down" ) 
    669669            ! 
    670670            CALL wheneq ( jpj*jpk, MIN(sjk(:,:,1), 1._wp), 1, 1., ndex  , ndim  )      ! Lat-Depth 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90

    r3704 r3862  
    323323            &          nit000-1, zjulian, zdt, nh_T, nid_T, domain_id=nidom, snc4chunks=snc4set ) 
    324324         CALL histvert( nid_T, "deptht", "Vertical T levels",      &  ! Vertical grid: gdept 
    325             &           "m", ipk, gdept_0, nz_T, "down" ) 
     325            &           "m", ipk, gdept_1d, nz_T, "down" ) 
    326326         !                                                            ! Index of ocean points 
    327327         CALL wheneq( jpi*jpj*ipk, tmask, 1, 1., ndex_T , ndim_T  )      ! volume 
     
    359359            &          nit000-1, zjulian, zdt, nh_U, nid_U, domain_id=nidom, snc4chunks=snc4set ) 
    360360         CALL histvert( nid_U, "depthu", "Vertical U levels",      &  ! Vertical grid: gdept 
    361             &           "m", ipk, gdept_0, nz_U, "down" ) 
     361            &           "m", ipk, gdept_1d, nz_U, "down" ) 
    362362         !                                                            ! Index of ocean points 
    363363         CALL wheneq( jpi*jpj*ipk, umask, 1, 1., ndex_U , ndim_U  )      ! volume 
     
    372372            &          nit000-1, zjulian, zdt, nh_V, nid_V, domain_id=nidom, snc4chunks=snc4set ) 
    373373         CALL histvert( nid_V, "depthv", "Vertical V levels",      &  ! Vertical grid : gdept 
    374             &          "m", ipk, gdept_0, nz_V, "down" ) 
     374            &          "m", ipk, gdept_1d, nz_V, "down" ) 
    375375         !                                                            ! Index of ocean points 
    376376         CALL wheneq( jpi*jpj*ipk, vmask, 1, 1., ndex_V , ndim_V  )      ! volume 
     
    385385            &          nit000-1, zjulian, zdt, nh_W, nid_W, domain_id=nidom, snc4chunks=snc4set ) 
    386386         CALL histvert( nid_W, "depthw", "Vertical W levels",      &  ! Vertical grid: gdepw 
    387             &          "m", ipk, gdepw_0, nz_W, "down" ) 
     387            &          "m", ipk, gdepw_1d, nz_W, "down" ) 
    388388 
    389389 
     
    811811          1, jpi, 1, jpj, nit000-1, zjulian, zdt, nh_i, id_i, domain_id=nidom, snc4chunks=snc4set ) ! Horizontal grid : glamt and gphit 
    812812      CALL histvert( id_i, "deptht", "Vertical T levels",   &    ! Vertical grid : gdept 
    813           "m", jpk, gdept_0, nz_i, "down") 
     813          "m", jpk, gdept_1d, nz_i, "down") 
    814814 
    815815      ! Declare all the output fields as NetCDF variables 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DOM/dom_oce.F90

    r3851 r3862  
    165165   !! z-coordinate with full steps (also used in the other cases as reference z-coordinate) 
    166166   !! =-----------------====------ 
    167    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:)   :: gdept_0, gdepw_0 !: reference depth of t- and w-points (m) 
    168    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:)   :: e3t_0  , e3w_0   !: reference vertical scale factors at T- and W-pts (m) 
    169    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e3tp   , e3wp    !: ocean bottom level thickness at T and W points 
     167   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:)   :: gdept_1d, gdepw_1d !: reference depth of t- and w-points (m) 
     168   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:)   :: e3t_1d  , e3w_1d   !: reference vertical scale factors at T- and W-pts (m) 
     169   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: e3tp    , e3wp     !: ocean bottom level thickness at T and W points 
    170170 
    171171   !! s-coordinate and hybrid z-s-coordinate 
     
    293293         &      hv(jpi,jpj) , hvr(jpi,jpj) , hv_0(jpi,jpj) , STAT=ierr(6) ) 
    294294         ! 
    295       ALLOCATE( gdept_0(jpk) , gdepw_0(jpk) ,                                     & 
    296          &      e3t_0  (jpk) , e3w_0  (jpk) , e3tp (jpi,jpj), e3wp(jpi,jpj) ,     & 
     295      ALLOCATE( gdept_1d(jpk) , gdepw_1d(jpk) ,                                   & 
     296         &      e3t_1d  (jpk) , e3w_1d  (jpk) , e3tp (jpi,jpj), e3wp(jpi,jpj) ,   & 
    297297         &      gsigt  (jpk) , gsigw  (jpk) , gsi3w(jpk)    ,                     & 
    298298         &      esigt  (jpk) , esigw  (jpk)                                 , STAT=ierr(7) ) 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DOM/domstp.F90

    r2715 r3862  
    9090 
    9191            DO jk = 1, jpk 
    92                IF( gdept_0(jk) <= rdth ) rdttra(jk) = rdtmin 
    93                IF( gdept_0(jk) >  rdth ) THEN 
     92               IF( gdept_1d(jk) <= rdth ) rdttra(jk) = rdtmin 
     93               IF( gdept_1d(jk) >  rdth ) THEN 
    9494                  rdttra(jk) = rdtmin + ( rdtmax - rdtmin )   & 
    95                                       * ( EXP( ( gdept_0(jk ) - rdth ) / rdth ) - 1. )   & 
    96                                       / ( EXP( ( gdept_0(jpk) - rdth ) / rdth ) - 1. ) 
     95                                      * ( EXP( ( gdept_1d(jk ) - rdth ) / rdth ) - 1. )   & 
     96                                      / ( EXP( ( gdept_1d(jpk) - rdth ) / rdth ) - 1. ) 
    9797               ENDIF 
    9898               IF(lwp) WRITE(numout,"(36x,f5.2,5x,i3)") rdttra(jk)/3600., jk 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DOM/domwri.F90

    r3680 r3862  
    236236         ENDIF 
    237237         ! 
    238          CALL iom_rstput( 0, 0, inum4, 'gdept_0', gdept_0 )     !    ! reference z-coord. 
    239          CALL iom_rstput( 0, 0, inum4, 'gdepw_0', gdepw_0 ) 
    240          CALL iom_rstput( 0, 0, inum4, 'e3t_0'  , e3t_0   ) 
    241          CALL iom_rstput( 0, 0, inum4, 'e3w_0'  , e3w_0   ) 
     238         CALL iom_rstput( 0, 0, inum4, 'gdept_1d', gdept_1d )   !    ! reference z-coord. 
     239         CALL iom_rstput( 0, 0, inum4, 'gdepw_1d', gdepw_1d ) 
     240         CALL iom_rstput( 0, 0, inum4, 'e3t_1d'  , e3t_1d   ) 
     241         CALL iom_rstput( 0, 0, inum4, 'e3w_1d'  , e3w_1d   ) 
    242242      ENDIF 
    243243       
    244244      IF( ln_zco ) THEN 
    245245         !                                                      ! z-coordinate - full steps 
    246          CALL iom_rstput( 0, 0, inum4, 'gdept_0', gdept_0 )     !    ! depth 
    247          CALL iom_rstput( 0, 0, inum4, 'gdepw_0', gdepw_0 ) 
    248          CALL iom_rstput( 0, 0, inum4, 'e3t_0'  , e3t_0   )     !    ! scale factors 
    249          CALL iom_rstput( 0, 0, inum4, 'e3w_0'  , e3w_0   ) 
     246         CALL iom_rstput( 0, 0, inum4, 'gdept_1d', gdept_1d )   !    ! depth 
     247         CALL iom_rstput( 0, 0, inum4, 'gdepw_1d', gdepw_1d ) 
     248         CALL iom_rstput( 0, 0, inum4, 'e3t_1d'  , e3t_1d   )     !    ! scale factors 
     249         CALL iom_rstput( 0, 0, inum4, 'e3w_1d'  , e3w_1d   ) 
    250250      ENDIF 
    251251      !                                     ! ============================ 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DOM/domzgr.F90

    r3764 r3862  
    8888      !!              vertical scale factors. 
    8989      !! 
    90       !! ** Method  : - reference 1D vertical coordinate (gdep._0, e3._0) 
     90      !! ** Method  : - reference 1D vertical coordinate (gdep._1d, e3._1d) 
    9191      !!              - read/set ocean depth and ocean levels (bathy, mbathy) 
    9292      !!              - vertical coordinate (gdep., e3.) depending on the  
     
    177177      !!      function the derivative of which gives the scale factors. 
    178178      !!        both depth and scale factors only depend on k (1d arrays). 
    179       !!              w-level: gdepw_0  = fsdep(k) 
    180       !!                       e3w_0(k) = dk(fsdep)(k)     = fse3(k) 
    181       !!              t-level: gdept_0  = fsdep(k+0.5) 
    182       !!                       e3t_0(k) = dk(fsdep)(k+0.5) = fse3(k+0.5) 
    183       !! 
    184       !! ** Action  : - gdept_0, gdepw_0 : depth of T- and W-point (m) 
    185       !!              - e3t_0  , e3w_0   : scale factors at T- and W-levels (m) 
     179      !!              w-level: gdepw_1d  = fsdep(k) 
     180      !!                       e3w_1d(k) = dk(fsdep)(k)     = fse3(k) 
     181      !!              t-level: gdept_1d  = fsdep(k+0.5) 
     182      !!                       e3t_1d(k) = dk(fsdep)(k+0.5) = fse3(k+0.5) 
     183      !! 
     184      !! ** Action  : - gdept_1d, gdepw_1d : depth of T- and W-point (m) 
     185      !!              - e3t_1d  , e3w_1d   : scale factors at T- and W-levels (m) 
    186186      !! 
    187187      !! Reference : Marti, Madec & Delecluse, 1992, JGR, 97, No8, 12,763-12,766. 
     
    255255            zw = FLOAT( jk ) 
    256256            zt = FLOAT( jk ) + 0.5_wp 
    257             gdepw_0(jk) = ( zw - 1 ) * za1 
    258             gdept_0(jk) = ( zt - 1 ) * za1 
    259             e3w_0  (jk) =  za1 
    260             e3t_0  (jk) =  za1 
     257            gdepw_1d(jk) = ( zw - 1 ) * za1 
     258            gdept_1d(jk) = ( zt - 1 ) * za1 
     259            e3w_1d  (jk) =  za1 
     260            e3t_1d  (jk) =  za1 
    261261         END DO 
    262262      ELSE                                ! Madec & Imbard 1996 function 
     
    265265               zw = REAL( jk , wp ) 
    266266               zt = REAL( jk , wp ) + 0.5_wp 
    267                gdepw_0(jk) = ( zsur + za0 * zw + za1 * zacr * LOG ( COSH( (zw-zkth) / zacr ) )  ) 
    268                gdept_0(jk) = ( zsur + za0 * zt + za1 * zacr * LOG ( COSH( (zt-zkth) / zacr ) )  ) 
    269                e3w_0  (jk) =          za0      + za1        * TANH(       (zw-zkth) / zacr   ) 
    270                e3t_0  (jk) =          za0      + za1        * TANH(       (zt-zkth) / zacr   ) 
     267               gdepw_1d(jk) = ( zsur + za0 * zw + za1 * zacr * LOG ( COSH( (zw-zkth) / zacr ) )  ) 
     268               gdept_1d(jk) = ( zsur + za0 * zt + za1 * zacr * LOG ( COSH( (zt-zkth) / zacr ) )  ) 
     269               e3w_1d  (jk) =          za0      + za1        * TANH(       (zw-zkth) / zacr   ) 
     270               e3t_1d  (jk) =          za0      + za1        * TANH(       (zt-zkth) / zacr   ) 
    271271            END DO 
    272272         ELSE 
     
    275275               zt = FLOAT( jk ) + 0.5_wp 
    276276               ! Double tanh function 
    277                gdepw_0(jk) = ( zsur + za0 * zw + za1 * zacr * LOG ( COSH( (zw-zkth ) / zacr  ) )    & 
    278                   &                            + za2 * zacr2* LOG ( COSH( (zw-zkth2) / zacr2 ) )  ) 
    279                gdept_0(jk) = ( zsur + za0 * zt + za1 * zacr * LOG ( COSH( (zt-zkth ) / zacr  ) )    & 
    280                   &                            + za2 * zacr2* LOG ( COSH( (zt-zkth2) / zacr2 ) )  ) 
    281                e3w_0  (jk) =          za0      + za1        * TANH(       (zw-zkth ) / zacr  )    & 
    282                   &                            + za2        * TANH(       (zw-zkth2) / zacr2 ) 
    283                e3t_0  (jk) =          za0      + za1        * TANH(       (zt-zkth ) / zacr  )    & 
    284                   &                            + za2        * TANH(       (zt-zkth2) / zacr2 ) 
     277               gdepw_1d(jk) = ( zsur + za0 * zw + za1 * zacr * LOG ( COSH( (zw-zkth ) / zacr  ) )    & 
     278                  &                             + za2 * zacr2* LOG ( COSH( (zw-zkth2) / zacr2 ) )  ) 
     279               gdept_1d(jk) = ( zsur + za0 * zt + za1 * zacr * LOG ( COSH( (zt-zkth ) / zacr  ) )    & 
     280                  &                             + za2 * zacr2* LOG ( COSH( (zt-zkth2) / zacr2 ) )  ) 
     281               e3w_1d  (jk) =          za0      + za1        * TANH(       (zw-zkth ) / zacr  )      & 
     282                  &                             + za2        * TANH(       (zw-zkth2) / zacr2 ) 
     283               e3t_1d  (jk) =          za0      + za1        * TANH(       (zt-zkth ) / zacr  )      & 
     284                  &                             + za2        * TANH(       (zt-zkth2) / zacr2 ) 
    285285            END DO 
    286286         ENDIF 
    287          gdepw_0(1) = 0._wp                    ! force first w-level to be exactly at zero 
     287         gdepw_1d(1) = 0._wp                    ! force first w-level to be exactly at zero 
    288288      ENDIF 
    289289 
    290290!!gm BUG in s-coordinate this does not work! 
    291291      ! deepest/shallowest W level Above/Below ~10m 
    292       zrefdep = 10._wp - 0.1_wp * MINVAL( e3w_0 )                    ! ref. depth with tolerance (10% of minimum layer thickness) 
    293       nlb10 = MINLOC( gdepw_0, mask = gdepw_0 > zrefdep, dim = 1 )  ! shallowest W level Below ~10m 
     292      zrefdep = 10._wp - 0.1_wp * MINVAL( e3w_1d )                   ! ref. depth with tolerance (10% of minimum layer thickness) 
     293      nlb10 = MINLOC( gdepw_1d, mask = gdepw_1d > zrefdep, dim = 1 ) ! shallowest W level Below ~10m 
    294294      nla10 = nlb10 - 1                                              ! deepest    W level Above ~10m 
    295295!!gm end bug 
     
    299299         WRITE(numout,*) '              Reference z-coordinate depth and scale factors:' 
    300300         WRITE(numout, "(9x,' level   gdept    gdepw     e3t      e3w  ')" ) 
    301          WRITE(numout, "(10x, i4, 4f9.2)" ) ( jk, gdept_0(jk), gdepw_0(jk), e3t_0(jk), e3w_0(jk), jk = 1, jpk ) 
     301         WRITE(numout, "(10x, i4, 4f9.2)" ) ( jk, gdept_1d(jk), gdepw_1d(jk), e3t_1d(jk), e3w_1d(jk), jk = 1, jpk ) 
    302302      ENDIF 
    303303      DO jk = 1, jpk                      ! control positivity 
    304          IF( e3w_0  (jk) <= 0._wp .OR. e3t_0  (jk) <= 0._wp )   CALL ctl_stop( 'dom:zgr_z: e3w or e3t =< 0 '    ) 
    305          IF( gdepw_0(jk) <  0._wp .OR. gdept_0(jk) <  0._wp )   CALL ctl_stop( 'dom:zgr_z: gdepw or gdept < 0 ' ) 
     304         IF( e3w_1d  (jk) <= 0._wp .OR. e3t_1d  (jk) <= 0._wp )   CALL ctl_stop( 'dom:zgr_z: e3w or e3t =< 0 '    ) 
     305         IF( gdepw_1d(jk) <  0._wp .OR. gdept_1d(jk) <  0._wp )   CALL ctl_stop( 'dom:zgr_z: gdepw or gdept < 0 ' ) 
    306306      END DO 
    307307      ! 
     
    369369            IF(lwp) WRITE(numout,*) '         bathymetry field: flat basin' 
    370370            idta(:,:) = jpkm1                            ! before last level 
    371             zdta(:,:) = gdepw_0(jpk)                     ! last w-point depth 
    372             h_oce     = gdepw_0(jpk) 
     371            zdta(:,:) = gdepw_1d(jpk)                     ! last w-point depth 
     372            h_oce     = gdepw_1d(jpk) 
    373373         ELSE                                         ! bump centered in the basin 
    374374            IF(lwp) WRITE(numout,*) 
     
    378378            r_bump  = 50000._wp                            ! bump radius (meters)        
    379379            h_bump  =  2700._wp                            ! bump height (meters) 
    380             h_oce   = gdepw_0(jpk)                         ! background ocean depth (meters) 
     380            h_oce   = gdepw_1d(jpk)                         ! background ocean depth (meters) 
    381381            IF(lwp) WRITE(numout,*) '            bump characteristics: ' 
    382382            IF(lwp) WRITE(numout,*) '               bump center (i,j)   = ', ii_bump, ii_bump 
     
    398398               idta(:,:) = jpkm1 
    399399               DO jk = 1, jpkm1 
    400                   WHERE( gdept_0(jk) < zdta(:,:) .AND. zdta(:,:) <= gdept_0(jk+1) )   idta(:,:) = jk 
     400                  WHERE( gdept_1d(jk) < zdta(:,:) .AND. zdta(:,:) <= gdept_1d(jk+1) )   idta(:,:) = jk 
    401401               END DO 
    402402            ENDIF 
     
    511511      IF ( .not. ln_sco ) THEN                                !==  set a minimum depth  ==! 
    512512         IF( rn_hmin < 0._wp ) THEN    ;   ik = - INT( rn_hmin )                                      ! from a nb of level 
    513          ELSE                          ;   ik = MINLOC( gdepw_0, mask = gdepw_0 > rn_hmin, dim = 1 )  ! from a depth 
     513         ELSE                          ;   ik = MINLOC( gdepw_1d, mask = gdepw_1d > rn_hmin, dim = 1 )  ! from a depth 
    514514         ENDIF 
    515          zhmin = gdepw_0(ik+1)                                                         ! minimum depth = ik+1 w-levels  
     515         zhmin = gdepw_1d(ik+1)                                                         ! minimum depth = ik+1 w-levels  
    516516         WHERE( bathy(:,:) <= 0._wp )   ;   bathy(:,:) = 0._wp                         ! min=0     over the lands 
    517517         ELSE WHERE                     ;   bathy(:,:) = MAX(  zhmin , bathy(:,:)  )   ! min=zhmin over the oceans 
     
    780780      ! 
    781781      DO jk = 1, jpk 
    782             gdept(:,:,jk) = gdept_0(jk) 
    783             gdepw(:,:,jk) = gdepw_0(jk) 
    784             gdep3w(:,:,jk) = gdepw_0(jk) 
    785             e3t (:,:,jk) = e3t_0(jk) 
    786             e3u (:,:,jk) = e3t_0(jk) 
    787             e3v (:,:,jk) = e3t_0(jk) 
    788             e3f (:,:,jk) = e3t_0(jk) 
    789             e3w (:,:,jk) = e3w_0(jk) 
    790             e3uw(:,:,jk) = e3w_0(jk) 
    791             e3vw(:,:,jk) = e3w_0(jk) 
     782            gdept(:,:,jk) = gdept_1d(jk) 
     783            gdepw(:,:,jk) = gdepw_1d(jk) 
     784            gdep3w(:,:,jk) = gdepw_1d(jk) 
     785            e3t (:,:,jk) = e3t_1d(jk) 
     786            e3u (:,:,jk) = e3t_1d(jk) 
     787            e3v (:,:,jk) = e3t_1d(jk) 
     788            e3f (:,:,jk) = e3t_1d(jk) 
     789            e3w (:,:,jk) = e3w_1d(jk) 
     790            e3uw(:,:,jk) = e3w_1d(jk) 
     791            e3vw(:,:,jk) = e3w_1d(jk) 
    792792      END DO 
    793793      ! 
     
    837837      !!      schemes. 
    838838      !! 
    839       !!         c a u t i o n : gdept_0, gdepw_0 and e3._0 are positives 
     839      !!         c a u t i o n : gdept_1d, gdepw_1d and e3._1d are positives 
    840840      !!         - - - - - - -   gdept, gdepw and e3. are positives 
    841841      !!       
     
    874874      ! bathymetry in level (from bathy_meter) 
    875875      ! =================== 
    876       zmax = gdepw_0(jpk) + e3t_0(jpk)          ! maximum depth (i.e. the last ocean level thickness <= 2*e3t_0(jpkm1) ) 
     876      zmax = gdepw_1d(jpk) + e3t_1d(jpk)        ! maximum depth (i.e. the last ocean level thickness <= 2*e3t_1d(jpkm1) ) 
    877877      bathy(:,:) = MIN( zmax ,  bathy(:,:) )    ! bounded value of bathy (min already set at the end of zgr_bat) 
    878878      WHERE( bathy(:,:) == 0._wp )   ;   mbathy(:,:) = 0       ! land  : set mbathy to 0 
     
    882882      ! Compute mbathy for ocean points (i.e. the number of ocean levels) 
    883883      ! find the number of ocean levels such that the last level thickness 
    884       ! is larger than the minimum of e3zps_min and e3zps_rat * e3t_0 (where 
    885       ! e3t_0 is the reference level thickness 
     884      ! is larger than the minimum of e3zps_min and e3zps_rat * e3t_1d (where 
     885      ! e3t_1d is the reference level thickness 
    886886      DO jk = jpkm1, 1, -1 
    887          zdepth = gdepw_0(jk) + MIN( e3zps_min, e3t_0(jk)*e3zps_rat ) 
     887         zdepth = gdepw_1d(jk) + MIN( e3zps_min, e3t_1d(jk)*e3zps_rat ) 
    888888         WHERE( 0._wp < bathy(:,:) .AND. bathy(:,:) <= zdepth )   mbathy(:,:) = jk-1 
    889889      END DO 
     
    891891      ! Scale factors and depth at T- and W-points 
    892892      DO jk = 1, jpk                        ! intitialization to the reference z-coordinate 
    893          gdept(:,:,jk) = gdept_0(jk) 
    894          gdepw(:,:,jk) = gdepw_0(jk) 
    895          e3t  (:,:,jk) = e3t_0  (jk) 
    896          e3w  (:,:,jk) = e3w_0  (jk) 
     893         gdept(:,:,jk) = gdept_1d(jk) 
     894         gdepw(:,:,jk) = gdepw_1d(jk) 
     895         e3t  (:,:,jk) = e3t_1d  (jk) 
     896         e3w  (:,:,jk) = e3w_1d  (jk) 
    897897      END DO 
    898898      !  
     
    904904               IF( ik == jpkm1 ) THEN 
    905905                  zdepwp = bathy(ji,jj) 
    906                   ze3tp  = bathy(ji,jj) - gdepw_0(ik) 
    907                   ze3wp = 0.5_wp * e3w_0(ik) * ( 1._wp + ( ze3tp/e3t_0(ik) ) ) 
     906                  ze3tp  = bathy(ji,jj) - gdepw_1d(ik) 
     907                  ze3wp = 0.5_wp * e3w_1d(ik) * ( 1._wp + ( ze3tp/e3t_1d(ik) ) ) 
    908908                  e3t(ji,jj,ik  ) = ze3tp 
    909909                  e3t(ji,jj,ik+1) = ze3tp 
     
    911911                  e3w(ji,jj,ik+1) = ze3tp 
    912912                  gdepw(ji,jj,ik+1) = zdepwp 
    913                   gdept(ji,jj,ik  ) = gdept_0(ik-1) + ze3wp 
     913                  gdept(ji,jj,ik  ) = gdept_1d(ik-1) + ze3wp 
    914914                  gdept(ji,jj,ik+1) = gdept(ji,jj,ik) + ze3tp 
    915915                  ! 
    916916               ELSE                         ! standard case 
    917                   IF( bathy(ji,jj) <= gdepw_0(ik+1) ) THEN   ;   gdepw(ji,jj,ik+1) = bathy(ji,jj) 
    918                   ELSE                                       ;   gdepw(ji,jj,ik+1) = gdepw_0(ik+1) 
     917                  IF( bathy(ji,jj) <= gdepw_1d(ik+1) ) THEN  ;   gdepw(ji,jj,ik+1) = bathy(ji,jj) 
     918                  ELSE                                       ;   gdepw(ji,jj,ik+1) = gdepw_1d(ik+1) 
    919919                  ENDIF 
    920 !gm Bug?  check the gdepw_0 
     920!gm Bug?  check the gdepw_1d 
    921921                  !       ... on ik 
    922                   gdept(ji,jj,ik) = gdepw_0(ik) + ( gdepw  (ji,jj,ik+1) - gdepw_0(ik) )   & 
    923                      &                          * ((gdept_0(      ik  ) - gdepw_0(ik) )   & 
    924                      &                          / ( gdepw_0(      ik+1) - gdepw_0(ik) )) 
    925                   e3t  (ji,jj,ik) = e3t_0  (ik) * ( gdepw  (ji,jj,ik+1) - gdepw_0(ik) )   &  
    926                      &                          / ( gdepw_0(      ik+1) - gdepw_0(ik) )  
    927                   e3w  (ji,jj,ik) = 0.5_wp * ( gdepw(ji,jj,ik+1) + gdepw_0(ik+1) - 2._wp * gdepw_0(ik) )   & 
    928                      &                     * ( e3w_0(ik) / ( gdepw_0(ik+1) - gdepw_0(ik) ) ) 
     922                  gdept(ji,jj,ik) = gdepw_1d(ik) + ( gdepw  (ji,jj,ik+1) - gdepw_1d(ik) )   & 
     923                     &                           * ((gdept_1d(     ik  ) - gdepw_1d(ik) )   & 
     924                     &                           / ( gdepw_1d(     ik+1) - gdepw_1d(ik) )) 
     925                  e3t  (ji,jj,ik) = e3t_1d  (ik) * ( gdepw  (ji,jj,ik+1) - gdepw_1d(ik)  )  &  
     926                     &                           / ( gdepw_1d(      ik+1) - gdepw_1d(ik) )  
     927                  e3w  (ji,jj,ik) = 0.5_wp * ( gdepw(ji,jj,ik+1) + gdepw_1d(ik+1) - 2._wp * gdepw_1d(ik) )   & 
     928                     &                     * ( e3w_1d(ik) / ( gdepw_1d(ik+1) - gdepw_1d(ik) ) ) 
    929929                  !       ... on ik+1 
    930930                  e3w  (ji,jj,ik+1) = e3t  (ji,jj,ik) 
     
    958958      ! Scale factors and depth at U-, V-, UW and VW-points 
    959959      DO jk = 1, jpk                        ! initialisation to z-scale factors 
    960          e3u (:,:,jk) = e3t_0(jk) 
    961          e3v (:,:,jk) = e3t_0(jk) 
    962          e3uw(:,:,jk) = e3w_0(jk) 
    963          e3vw(:,:,jk) = e3w_0(jk) 
     960         e3u (:,:,jk) = e3t_1d(jk) 
     961         e3v (:,:,jk) = e3t_1d(jk) 
     962         e3uw(:,:,jk) = e3w_1d(jk) 
     963         e3vw(:,:,jk) = e3w_1d(jk) 
    964964      END DO 
    965965      DO jk = 1,jpk                         ! Computed as the minimum of neighbooring scale factors 
     
    977977      ! 
    978978      DO jk = 1, jpk                        ! set to z-scale factor if zero (i.e. along closed boundaries) 
    979          WHERE( e3u (:,:,jk) == 0._wp )   e3u (:,:,jk) = e3t_0(jk) 
    980          WHERE( e3v (:,:,jk) == 0._wp )   e3v (:,:,jk) = e3t_0(jk) 
    981          WHERE( e3uw(:,:,jk) == 0._wp )   e3uw(:,:,jk) = e3w_0(jk) 
    982          WHERE( e3vw(:,:,jk) == 0._wp )   e3vw(:,:,jk) = e3w_0(jk) 
     979         WHERE( e3u (:,:,jk) == 0._wp )   e3u (:,:,jk) = e3t_1d(jk) 
     980         WHERE( e3v (:,:,jk) == 0._wp )   e3v (:,:,jk) = e3t_1d(jk) 
     981         WHERE( e3uw(:,:,jk) == 0._wp )   e3uw(:,:,jk) = e3w_1d(jk) 
     982         WHERE( e3vw(:,:,jk) == 0._wp )   e3vw(:,:,jk) = e3w_1d(jk) 
    983983      END DO 
    984984       
    985985      ! Scale factor at F-point 
    986986      DO jk = 1, jpk                        ! initialisation to z-scale factors 
    987          e3f(:,:,jk) = e3t_0(jk) 
     987         e3f(:,:,jk) = e3t_1d(jk) 
    988988      END DO 
    989989      DO jk = 1, jpk                        ! Computed as the minimum of neighbooring V-scale factors 
     
    997997      ! 
    998998      DO jk = 1, jpk                        ! set to z-scale factor if zero (i.e. along closed boundaries) 
    999          WHERE( e3f(:,:,jk) == 0._wp )   e3f(:,:,jk) = e3t_0(jk) 
     999         WHERE( e3f(:,:,jk) == 0._wp )   e3f(:,:,jk) = e3t_1d(jk) 
    10001000      END DO 
    10011001!!gm  bug ? :  must be a do loop with mj0,mj1 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DOM/dtatsd.F90

    r3294 r3862  
    223223               DO jk = 1, jpk                        ! determines the intepolated T-S profiles at each (i,j) points 
    224224                  zl = fsdept_0(ji,jj,jk) 
    225                   IF(     zl < gdept_0(1  ) ) THEN          ! above the first level of data 
     225                  IF(     zl < gdept_1d(1  ) ) THEN          ! above the first level of data 
    226226                     ztp(jk) =  ptsd(ji,jj,1    ,jp_tem) 
    227227                     zsp(jk) =  ptsd(ji,jj,1    ,jp_sal) 
    228                   ELSEIF( zl > gdept_0(jpk) ) THEN          ! below the last level of data 
     228                  ELSEIF( zl > gdept_1d(jpk) ) THEN          ! below the last level of data 
    229229                     ztp(jk) =  ptsd(ji,jj,jpkm1,jp_tem) 
    230230                     zsp(jk) =  ptsd(ji,jj,jpkm1,jp_sal) 
    231231                  ELSE                                      ! inbetween : vertical interpolation between jkk & jkk+1 
    232232                     DO jkk = 1, jpkm1                                  ! when  gdept(jkk) < zl < gdept(jkk+1) 
    233                         IF( (zl-gdept_0(jkk)) * (zl-gdept_0(jkk+1)) <= 0._wp ) THEN 
    234                            zi = ( zl - gdept_0(jkk) ) / (gdept_0(jkk+1)-gdept_0(jkk)) 
     233                        IF( (zl-gdept_1d(jkk)) * (zl-gdept_1d(jkk+1)) <= 0._wp ) THEN 
     234                           zi = ( zl - gdept_1d(jkk) ) / (gdept_1d(jkk+1)-gdept_1d(jkk)) 
    235235                           ztp(jk) = ptsd(ji,jj,jkk,jp_tem) + ( ptsd(ji,jj,jkk+1,jp_tem) - ptsd(ji,jj,jkk,jp_tem) ) * zi  
    236236                           zsp(jk) = ptsd(ji,jj,jkk,jp_sal) + ( ptsd(ji,jj,jkk+1,jp_sal) - ptsd(ji,jj,jkk,jp_sal) ) * zi 
     
    260260                  ik = mbkt(ji,jj)  
    261261                  IF( ik > 1 ) THEN 
    262                      zl = ( gdept_0(ik) - fsdept_0(ji,jj,ik) ) / ( gdept_0(ik) - gdept_0(ik-1) ) 
     262                     zl = ( gdept_1d(ik) - fsdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) ) 
    263263                     ptsd(ji,jj,ik,jp_tem) = (1.-zl) * ptsd(ji,jj,ik,jp_tem) + zl * ptsd(ji,jj,ik-1,jp_tem) 
    264264                     ptsd(ji,jj,ik,jp_sal) = (1.-zl) * ptsd(ji,jj,ik,jp_sal) + zl * ptsd(ji,jj,ik-1,jp_sal) 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90

    r3764 r3862  
    218218            IF(lwp) WRITE(numout,*) '~~~~~~~~~~' 
    219219            ! 
    220             zh1 = gdept_0(  1  ) 
    221             zh2 = gdept_0(jpkm1) 
     220            zh1 = gdept_1d(  1  ) 
     221            zh2 = gdept_1d(jpkm1) 
    222222            ! 
    223223            zslope = ( zt1 - zt2 ) / ( zh1 - zh2 ) 
     
    399399         WRITE(numout,*) 
    400400         WRITE(numout,*) '              Initial temperature and salinity profiles:' 
    401          WRITE(numout, "(9x,' level   gdept_0   temperature   salinity   ')" ) 
    402          WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_0(jk), tsn(2,2,jk,jp_tem), tsn(2,2,jk,jp_sal), jk = 1, jpk ) 
     401         WRITE(numout, "(9x,' level   gdept_1d   temperature   salinity   ')" ) 
     402         WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_1d(jk), tsn(2,2,jk,jp_tem), tsn(2,2,jk,jp_sal), jk = 1, jpk ) 
    403403      ENDIF 
    404404 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/IOM/iom.F90

    r3771 r3862  
    118118 
    119119      ! vertical grid definition 
    120       CALL iom_set_axis_attr( "deptht", gdept_0 ) 
    121       CALL iom_set_axis_attr( "depthu", gdept_0 ) 
    122       CALL iom_set_axis_attr( "depthv", gdept_0 ) 
    123       CALL iom_set_axis_attr( "depthw", gdepw_0 ) 
     120      CALL iom_set_axis_attr( "deptht", gdept_1d ) 
     121      CALL iom_set_axis_attr( "depthu", gdept_1d ) 
     122      CALL iom_set_axis_attr( "depthv", gdept_1d ) 
     123      CALL iom_set_axis_attr( "depthw", gdepw_1d ) 
    124124# if defined key_floats 
    125125      CALL iom_set_axis_attr( "nfloat", (ji, ji=1,nfloat) ) 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/IOM/iom_ioipsl.F90

    r2715 r3862  
    406406               CALL flioputv( ioipslid, 'nav_lon'     , glamt(ix1:ix2, iy1:iy2) ) 
    407407               CALL flioputv( ioipslid, 'nav_lat'     , gphit(ix1:ix2, iy1:iy2) ) 
    408                CALL flioputv( ioipslid, 'nav_lev'     , gdept_0 ) 
     408               CALL flioputv( ioipslid, 'nav_lev'     , gdept_1d ) 
    409409               ! +++ WRONG VALUE: to be improved but not really useful... 
    410410               CALL flioputv( ioipslid, 'time_counter', kt ) 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/IOM/iom_nf90.F90

    r2715 r3862  
    528528               CALL iom_nf90_check(NF90_PUT_VAR( if90id, idmy, gphit(ix1:ix2, iy1:iy2) ), clinfo) 
    529529               CALL iom_nf90_check(NF90_INQ_VARID( if90id, 'nav_lev'     , idmy ), clinfo) 
    530                CALL iom_nf90_check(NF90_PUT_VAR( if90id, idmy, gdept_0                 ), clinfo) 
     530               CALL iom_nf90_check(NF90_PUT_VAR( if90id, idmy, gdept_1d                ), clinfo) 
    531531               ! +++ WRONG VALUE: to be improved but not really useful... 
    532532               CALL iom_nf90_check(NF90_INQ_VARID( if90id, 'time_counter', idmy ), clinfo) 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/LDF/ldfdyn.F90

    r3634 r3862  
    184184      !!---------------------------------------------------------------------- 
    185185 
    186       zm00 = TANH( ( pdam - gdept_0(1    ) ) / pwam ) 
    187       zm01 = TANH( ( pdam - gdept_0(jpkm1) ) / pwam ) 
     186      zm00 = TANH( ( pdam - gdept_1d(1    ) ) / pwam ) 
     187      zm01 = TANH( ( pdam - gdept_1d(jpkm1) ) / pwam ) 
    188188      zmhs = zm00 / zm01 
    189189      zmhb = ( 1.e0 - pbot ) / ( 1.e0 - zmhs ) / zm01 
     
    225225      !!---------------------------------------------------------------------- 
    226226 
    227       zm00 = TANH( ( pdam - gdept_0(1    ) ) / pwam ) 
    228       zm01 = TANH( ( pdam - gdept_0(jpkm1) ) / pwam ) 
     227      zm00 = TANH( ( pdam - gdept_1d(1    ) ) / pwam ) 
     228      zm01 = TANH( ( pdam - gdept_1d(jpkm1) ) / pwam ) 
    229229      zmhs = zm00 / zm01 
    230230      zmhb = ( 1.e0 - pbot ) / ( 1.e0 - zmhs ) / zm01 
     
    267267      !!---------------------------------------------------------------------- 
    268268 
    269       zm00 = TANH( ( pdam - gdept_0(1    ) ) / pwam )    
    270       zm01 = TANH( ( pdam - gdept_0(jpkm1) ) / pwam ) 
     269      zm00 = TANH( ( pdam - gdept_1d(1    ) ) / pwam )    
     270      zm01 = TANH( ( pdam - gdept_1d(jpkm1) ) / pwam ) 
    271271      zmhs = zm00 / zm01 
    272272      zmhb = ( 1.e0 - pbot ) / ( 1.e0 - zmhs ) / zm01 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/LDF/ldfdyn_c3d.h90

    r3294 r3862  
    386386 
    387387      DO jk=1, jpk 
    388          zcoef(jk) = 1.0_wp + NINT(9.0_wp*(gdept_0(jk)-800.0_wp)/(3000.0_wp-800.0_wp)) 
     388         zcoef(jk) = 1.0_wp + NINT(9.0_wp*(gdept_1d(jk)-800.0_wp)/(3000.0_wp-800.0_wp)) 
    389389         zcoef(jk) = MIN(10.0_wp, MAX(1.0_wp, zcoef(jk))) 
    390390         IF(lwp) WRITE(numout,'(4x,i3,6x,f7.3)') jk,zcoef(jk) 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/OBS/diaobs.F90

    r3651 r3862  
    10171017      USE dom_oce, ONLY : &             ! Ocean space and time domain variables 
    10181018         & rdt,           &                        
    1019          & gdept_0,       &              
     1019         & gdept_1d,       &              
    10201020         & tmask, umask, vmask                             
    10211021      USE phycst, ONLY : &              ! Physical constants 
     
    10791079                  &              kstp, jpi, jpj, jpk, nit000, idaystp,   & 
    10801080                  &              tsn(:,:,:,jp_tem), tsn(:,:,:,jp_sal),   & 
    1081                   &              gdept_0, tmask, n1dint, n2dint,         & 
     1081                  &              gdept_1d, tmask, n1dint, n2dint,        & 
    10821082                  &              kdailyavtypes = endailyavtypes ) 
    10831083            ELSE 
     
    10851085                  &              kstp, jpi, jpj, jpk, nit000, idaystp,   & 
    10861086                  &              tsn(:,:,:,jp_tem), tsn(:,:,:,jp_sal),   & 
    1087                   &              gdept_0, tmask, n1dint, n2dint               ) 
     1087                  &              gdept_1d, tmask, n1dint, n2dint              ) 
    10881088            ENDIF 
    10891089         END DO 
     
    11291129           ! zonal component of velocity 
    11301130           CALL obs_vel_opt( veldatqc(jveloset), kstp, jpi, jpj, jpk, & 
    1131               &              nit000, idaystp, un, vn, gdept_0, umask, vmask, & 
     1131              &              nit000, idaystp, un, vn, gdept_1d, umask, vmask, & 
    11321132                             n1dint, n2dint, ld_velav(jveloset) ) 
    11331133         END DO 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/OBS/obs_prep.F90

    r2715 r3862  
    7575         & glamt,   & 
    7676         & gphit,   & 
    77          & gdept_0, & 
     77         & gdept_1d,& 
    7878         & tmask,   & 
    7979         & nproc 
     
    193193         &                 profdata%var(1)%vdep,                        & 
    194194         &                 glamt,                 gphit,                & 
    195          &                 gdept_0,               tmask,                & 
     195         &                 gdept_1d,              tmask,                & 
    196196         &                 profdata%nqc,          profdata%var(1)%nvqc, & 
    197197         &                 iosdtobs,              ilantobs,             & 
     
    213213         &                 profdata%var(2)%vdep,                        & 
    214214         &                 glamt,                 gphit,                & 
    215          &                 gdept_0,               tmask,                & 
     215         &                 gdept_1d,              tmask,                & 
    216216         &                 profdata%nqc,          profdata%var(2)%nvqc, & 
    217217         &                 iosdsobs,              ilansobs,             & 
     
    916916         & glamt, glamu, glamv,    & 
    917917         & gphit, gphiu, gphiv,    & 
    918          & gdept_0, & 
     918         & gdept_1d,            & 
    919919         & tmask, umask, vmask,  & 
    920920         & nproc 
     
    10321032         &                 profdata%var(1)%vdep,                        & 
    10331033         &                 glamu,                 gphiu,                & 
    1034          &                 gdept_0,               umask,                & 
     1034         &                 gdept_1d,              umask,                & 
    10351035         &                 profdata%nqc,          profdata%var(1)%nvqc, & 
    10361036         &                 iosduobs,              ilanuobs,             & 
     
    10521052         &                 profdata%var(2)%vdep,                        & 
    10531053         &                 glamv,                 gphiv,                & 
    1054          &                 gdept_0,               vmask,                & 
     1054         &                 gdept_1d,              vmask,                & 
    10551055         &                 profdata%nqc,          profdata%var(2)%nvqc, & 
    10561056         &                 iosdvobs,              ilanvobs,             & 
     
    17091709      !! * Modules used 
    17101710      USE dom_oce, ONLY : &       ! Geographical information 
    1711          & gdepw_0                         
     1711         & gdepw_1d                         
    17121712 
    17131713      !! * Arguments 
     
    18261826               &  .OR. ( pobsphi(jobs) >   90.         )       & 
    18271827               &  .OR. ( pobsdep(jobsp) < 0.0          )       & 
    1828                &  .OR. ( pobsdep(jobsp) > gdepw_0(kpk) ) ) THEN 
     1828               &  .OR. ( pobsdep(jobsp) > gdepw_1d(kpk)) ) THEN 
    18291829               kobsqc(jobsp) = kobsqc(jobsp) + 11 
    18301830               kosdobs = kosdobs + 1 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/OBS/obs_read_prof.F90

    r2715 r3862  
    793793      !----------------------------------------------------------------------- 
    794794      IF ( ldt3d ) THEN 
    795          CALL obs_level_search( jpk, gdept_0, & 
     795         CALL obs_level_search( jpk, gdept_1d, & 
    796796            & profdata%nvprot(1), profdata%var(1)%vdep, & 
    797797            & profdata%var(1)%mvk ) 
    798798      ENDIF 
    799799      IF ( lds3d ) THEN 
    800          CALL obs_level_search( jpk, gdept_0, & 
     800         CALL obs_level_search( jpk, gdept_1d, & 
    801801            & profdata%nvprot(2), profdata%var(2)%vdep, & 
    802802            & profdata%var(2)%mvk ) 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/OBS/obs_read_vel.F90

    r2715 r3862  
    614614      ! Model level search 
    615615      !----------------------------------------------------------------------- 
    616       CALL obs_level_search( jpk, gdept_0, & 
     616      CALL obs_level_search( jpk, gdept_1d,          & 
    617617         & profdata%nvprot(1), profdata%var(1)%vdep, & 
    618618         & profdata%var(1)%mvk ) 
    619       CALL obs_level_search( jpk, gdept_0, & 
     619      CALL obs_level_search( jpk, gdept_1d,          & 
    620620         & profdata%nvprot(2), profdata%var(2)%vdep, & 
    621621         & profdata%var(2)%mvk ) 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/SBC/sbcrnf.F90

    r3832 r3862  
    390390         IF( rn_hrnf > 0._wp ) THEN 
    391391            nkrnf = 2 
    392             DO WHILE( nkrnf /= jpkm1 .AND. gdepw_0(nkrnf+1) < rn_hrnf )   ;   nkrnf = nkrnf + 1   ;   END DO 
     392            DO WHILE( nkrnf /= jpkm1 .AND. gdepw_1d(nkrnf+1) < rn_hrnf )   ;   nkrnf = nkrnf + 1   ;   END DO 
    393393            IF( ln_sco )   CALL ctl_warn( 'sbc_rnf: number of levels over which Kz is increased is computed for zco...' ) 
    394394         ENDIF 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/TRA/tradmp.F90

    r3294 r3862  
    759759      clname = 'dist.coast' 
    760760      itime  = 0 
    761       CALL ymds2ju( 0     , 1      , 1     , 0._wp , zdate0 ) 
    762       CALL restini( 'NONE', jpi    , jpj   , glamt, gphit ,   & 
    763          &          jpk   , gdept_0, clname, itime, zdate0,   & 
     761      CALL ymds2ju( 0     , 1       , 1     , 0._wp , zdate0 ) 
     762      CALL restini( 'NONE', jpi     , jpj   , glamt, gphit ,   & 
     763         &          jpk   , gdept_1d, clname, itime, zdate0,   & 
    764764         &          rdt   , icot                         ) 
    765765      CALL restput( icot, 'Tcoast', jpi, jpj, jpk, 0, pdct ) 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/TRA/traqsr.F90

    r3680 r3862  
    407407            ENDIF 
    408408 
    409             IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksr, ' ref depth = ', gdepw_0(nksr+1), ' m' 
     409            IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksr, ' ref depth = ', gdepw_1d(nksr+1), ' m' 
    410410            ! 
    411411            IF( nn_chldta == 1 ) THEN           !* Chl data : set sf_chl structure 
     
    477477            IF(lwp) THEN 
    478478               WRITE(numout,*) 
    479             IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksr, ' ref depth = ', gdepw_0(nksr+1), ' m' 
     479            IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksr, ' ref depth = ', gdepw_1d(nksr+1), ' m' 
    480480            ENDIF 
    481481            ! 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfini.F90

    r3680 r3862  
    147147      ELSE                                  ! Background profile of avt (fit a theoretical/observational profile (Krauss 1990) 
    148148         avmb(:) = rn_avm0 
    149          avtb(:) = rn_avt0 + ( 3.e-4_wp - 2._wp * rn_avt0 ) * 1.e-4_wp * gdepw_0(:)   ! m2/s 
     149         avtb(:) = rn_avt0 + ( 3.e-4_wp - 2._wp * rn_avt0 ) * 1.e-4_wp * gdepw_1d(:)   ! m2/s 
    150150         IF(ln_sco .AND. lwp)   CALL ctl_warn( 'avtb profile not valid in sco' ) 
    151151      ENDIF 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/SAS_SRC/diawri.F90

    r3331 r3862  
    220220            &          nit000-1, zjulian, zdt, nh_T, nid_T, domain_id=nidom, snc4chunks=snc4set ) 
    221221         CALL histvert( nid_T, "deptht", "Vertical T levels",      &  ! Vertical grid: gdept 
    222             &           "m", ipk, gdept_0, nz_T, "down" ) 
     222            &           "m", ipk, gdept_1d, nz_T, "down" ) 
    223223         !                                                            ! Index of ocean points 
    224224         CALL wheneq( jpi*jpj    , tmask, 1, 1., ndex_hT, ndim_hT )      ! surface 
     
    232232            &          nit000-1, zjulian, zdt, nh_U, nid_U, domain_id=nidom, snc4chunks=snc4set ) 
    233233         CALL histvert( nid_U, "depthu", "Vertical U levels",      &  ! Vertical grid: gdept 
    234             &           "m", ipk, gdept_0, nz_U, "down" ) 
     234            &           "m", ipk, gdept_1d, nz_U, "down" ) 
    235235         !                                                            ! Index of ocean points 
    236236         CALL wheneq( jpi*jpj    , umask, 1, 1., ndex_hU, ndim_hU )      ! surface 
     
    244244            &          nit000-1, zjulian, zdt, nh_V, nid_V, domain_id=nidom, snc4chunks=snc4set ) 
    245245         CALL histvert( nid_V, "depthv", "Vertical V levels",      &  ! Vertical grid : gdept 
    246             &          "m", ipk, gdept_0, nz_V, "down" ) 
     246            &          "m", ipk, gdept_1d, nz_V, "down" ) 
    247247         !                                                            ! Index of ocean points 
    248248         CALL wheneq( jpi*jpj    , vmask, 1, 1., ndex_hV, ndim_hV )      ! surface 
     
    391391          1, jpi, 1, jpj, nit000-1, zjulian, zdt, nh_i, id_i, domain_id=nidom, snc4chunks=snc4set ) ! Horizontal grid : glamt and gphit 
    392392      CALL histvert( id_i, "deptht", "Vertical T levels",   &    ! Vertical grid : gdept 
    393           "m", jpk, gdept_0, nz_i, "down") 
     393          "m", jpk, gdept_1d, nz_i, "down") 
    394394 
    395395      ! Declare all the output fields as NetCDF variables 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/TOP_SRC/PISCES/P4Z/p4zopt.F90

    r3475 r3862  
    370370      nksrp = trc_oce_ext_lev( r_si2, 0.33e2 )     ! max level of light extinction (Blue Chl=0.01) 
    371371      ! 
    372       IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksrp, ' ref depth = ', gdepw_0(nksrp+1), ' m' 
     372      IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksrp, ' ref depth = ', gdepw_1d(nksrp+1), ' m' 
    373373      ! 
    374374                         etot (:,:,:) = 0._wp 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/TOP_SRC/PISCES/SED/sed.F90

    r3443 r3862  
    1919   USE dom_oce , ONLY :   glamt     =>   glamt          !: longitude of t-point (degre) 
    2020   USE dom_oce , ONLY :   gphit     =>   gphit          !: latitude  of t-point (degre) 
    21    USE dom_oce , ONLY :   e3t_0     =>   e3t_0          !: reference depth of t-points (m) 
     21   USE dom_oce , ONLY :   e3t_1d    =>   e3t_1d         !: reference depth of t-points (m) 
    2222   USE dom_oce , ONLY :   mbkt      =>   mbkt           !: vertical index of the bottom last T- ocean level 
    2323   USE dom_oce , ONLY :   tmask     =>   tmask          !: land/ocean mask at t-points 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/TOP_SRC/PISCES/SED/sedini.F90

    r3443 r3862  
    142142         DO ji = 1, jpi 
    143143            ikt = mbkt(ji,jj)  
    144             IF( tmask(ji,jj,ikt) == 1 ) epkbot(ji,jj) = e3t_0(ikt) 
     144            IF( tmask(ji,jj,ikt) == 1 ) epkbot(ji,jj) = e3t_1d(ikt) 
    145145         ENDDO 
    146146      ENDDO 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/TOP_SRC/trcdia.F90

    r3294 r3862  
    187187 
    188188         ! Vertical grid for tracer : gdept 
    189          CALL histvert( nit5, 'deptht', 'Vertical T levels', 'm', ipk, gdept_0, ndepit5) 
     189         CALL histvert( nit5, 'deptht', 'Vertical T levels', 'm', ipk, gdept_1d, ndepit5) 
    190190 
    191191         ! Index of ocean points in 3D and 2D (surface) 
     
    308308         ! Vertical grid for 2d and 3d arrays 
    309309 
    310          CALL histvert( nitd, 'deptht', 'Vertical T levels','m', ipk, gdept_0, ndepitd) 
     310         CALL histvert( nitd, 'deptht', 'Vertical T levels','m', ipk, gdept_1d, ndepitd) 
    311311 
    312312         ! Declare all the output fields as NETCDF variables 
     
    439439            &    iiter, zjulian, zdt, nhoritb, nitb , domain_id=nidom, snc4chunks=snc4set ) 
    440440         ! Vertical grid for biological trends 
    441          CALL histvert(nitb, 'deptht', 'Vertical T levels', 'm', ipk, gdept_0, ndepitb) 
     441         CALL histvert(nitb, 'deptht', 'Vertical T levels', 'm', ipk, gdept_1d, ndepitb) 
    442442 
    443443         ! Declare all the output fields as NETCDF variables 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/TOP_SRC/trcdta.F90

    r3827 r3862  
    195195                     DO jk = 1, jpk                        ! determines the intepolated T-S profiles at each (i,j) points 
    196196                        zl = fsdept_0(ji,jj,jk) 
    197                         IF(     zl < gdept_0(1  ) ) THEN          ! above the first level of data 
     197                        IF(     zl < gdept_1d(1  ) ) THEN         ! above the first level of data 
    198198                           ztp(jk) =  ptrc(ji,jj,1    ,jn) 
    199                         ELSEIF( zl > gdept_0(jpk) ) THEN          ! below the last level of data 
     199                        ELSEIF( zl > gdept_1d(jpk) ) THEN         ! below the last level of data 
    200200                           ztp(jk) =  ptrc(ji,jj,jpkm1,jn) 
    201201                        ELSE                                      ! inbetween : vertical interpolation between jkk & jkk+1 
    202202                           DO jkk = 1, jpkm1                                  ! when  gdept(jkk) < zl < gdept(jkk+1) 
    203                               IF( (zl-gdept_0(jkk)) * (zl-gdept_0(jkk+1)) <= 0._wp ) THEN 
    204                                  zi = ( zl - gdept_0(jkk) ) / (gdept_0(jkk+1)-gdept_0(jkk)) 
     203                              IF( (zl-gdept_1d(jkk)) * (zl-gdept_1d(jkk+1)) <= 0._wp ) THEN 
     204                                 zi = ( zl - gdept_1d(jkk) ) / (gdept_1d(jkk+1)-gdept_1d(jkk)) 
    205205                                 ztp(jk) = ptrc(ji,jj,jkk,jn) + ( ptrc(ji,jj,jkk+1,jn) - ptrc(ji,jj,jkk,jn) ) * zi  
    206206                              ENDIF 
     
    226226                        ik = mbkt(ji,jj)  
    227227                        IF( ik > 1 ) THEN 
    228                            zl = ( gdept_0(ik) - fsdept_0(ji,jj,ik) ) / ( gdept_0(ik) - gdept_0(ik-1) ) 
     228                           zl = ( gdept_1d(ik) - fsdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) ) 
    229229                           ptrc(ji,jj,ik,jn) = (1.-zl) * ptrc(ji,jj,ik,jn) + zl * ptrc(ji,jj,ik-1,jn) 
    230230                        ENDIF 
  • branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/TOOLS/COMPILE/Fcheck_archfile.sh

    r3294 r3862  
    6969   fi 
    7070else 
    71    myfile=$( find ${MAIN_DIR}/ARCH -name arch-${2}.fcm -print ) 
     71   myfile=$( find -L ${MAIN_DIR}/ARCH -name arch-${2}.fcm -print ) 
    7272   if [ ${#myfile} -gt 0 ]; then 
    7373      ln -sf  ${myfile} ${COMPIL_DIR}/$1 
     
    7777      echo "Try makenemo -h for help" 
    7878      echo "EXITING..." 
     79      pwd 
     80      echo "find ${MAIN_DIR}/ARCH -name arch-${2}.fcm -print" 
    7981      exit 1        
    8082   fi    
Note: See TracChangeset for help on using the changeset viewer.