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

Changeset 15131


Ignore:
Timestamp:
2021-07-21T16:37:16+02:00 (3 years ago)
Author:
dbruciaferri
Message:

introduce mes_init subroutine

Location:
NEMO/branches/UKMO/tools_r4.0-HEAD_dev_MEs/DOMAINcfg
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/UKMO/tools_r4.0-HEAD_dev_MEs/DOMAINcfg/namelist_cfg

    r15129 r15131  
    4040&namzgr_mes    !   MEs-coordinate                      
    4141!----------------------------------------------------------------------- 
    42    ln_envl     =   .TRUE. , .TRUE. , .TRUE. , .TRUE., .TRUE.   ! (T/F) If the envelope is used 
     42   ln_envl     =   .TRUE. , .TRUE. , .TRUE. , .TRUE., .FALSE.  ! (T/F) If the envelope is used 
    4343   nn_strt     =     1    ,   1    ,    1   ,   1    ,   1     ! Stretch. funct.: Madec 1996 (0) or 
    4444                                                               ! Song & Haidvogel 1994 (1) or                     
  • NEMO/branches/UKMO/tools_r4.0-HEAD_dev_MEs/DOMAINcfg/src/lib_mpp.f90

    r9598 r15131  
    7676   PUBLIC   mpp_lnk_2d_9 , mpp_lnk_2d_multiple  
    7777   PUBLIC   mpp_lnk_sum_3d, mpp_lnk_sum_2d 
    78    PUBLIC   mppscatter, mppgather 
     78   PUBLIC   mppscatter, mppgather, mppbcast_a_real 
    7979   PUBLIC   mpp_ini_ice, mpp_ini_znl 
    8080   PUBLIC   mppsize 
     
    18471847 
    18481848 
     1849   SUBROUTINE mppbcast_a_real( kvals, kno, kroot ) 
     1850      !!---------------------------------------------------------------------- 
     1851      !!                  ***  routine mppbcast_a_real  *** 
     1852      !! 
     1853      !! ** Purpose : Send array kvals to all processors 
     1854      !! 
     1855      !! ** Method  : MPI broadcast 
     1856      !! 
     1857      !!----------------------------------------------------------------------- 
     1858      INTEGER                 , INTENT(in   ) :: kno     ! Number of elements in array 
     1859      INTEGER                 , INTENT(in   ) :: kroot   ! Processor to send data 
     1860      REAL(wp), DIMENSION(kno), INTENT(inout) :: kvals   ! Array to send on kroot, receive for non-kroot 
     1861      !! 
     1862      INTEGER                   ::   ierr    ! temporary integer 
     1863      INTEGER                   ::   localcomm 
     1864      !!----------------------------------------------------------------------- 
     1865      ! 
     1866      localcomm = mpi_comm_opa 
     1867      CALL mpi_bcast( kvals, kno, mpi_double_precision, kroot, localcomm, ierr ) 
     1868      ! 
     1869   END SUBROUTINE mppbcast_a_real 
     1870 
     1871 
    18491872   SUBROUTINE mppmax_a_int( ktab, kdim, kcom ) 
    18501873      !!---------------------------------------------------------------------- 
  • NEMO/branches/UKMO/tools_r4.0-HEAD_dev_MEs/DOMAINcfg/src/mes.F90

    r15129 r15131  
    3535   INTEGER            :: nn_strt(max_nn_env)     ! Array specifing the stretching function for each envelope: 
    3636                                                 ! Madec 1996 (0), Song and Haidvogel 1994 (1)  
    37    REAL(wp)           :: rn_ebot_max(max_nn_env) ! maximum depth of envelopes (= envelopes depths) (>0) (m) 
     37   REAL(wp)           :: rn_ebot_max(max_nn_env) ! depth of envelopes in one point (for checking 1D transform.) 
     38   REAL(wp)           :: max_dep_env(max_nn_env) ! global maximum depth of envelopes 
     39   REAL(wp)           :: min_dep_env(max_nn_env) ! global minimum depth of envelopes 
    3840   INTEGER            :: nn_slev(max_nn_env)     ! Array specifing number of levels of each enveloped vertical zone 
    3941   REAL(wp)           :: rn_e_hc(max_nn_env)     ! Array specifing critical depth for transition to stretched  
     
    161163      REAL(wp), POINTER, DIMENSION(:,:  ) :: env0, env1, env2, env3     ! for loops over envelopes 
    162164                                                                        ! for cubic splines 
    163       INTEGER                             :: gst_envl(max_nn_env)       ! Array to deal with a ghost last envelope  
    164  
    165165      !!---------------------------------------------------------------------- 
    166166      ! 
     
    171171      CALL wrk_alloc( jpi, jpj, env_up, env_dw, env0, env1, env2, env3 ) 
    172172      CALL wrk_alloc( jpk, z_gsigw1, z_gsigt1, z_esigw1, z_esigt1 ) 
    173       ! 
    174       ! Determining if there is a "ghost" envelope:  
    175       gst_envl(:) = 0 
    176       IF ( nn_slev(tot_env) == 0 ) gst_envl(tot_env) = 1 
    177173      ! 
    178174      ! Initializing to 0.0 some arrays: 
     
    227223 
    228224         IF ( isodd(je) ) THEN 
    229             IF ( gst_envl(je) == 0 ) THEN 
    230                DO jk = 1, num_s 
    231                   ind = s_1st+jk-1 
    232                   ! W-GRID 
    233                   z_gsigw1(ind) = -s_coord(sigma(1,jk,num_s), coefa, & 
    234                                            coefb, alpha, num_s, nn_strt(je)) 
    235                   ! T-GRID 
    236                   z_gsigt1(ind) = -s_coord(sigma(0,jk,num_s), coefa, & 
    237                                            coefb, alpha, num_s, nn_strt(je)) 
    238                   zcoefw = z_gsigw1(ind) 
    239                   zcoeft = z_gsigt1(ind) 
    240                                
    241                   IF (nn_strt(je) /= 2) THEN 
    242                      gdept_1d(ind) = z_mes(dep_up, dep_dw, zcoeft, & 
    243                                            -sigma(0,jk,num_s), env_hc) 
    244                      gdepw_1d(ind) = z_mes(dep_up, dep_dw, zcoefw, & 
    245                                            -sigma(1,jk,num_s), env_hc) 
    246                   ELSE 
    247                      gdept_1d(ind) = z_mes(dep_up, dep_dw, zcoeft, & 
    248                                            -sigma(0,jk,num_s), 0._wp) 
    249                      gdepw_1d(ind) = z_mes(dep_up, dep_dw, zcoefw, & 
    250                                            -sigma(1,jk,num_s), 0._wp) 
    251                   END IF 
    252                ENDDO 
    253             ENDIF 
    254  
     225            DO jk = 1, num_s 
     226               ind = s_1st+jk-1 
     227               ! W-GRID 
     228               z_gsigw1(ind) = -s_coord(sigma(1,jk,num_s), coefa, & 
     229                                        coefb, alpha, num_s, nn_strt(je)) 
     230               ! T-GRID 
     231               z_gsigt1(ind) = -s_coord(sigma(0,jk,num_s), coefa, & 
     232                                        coefb, alpha, num_s, nn_strt(je)) 
     233               zcoefw = z_gsigw1(ind) 
     234               zcoeft = z_gsigt1(ind) 
     235               !                   
     236               IF (nn_strt(je) /= 2) THEN 
     237                   gdept_1d(ind) = z_mes(dep_up, dep_dw, zcoeft, & 
     238                                        -sigma(0,jk,num_s), env_hc) 
     239                   gdepw_1d(ind) = z_mes(dep_up, dep_dw, zcoefw, & 
     240                                        -sigma(1,jk,num_s), env_hc) 
     241               ELSE 
     242                   gdept_1d(ind) = z_mes(dep_up, dep_dw, zcoeft, & 
     243                                        -sigma(0,jk,num_s), 0._wp) 
     244                   gdepw_1d(ind) = z_mes(dep_up, dep_dw, zcoefw, & 
     245                                        -sigma(1,jk,num_s), 0._wp) 
     246               END IF 
     247            ENDDO 
     248            ! 
    255249         ELSE 
    256  
    257250            ! Ensure that for splines we do not use SF12 
    258251            !IF (nn_strt(je) > 1) nn_strt(je) = 1 
    259  
     252            ! 
    260253            ! Interval's endpoints TAU are 0 and 1. 
    261254            tau(1) = 0._wp 
     
    269262            dep1   = dep_up 
    270263            dep2   = dep_dw 
    271             dep3  = rn_ebot_max(je+1) 
    272  
     264            IF ( je < tot_env ) dep3 = rn_ebot_max(je+1) 
     265            ! 
    273266            ! Computing derivative analytically 
    274267            ! 
     
    297290            ! 
    298291            ! 2. Derivative for deeper envelope 
    299             IF ( gst_envl(je+1) == 0 ) THEN 
     292            IF ( je < tot_env ) THEN 
    300293               ibcend = 1 ! Bound. cond for the 1st derivative           
    301294               IF (nn_strt(je+1) == 2) THEN 
     
    405398      jk = 1 
    406399      e3w_1d(jk) = 2.0_wp * (gdept_1d(1) - gdepw_1d(1)) 
    407  
    408400  
    409401      ! Control print 
     
    463455 
    464456         IF ( isodd(je) ) THEN 
    465  
    466             IF ( gst_envl(je) == 0 ) THEN 
    467                DO jk = 1, num_s 
    468                   ind = s_1st+jk-1 
    469                   DO jj = 1,jpj 
    470                      DO ji = 1,jpi 
    471                         IF (env_dw(ji,jj) < env_hc) THEN    
    472                            ! shallow water, uniform sigma 
    473                            zcoefw = -sigma(1, jk, num_s)  
    474                            zcoeft = -sigma(0, jk, num_s) 
    475                            !IF (nn_strt(je) == 2) THEN 
    476                            !   ! shallow water, z coordinates 
    477                            !   ! SF12 in AMM* configurations uses  
    478                            !   ! this (ln_sigcrit=.false.) 
    479                            !   zcoefw = zcoefw*(env_hc/(env_dw(ji,jj)-env_up(ji,jj))) 
    480                            !   zcoeft = zcoeft*(env_hc/(env_dw(ji,jj)-env_up(ji,jj))) 
    481                            !ENDIF 
    482                         ELSE                               ! deep water, stretched s 
    483                            IF (nn_strt(je) == 2) THEN 
    484                               coefa = rn_e_th(je) /  (env_dw(ji,jj)-env_up(ji,jj)) 
    485                               coefb = rn_e_bb(je) + ((env_dw(ji,jj)-env_up(ji,jj))*rn_e_ba(je)) 
    486                               coefb = 1.0_wp-(coefb/(env_dw(ji,jj)-env_up(ji,jj))) 
    487                               ! W-GRID 
    488                               zcoefw = -s_coord( sigma(1,jk,num_s), coefa, & 
    489                                                  coefb, alpha, num_s, nn_strt(je) ) 
    490                               ! T-GRID 
    491                               zcoeft = -s_coord( sigma(0,jk,num_s), coefa, & 
    492                                                  coefb, alpha, num_s, nn_strt(je) ) 
    493                            ELSE 
    494                               zcoefw = z_gsigw1(ind) 
    495                               zcoeft = z_gsigt1(ind) 
    496                            END IF 
    497                         ENDIF 
    498  
    499                         IF (nn_strt(je) /= 2) THEN 
    500                            gdept_0(ji,jj,ind) = z_mes(env_up(ji,jj), env_dw(ji,jj), zcoeft, & 
    501                                                       -sigma(0, jk, num_s), env_hc) 
     457            ! 
     458            DO jk = 1, num_s 
     459               ind = s_1st+jk-1 
     460               DO jj = 1,jpj 
     461                  DO ji = 1,jpi 
     462                     IF (env_dw(ji,jj) < env_hc) THEN    
     463                        ! shallow water, uniform sigma 
     464                        zcoefw = -sigma(1, jk, num_s)  
     465                        zcoeft = -sigma(0, jk, num_s) 
     466                        !IF (nn_strt(je) == 2) THEN 
     467                        !   ! shallow water, z coordinates 
     468                        !   ! SF12 in AMM* configurations uses  
     469                        !   ! this (ln_sigcrit=.false.) 
     470                        !   zcoefw = zcoefw*(env_hc/(env_dw(ji,jj)-env_up(ji,jj))) 
     471                        !   zcoeft = zcoeft*(env_hc/(env_dw(ji,jj)-env_up(ji,jj))) 
     472                        !ENDIF 
     473                     ELSE                               ! deep water, stretched s 
     474                        IF (nn_strt(je) == 2) THEN 
     475                           coefa = rn_e_th(je) /  (env_dw(ji,jj)-env_up(ji,jj)) 
     476                           coefb = rn_e_bb(je) + ((env_dw(ji,jj)-env_up(ji,jj))*rn_e_ba(je)) 
     477                           coefb = 1.0_wp-(coefb/(env_dw(ji,jj)-env_up(ji,jj))) 
     478                           ! W-GRID 
     479                           zcoefw = -s_coord( sigma(1,jk,num_s), coefa, & 
     480                                              coefb, alpha, num_s, nn_strt(je) ) 
     481                           ! T-GRID 
     482                           zcoeft = -s_coord( sigma(0,jk,num_s), coefa, & 
     483                                              coefb, alpha, num_s, nn_strt(je) ) 
     484                        ELSE 
     485                           zcoefw = z_gsigw1(ind) 
     486                           zcoeft = z_gsigt1(ind) 
     487                        END IF 
     488                     ENDIF 
     489                     ! 
     490                     IF (nn_strt(je) /= 2) THEN 
     491                        gdept_0(ji,jj,ind) = z_mes(env_up(ji,jj), env_dw(ji,jj), zcoeft, & 
     492                                                   -sigma(0, jk, num_s), env_hc) 
    502493                                                 
    503                            gdepw_0(ji,jj,ind) = z_mes(env_up(ji,jj), env_dw(ji,jj), zcoefw, & 
    504                                                       -sigma(1, jk, num_s), env_hc) 
    505                         ELSE 
    506                            gdept_0(ji,jj,ind) = z_mes(env_up(ji,jj), env_dw(ji,jj), zcoeft, & 
    507                                                       -sigma(0, jk, num_s), 0._wp) 
    508  
    509                            gdepw_0(ji,jj,ind) = z_mes(env_up(ji,jj), env_dw(ji,jj), zcoefw, & 
    510                                                       -sigma(1, jk, num_s), 0._wp) 
    511                         END IF 
    512                      END DO 
     494                        gdepw_0(ji,jj,ind) = z_mes(env_up(ji,jj), env_dw(ji,jj), zcoefw, & 
     495                                                   -sigma(1, jk, num_s), env_hc) 
     496                     ELSE 
     497                        gdept_0(ji,jj,ind) = z_mes(env_up(ji,jj), env_dw(ji,jj), zcoeft, & 
     498                                                   -sigma(0, jk, num_s), 0._wp) 
     499  
     500                        gdepw_0(ji,jj,ind) = z_mes(env_up(ji,jj), env_dw(ji,jj), zcoefw, & 
     501                                                   -sigma(1, jk, num_s), 0._wp) 
     502                     END IF 
    513503                  END DO 
    514504               END DO 
    515             END IF 
    516  
     505            END DO 
     506            ! 
    517507         ELSE 
    518  
     508            ! 
    519509            ! Interval's endpoints TAU are 0 and 1. 
    520510            tau(1) = 0._wp 
    521511            tau(n) = 1._wp 
    522  
     512            ! 
    523513            IF ( je == 2 ) THEN 
    524514               env0 = scosrf 
     
    528518            env1   = env_up 
    529519            env2   = env_dw 
    530             env3  = envlt(:,:,je+1) 
    531  
     520            IF ( je < tot_env ) env3 = envlt(:,:,je+1) 
     521            !             
    532522            DO jj = 1,jpj 
    533523               DO ji = 1,jpi 
     
    564554                  ! 
    565555                  ! 2. Derivative for lower envelope 
    566                   IF ( gst_envl(je+1) == 0 ) THEN 
     556                  IF ( je < tot_env ) THEN 
    567557                     ibcend = 1 ! Bound. cond for the 1st derivative 
    568558                     IF (nn_strt(je+1) == 2) THEN 
     
    742732            s_1st  = s_1st + num_s - 1 
    743733            num_s  = nn_slev(je) + 1 
    744             max_env_up = MAXVAL(envlt(:,:,je-1)) 
    745             min_env_up = MINVAL(envlt(:,:,je-1)) 
    746             IF( lk_mpp ) CALL mpp_max( max_env_up ) 
    747             IF( lk_mpp ) CALL mpp_min( min_env_up ) 
     734            max_env_up = max_dep_env(je-1) 
     735            min_env_up = min_dep_env(je-1) 
    748736         ENDIF 
    749737 
    750          max_env_dw  = MAXVAL(envlt(:,:,je)) 
    751          min_env_dw  = MINVAL(envlt(:,:,je)) 
    752          IF( lk_mpp ) CALL mpp_max( max_env_dw ) 
    753          IF( lk_mpp ) CALL mpp_min( min_env_dw ) 
     738         max_env_dw  = max_dep_env(je) 
     739         min_env_dw  = min_dep_env(je) 
    754740 
    755741         IF ( max_env_up == min_env_up .AND. max_env_dw == min_env_dw ) THEN 
     
    891877   SUBROUTINE mes_init 
    892878      CHARACTER(lc)                :: env_name   ! name of the externally defined envelope 
    893       INTEGER                      :: ji, jj, je, iiemax, ijemax 
     879      INTEGER                      :: ji, jj, je, eii, eij 
    894880      INTEGER                      :: cor_env, ios, inum 
    895881      REAL(wp)                     :: rn_bot_min ! minimum depth of ocean bottom (>0) (m) 
    896882      REAL(wp)                     :: rn_bot_max ! maximum depth of ocean bottom (= ocean depth) (>0) (m) 
    897       REAL(wp), DIMENSION(jpi,jpj) :: pmsk       ! for loops over envelopes 
     883      REAL(wp)                     :: mindep 
     884      REAL(wp), DIMENSION(jpi,jpj) :: pmsk, ewrk ! working arrays 
    898885 
    899886      NAMELIST/namzgr_mes/rn_bot_min , rn_bot_max , ln_envl, & 
     
    906893      ! 
    907894      CALL wrk_alloc( jpi, jpj, max_nn_env , envlt ) 
    908       pmsk(:,:) = 1._wp 
     895      pmsk(:,:) = 0.0 
     896      ! 
     897      ! Initialising some arrays 
     898      rn_ebot_max(:) = 0.0 
     899      max_dep_env(:) = 0.0 
     900      min_dep_env(:) = 0.0 
     901      ! 
     902      IF(lwp) THEN                           ! control print 
     903         WRITE(numout,*) '' 
     904         WRITE(numout,*) 'mes_init: Setting a Multi-Envelope s-ccordinate system (Bruciaferri et al. 2018)' 
     905         WRITE(numout,*) '~~~~~~~~' 
     906      END IF 
    909907      ! 
    910908      ! Namelist namzgr_mes in reference namelist : envelopes and 
     
    931929      WRITE(ctlmes,*) 'num. of REQUESTED env. and num. of CORRECTLY defined env. are DIFFERENT' 
    932930      IF ( tot_env /= cor_env ) CALL ctl_stop( ctlmes ) 
     931      IF(lwp) THEN 
     932        WRITE(numout,*) '          Num. or requested envelopes: ', tot_env 
     933        WRITE(numout,*) '' 
     934      END IF 
    933935      ! 
    934936      ! 2) Checking consistency of user defined parameters 
     
    938940      ! 3) Reading Bathymetry and envelopes 
    939941      IF( ntopo == 1 ) THEN 
     942        IF(lwp) THEN 
     943           WRITE(numout,*) '          Reading bathymetry and envelopes ...' 
     944           WRITE(numout,*) '' 
     945        END IF 
     946 
    940947        CALL iom_open ( 'bathy_meter.nc', inum ) 
    941948        CALL iom_get  ( inum, jpdom_data, 'Bathymetry' , bathy, lrowattr=ln_use_jattr  ) 
     
    956963         IF (MAXVAL(envlt(:,:,je+1)) < MAXVAL(envlt(:,:,je))) CALL ctl_stop( ctlmes ) 
    957964      ENDDO 
    958       ! 5) Computing max depth of envelopes in the deepest 
    959       !    point of the domain for first check of monotonicity 
     965      ! 5) Computing max depth of envelopes in the shallowest 
     966      !    point of the deeper envelope for first check of monotonicity 
    960967      !    of transformation and computing 1D MEs-levels depths 
    961968      !    (used by diawri.F90) 
    962       CALL mpp_maxloc( envlt(:,:,tot_env), pmsk(:,:), rn_ebot_max(tot_env), iiemax, ijemax ) 
     969      WHERE ( bathy > 0.0 ) pmsk = 1.0 
     970      ewrk(:,:) = envlt(:,:,tot_env) 
     971      CALL mpp_minloc( ewrk(:,:), pmsk(:,:), mindep, eii, eij ) 
     972      IF ((mi0(eii)>1 .AND. mi0(eii)<jpi) .AND. (mj0(eij)>1 .AND. mj0(eij)<jpj)) THEN         
     973         inum = mpprank 
     974         DO je = 1, tot_env 
     975            rn_ebot_max(je) = envlt(mi0(eii),mj0(eij),je) 
     976         END DO 
     977      ELSE 
     978         inum = -1 
     979      END IF 
     980      IF( lk_mpp ) CALL mppsync 
     981      IF( lk_mpp ) CALL mpp_max(inum) 
     982      IF( lk_mpp ) CALL mppbcast_a_real(rn_ebot_max, max_nn_env, inum ) 
    963983      DO je = 1, tot_env 
    964          rn_ebot_max(je) = envlt(iiemax,ijemax,je) 
     984         max_dep_env(je) = MAXVAL(envlt(:,:,je)) 
     985         min_dep_env(je) = MINVAL(envlt(:,:,je)) 
     986         IF( lk_mpp ) CALL mpp_max( max_dep_env(je) ) 
     987         IF( lk_mpp ) CALL mpp_min( min_dep_env(je) ) 
    965988      END DO 
     989      IF( lk_mpp ) CALL mppsync 
    966990      ! 
    967991      ! 6) Set maximum and minimum ocean depth 
     
    975999      IF(lwp) THEN                           ! control print 
    9761000         WRITE(numout,*) 
    977          WRITE(numout,*) 'domzgr_mes : Multi Enveloped S-coordinate (Bruciaferri, Shapiro and Wobus 2017)' 
    978          WRITE(numout,*) '~~~~~~~~~~~' 
    979          WRITE(numout,*) '   Namelist namzgr_mes' 
     1001         WRITE(numout,*) '          GEOMETRICAL FEATURES OF THE REQUESTED MEs GRID:' 
     1002         WRITE(numout,*) '          -----------------------------------------------' 
     1003         WRITE(numout,*) '          Minimum depth of the ocean   rn_bot_min = ', rn_bot_min 
     1004         WRITE(numout,*) '          Maximum depth of the ocean   rn_bot_max = ', rn_bot_max 
    9801005         WRITE(numout,*) '' 
    981          WRITE(numout,*) '   Minimum depth of the ocean   rn_bot_min, ', rn_bot_min 
    982          WRITE(numout,*) '   Maximum depth of the ocean   rn_bot_max, ', rn_bot_max 
     1006         WRITE(numout,*) '          ENVELOPES:' 
     1007         WRITE(numout,*) '          ========= ' 
    9831008         WRITE(numout,*) '' 
    984          WRITE(numout,*) '-------------------------------------------------------------------------------' 
    985          DO je = 1, max_nn_env 
    986             WRITE(numout,*) 'SUBDOMAIN ', je,':' 
     1009         DO je = 1, tot_env 
     1010            WRITE(numout,*) '             * envelope ', je, ':   min depth = ', min_dep_env(je) 
     1011            WRITE(numout,*) '                               max depth = ', max_dep_env(je) 
     1012            WRITE(numout,*) '' 
     1013         END DO 
     1014         WRITE(numout,*) '' 
     1015         WRITE(numout,*) '          SUBDOMAINS:' 
     1016         WRITE(numout,*) '          ========== ' 
     1017         WRITE(numout,*) '' 
     1018         DO je = 1, tot_env 
     1019            WRITE(numout,*) '             * subdomain ', je,':' 
    9871020            IF ( je == 1) THEN 
    988                WRITE(numout,*) '   Envelope up  : envelope 0, free surface' 
    989                WRITE(numout,*) '   Envelope down: envelope ',je,',ln_envl(',je,') = ',ln_envl(je) 
     1021               WRITE(numout,*) '                  Shallower envelope: free surface' 
    9901022            ELSE 
    991                WRITE(numout,*) '   Envelope up  : envelope ',je-1,',ln_envl(',je,') = ',ln_envl(je-1) 
    992                WRITE(numout,*) '   Envelope down: envelope ',je  ,',ln_envl(',je,') = ',ln_envl(je) 
     1023               WRITE(numout,*) '                  Shallower envelope: envelope ',je-1 
    9931024            END IF 
    994             WRITE(numout,*) '   max dep of envlp down rn_ebot_max(',je,') =',rn_ebot_max(je) 
    995             WRITE(numout,*) '   num. of MEs-lev.          nn_slev(',je,') = ',nn_slev(je) 
     1025            WRITE(numout,*) '                  Deeper    envelope: envelope ',je 
     1026            WRITE(numout,*) '                  Number of MEs-levs: nn_slev(',je,') = ',nn_slev(je) 
    9961027            IF ( isodd(je) ) THEN 
    997                WRITE(numout,*) '   Stretched s-coordinates: ' 
     1028               WRITE(numout,*) '                  Stretched s-coordinates: ' 
    9981029            ELSE 
    999                WRITE(numout,*) '   Stretched CUBIC SPLINES: ' 
     1030               WRITE(numout,*) '                  Stretched CUBIC SPLINES: ' 
    10001031            END IF 
    1001             IF (nn_strt(je) == 0) WRITE(numout,*) '     M96  stretching function' 
    1002             IF (nn_strt(je) == 1) WRITE(numout,*) '     SH94 stretching function' 
    1003             IF (nn_strt(je) == 2) WRITE(numout,*) '     SF12 stretching function' 
    1004             WRITE(numout,*) '     critical depth        rn_e_hc(',je,') = ',rn_e_hc(je) 
    1005             WRITE(numout,*) '     surface stretc. coef. rn_e_th(',je,') = ',rn_e_th(je) 
     1032            IF (nn_strt(je) == 0) WRITE(numout,*) '                    M96  stretching function' 
     1033            IF (nn_strt(je) == 1) WRITE(numout,*) '                    SH94 stretching function' 
     1034            IF (nn_strt(je) == 2) WRITE(numout,*) '                    SF12 stretching function' 
     1035            WRITE(numout,*) '                    critical depth        rn_e_hc(',je,') = ',rn_e_hc(je) 
     1036            WRITE(numout,*) '                    surface stretc. coef. rn_e_th(',je,') = ',rn_e_th(je) 
    10061037            IF (nn_strt(je) == 2) THEN 
    1007                WRITE(numout,*) '     bottom  stretc. coef. rn_e_ba(',je,') = ',rn_e_ba(je) 
     1038               WRITE(numout,*) '                    bottom  stretc. coef. rn_e_ba(',je,') = ',rn_e_ba(je) 
    10081039            END IF 
    1009             WRITE(numout,*) '     bottom  stretc. coef. rn_e_bb(',je,') = ',rn_e_bb(je) 
     1040            WRITE(numout,*) '                    bottom  stretc. coef. rn_e_bb(',je,') = ',rn_e_bb(je) 
    10101041            IF (nn_strt(je) == 2) THEN 
    1011                WRITE(numout,*) '     bottom  stretc. coef. rn_e_al(',je,') = ',rn_e_al(je) 
     1042               WRITE(numout,*) '                 bottom  stretc. coef. rn_e_al(',je,') = ',rn_e_al(je) 
    10121043            END IF 
    1013             WRITE(numout,*) '-------------------------------------------------------------------------------' 
     1044            WRITE(numout,*) '          ------------------------------------------------------------------' 
    10141045         ENDDO 
    10151046      ENDIF 
Note: See TracChangeset for help on using the changeset viewer.