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 6667 for branches/2016/dev_r6409_SIMPLIF_2_usrdef/NEMOGCM/NEMO/OPA_SRC/DOM/dommsk.F90 – NEMO

Ignore:
Timestamp:
2016-06-06T07:57:00+02:00 (8 years ago)
Author:
gm
Message:

#1692 - branch SIMPLIF_2_usrdef: reduced domain_cfg.nc file: GYRE OK using usrdef or reading file

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2016/dev_r6409_SIMPLIF_2_usrdef/NEMOGCM/NEMO/OPA_SRC/DOM/dommsk.F90

    r6140 r6667  
    99   !!             -   ! 1996-05  (G. Madec)  mask computed from tmask 
    1010   !!            8.0  ! 1997-02  (G. Madec)  mesh information put in domhgr.F 
    11    !!            8.1  ! 1997-07  (G. Madec)  modification of mbathy and fmask 
     11   !!            8.1  ! 1997-07  (G. Madec)  modification of kbat and fmask 
    1212   !!             -   ! 1998-05  (G. Roullet)  free surface 
    1313   !!            8.2  ! 2000-03  (G. Madec)  no slip accurate 
     
    5050CONTAINS 
    5151 
    52    SUBROUTINE dom_msk 
     52   SUBROUTINE dom_msk( k_top, k_bot ) 
    5353      !!--------------------------------------------------------------------- 
    5454      !!                 ***  ROUTINE dom_msk  *** 
     
    5757      !!      zontal velocity points (u & v), vorticity points (f) points. 
    5858      !! 
    59       !! ** Method  :   The ocean/land mask is computed from the basin bathy- 
    60       !!      metry in level (mbathy) which is defined or read in dommba. 
    61       !!      mbathy equals 0 over continental T-point  
    62       !!      and the number of ocean level over the ocean. 
    63       !! 
    64       !!      At a given position (ji,jj,jk) the ocean/land mask is given by: 
    65       !!      t-point : 0. IF mbathy( ji ,jj) =< 0 
    66       !!                1. IF mbathy( ji ,jj) >= jk 
    67       !!      u-point : 0. IF mbathy( ji ,jj)  or mbathy(ji+1, jj ) =< 0 
    68       !!                1. IF mbathy( ji ,jj) and mbathy(ji+1, jj ) >= jk. 
    69       !!      v-point : 0. IF mbathy( ji ,jj)  or mbathy( ji ,jj+1) =< 0 
    70       !!                1. IF mbathy( ji ,jj) and mbathy( ji ,jj+1) >= jk. 
    71       !!      f-point : 0. IF mbathy( ji ,jj)  or mbathy( ji ,jj+1) 
    72       !!                   or mbathy(ji+1,jj)  or mbathy(ji+1,jj+1) =< 0 
    73       !!                1. IF mbathy( ji ,jj) and mbathy( ji ,jj+1) 
    74       !!                  and mbathy(ji+1,jj) and mbathy(ji+1,jj+1) >= jk. 
    75       !!      tmask_i : interior ocean mask at t-point, i.e. excluding duplicated 
    76       !!                rows/lines due to cyclic or North Fold boundaries as well 
    77       !!                as MPP halos. 
    78       !! 
    79       !!        The lateral friction is set through the value of fmask along 
    80       !!      the coast and topography. This value is defined by rn_shlat, a 
    81       !!      namelist parameter: 
     59      !! ** Method  :   The ocean/land mask  at t-point is deduced from ko_top  
     60      !!      and ko_bot, the indices of the fist and last ocean t-levels which  
     61      !!      are either defined in usrdef_zgr or read in zgr_read. 
     62      !!                The velocity masks (umask, vmask, wmask, wumask, wvmask)  
     63      !!      are deduced from a product of the two neighboring tmask. 
     64      !!                The vorticity mask (fmask) is deduced from tmask taking 
     65      !!      into account the choice of lateral boundary condition (rn_shlat) : 
    8266      !!         rn_shlat = 0, free slip  (no shear along the coast) 
    8367      !!         rn_shlat = 2, no slip  (specified zero velocity at the coast) 
     
    8569      !!         2 < rn_shlat, strong slip        | in the lateral boundary layer 
    8670      !! 
    87       !!      N.B. If nperio not equal to 0, the land/ocean mask arrays 
    88       !!      are defined with the proper value at lateral domain boundaries. 
     71      !!      tmask_i : interior ocean mask at t-point, i.e. excluding duplicated 
     72      !!                rows/lines due to cyclic or North Fold boundaries as well 
     73      !!                as MPP halos. 
     74      !!      tmask_h : halo mask at t-point, i.e. excluding duplicated rows/lines 
     75      !!                due to cyclic or North Fold boundaries as well 
     76      !!                as MPP halos. 
    8977      !! 
    9078      !!      In case of open boundaries (lk_bdy=T): 
    91       !!        - tmask is set to 1 on the points to be computed bay the open 
     79      !!        - tmask is set to 1 on the points to be computed by the open 
    9280      !!          boundaries routines. 
    9381      !! 
    94       !! ** Action :   tmask    : land/ocean mask at t-point (=0. or 1.) 
    95       !!               umask    : land/ocean mask at u-point (=0. or 1.) 
    96       !!               vmask    : land/ocean mask at v-point (=0. or 1.) 
    97       !!               fmask    : land/ocean mask at f-point (=0. or 1.) 
    98       !!                          =rn_shlat along lateral boundaries 
    99       !!               tmask_i  : interior ocean mask 
     82      !! ** Action :   tmask, umask, vmask, wmask, wumask, wvmask : land/ocean mask  
     83      !!                         at t-, u-, v- w, wu-, and wv-points (=0. or 1.) 
     84      !!               fmask   : land/ocean mask at f-point (=0., or =1., or  
     85      !!                         =rn_shlat along lateral boundaries) 
     86      !!               tmask_i : interior ocean mask  
     87      !!               tmask_h : halo mask 
     88      !!               ssmask , ssumask, ssvmask, ssfmask : 2D ocean mask 
    10089      !!---------------------------------------------------------------------- 
     90      INTEGER, DIMENSION(:,:), INTENT(in) ::   k_top, k_bot   ! first and last ocean level 
     91      ! 
    10192      INTEGER  ::   ji, jj, jk               ! dummy loop indices 
    10293      INTEGER  ::   iif, iil, ii0, ii1, ii   ! local integers 
    10394      INTEGER  ::   ijf, ijl, ij0, ij1       !   -       - 
     95      INTEGER  ::   iktop, ikbot             !   -       - 
    10496      INTEGER  ::   ios 
    10597      INTEGER  ::   isrow                    ! index for ORCA1 starting row 
    106       INTEGER , POINTER, DIMENSION(:,:) ::  imsk 
    10798      REAL(wp), POINTER, DIMENSION(:,:) ::  zwf 
    10899      !! 
     
    111102      ! 
    112103      IF( nn_timing == 1 )  CALL timing_start('dom_msk') 
    113       ! 
    114       CALL wrk_alloc( jpi, jpj, imsk ) 
    115       CALL wrk_alloc( jpi, jpj, zwf  ) 
    116104      ! 
    117105      REWIND( numnam_ref )              ! Namelist namlbc in reference namelist : Lateral momentum boundary condition 
     
    142130      ENDIF 
    143131 
    144       ! 1. Ocean/land mask at t-point (computed from mbathy) 
    145       ! ----------------------------- 
    146       ! N.B. tmask has already the right boundary conditions since mbathy is ok 
    147       ! 
    148       tmask(:,:,:) = 0._wp 
    149       DO jk = 1, jpk 
    150          DO jj = 1, jpj 
    151             DO ji = 1, jpi 
    152                IF( REAL( mbathy(ji,jj) - jk, wp ) + 0.1_wp >= 0._wp )   tmask(ji,jj,jk) = 1._wp 
    153             END DO   
     132 
     133      !  Ocean/land mask at t-point  (computed from ko_top and ko_bot) 
     134      ! ---------------------------- 
     135      ! 
     136      DO jj = 1, jpj 
     137         DO ji = 1, jpi 
     138            iktop = k_top(ji,jj) 
     139            ikbot = k_bot(ji,jj) 
     140            tmask(ji,jj,      1:iktop-1) = 0._wp      ! mask the iceshelves 
     141            tmask(ji,jj,iktop  :ikbot  ) = 1._wp 
     142            tmask(ji,jj,ikbot+1:jpkglo ) = 0._wp      ! mask the ocean topography 
    154143         END DO   
    155144      END DO   
     145 
     146      ! 2D ocean mask (=1 if at least one level of the water column is ocean, =0 otherwise) 
     147      WHERE( k_bot(:,:) > 0 )   ;   ssmask(:,:) = 1._wp 
     148      ELSEWHERE                 ;   ssmask(:,:) = 0._wp 
     149      END WHERE 
    156150       
    157       ! (ISF) define barotropic mask and mask the ice shelf point 
    158       ssmask(:,:)=tmask(:,:,1) ! at this stage ice shelf is not masked 
    159151       
    160       DO jk = 1, jpk 
    161          DO jj = 1, jpj 
    162             DO ji = 1, jpi 
    163                IF( REAL( misfdep(ji,jj) - jk, wp ) - 0.1_wp >= 0._wp )   THEN 
    164                   tmask(ji,jj,jk) = 0._wp 
    165                END IF 
    166             END DO   
    167          END DO   
    168       END DO   
    169  
    170       ! Interior domain mask (used for global sum) 
     152      ! Interior domain mask  (used for global sum) 
    171153      ! -------------------- 
    172       tmask_i(:,:) = ssmask(:,:)            ! (ISH) tmask_i = 1 even on the ice shelf 
    173  
    174       tmask_h(:,:) = 1._wp                 ! 0 on the halo and 1 elsewhere 
    175       iif = jpreci                         ! ??? 
    176       iil = nlci - jpreci + 1 
    177       ijf = jprecj                         ! ??? 
    178       ijl = nlcj - jprecj + 1 
    179  
     154      ! 
     155      iif = jpreci   ;   iil = nlci - jpreci + 1 
     156      ijf = jprecj   ;   ijl = nlcj - jprecj + 1 
     157      ! 
     158      !                          ! halo mask : 0 on the halo and 1 elsewhere 
     159      tmask_h(:,:) = 1._wp                   
    180160      tmask_h( 1 :iif,   :   ) = 0._wp      ! first columns 
    181161      tmask_h(iil:jpi,   :   ) = 0._wp      ! last  columns (including mpp extra columns) 
    182162      tmask_h(   :   , 1 :ijf) = 0._wp      ! first rows 
    183163      tmask_h(   :   ,ijl:jpj) = 0._wp      ! last  rows (including mpp extra rows) 
    184  
    185       ! north fold mask 
    186       ! --------------- 
     164      ! 
     165      !                          ! north fold mask 
    187166      tpol(1:jpiglo) = 1._wp  
    188167      fpol(1:jpiglo) = 1._wp 
     
    190169         tpol(jpiglo/2+1:jpiglo) = 0._wp 
    191170         fpol(     1    :jpiglo) = 0._wp 
    192          IF( mjg(nlej) == jpjglo ) THEN                  ! only half of the nlcj-1 row 
     171         IF( mjg(nlej) == jpjglo ) THEN                  ! only half of the nlcj-1 row for tmask_h 
    193172            DO ji = iif+1, iil-1 
    194173               tmask_h(ji,nlej-1) = tmask_h(ji,nlej-1) * tpol(mig(ji)) 
     
    196175         ENDIF 
    197176      ENDIF 
    198       
    199       tmask_i(:,:) = tmask_i(:,:) * tmask_h(:,:) 
    200  
     177      ! 
    201178      IF( jperio == 5 .OR. jperio == 6 ) THEN      ! F-point pivot 
    202179         tpol(     1    :jpiglo) = 0._wp 
    203180         fpol(jpiglo/2+1:jpiglo) = 0._wp 
    204181      ENDIF 
    205  
    206       ! 2. Ocean/land mask at u-,  v-, and z-points (computed from tmask) 
    207       ! ------------------------------------------- 
     182      ! 
     183      !                          ! interior mask : 2D ocean mask x halo mask  
     184      tmask_i(:,:) = ssmask(:,:) * tmask_h(:,:) 
     185 
     186 
     187      !  Ocean/land mask at u-, v-, and z-points (computed from tmask) 
     188      ! ---------------------------------------- 
    208189      DO jk = 1, jpk 
    209190         DO jj = 1, jpjm1 
     
    221202      DO jj = 1, jpjm1 
    222203         DO ji = 1, fs_jpim1   ! vector loop 
     204!!gm  simpler : 
     205!            ssumask(ji,jj)  = MIN(  1._wp , SUM( umask(ji,jj,:) )  ) 
     206!            ssvmask(ji,jj)  = MIN(  1._wp , SUM( vmask(ji,jj,:) )  ) 
     207!!gm 
     208!!gm  faster : 
     209!         ssumask(ji,jj) = ssmask(ji,jj) * tmask(ji+1,jj  ) 
     210!         ssvmask(ji,jj) = ssmask(ji,jj) * tmask(ji  ,jj+1) 
     211!!gm 
    223212            ssumask(ji,jj)  = ssmask(ji,jj) * ssmask(ji+1,jj  )  * MIN(1._wp,SUM(umask(ji,jj,:))) 
    224213            ssvmask(ji,jj)  = ssmask(ji,jj) * ssmask(ji  ,jj+1)  * MIN(1._wp,SUM(vmask(ji,jj,:))) 
     214!!end 
    225215         END DO 
    226216         DO ji = 1, jpim1      ! NO vector opt. 
     217!!gm faster 
     218!            ssfmask(ji,jj) =  ssmask(ji,jj  ) * ssmask(ji+1,jj  )   & 
     219!               &            * ssmask(ji,jj+1) * ssmask(ji+1,jj+1) 
     220!!gm  
    227221            ssfmask(ji,jj) =  ssmask(ji,jj  ) * ssmask(ji+1,jj  )   & 
    228222               &            * ssmask(ji,jj+1) * ssmask(ji+1,jj+1) * MIN(1._wp,SUM(fmask(ji,jj,:))) 
     223!!gm 
    229224         END DO 
    230225      END DO 
    231226      CALL lbc_lnk( umask  , 'U', 1._wp )      ! Lateral boundary conditions 
    232227      CALL lbc_lnk( vmask  , 'V', 1._wp ) 
    233       CALL lbc_lnk( fmask  , 'F', 1._wp ) 
    234       CALL lbc_lnk( ssumask, 'U', 1._wp )      ! Lateral boundary conditions 
     228!      CALL lbc_lnk( fmask  , 'F', 1._wp )         ! applied after the specification of lateral b.c. 
     229      CALL lbc_lnk( ssumask, 'U', 1._wp ) 
    235230      CALL lbc_lnk( ssvmask, 'V', 1._wp ) 
    236231      CALL lbc_lnk( ssfmask, 'F', 1._wp ) 
    237232 
    238       ! 3. Ocean/land mask at wu-, wv- and w points  
     233 
     234      ! Ocean/land mask at wu-, wv- and w points  
    239235      !---------------------------------------------- 
    240236      wmask (:,:,1) = tmask(:,:,1)     ! surface 
     
    247243      END DO 
    248244 
     245 
    249246      ! Lateral boundary conditions on velocity (modify fmask) 
    250247      ! ---------------------------------------      
     248      CALL wrk_alloc( jpi,jpj,   zwf  ) 
     249      ! 
    251250      DO jk = 1, jpk 
    252251         zwf(:,:) = fmask(:,:,jk)          
     
    277276      END DO 
    278277      ! 
     278      CALL wrk_dealloc( jpi,jpj,   zwf  ) 
     279      ! 
    279280      IF( cp_cfg == "orca" .AND. jp_cfg == 2 ) THEN   ! ORCA_R2 configuration 
    280281         !                                                 ! Increased lateral friction near of some straits 
     
    349350      ! CAUTION : The fmask may be further modified in dyn_vor_init ( dynvor.F90 ) 
    350351      ! 
    351       CALL wrk_dealloc( jpi, jpj, imsk ) 
    352       CALL wrk_dealloc( jpi, jpj, zwf  ) 
    353       ! 
    354352      IF( nn_timing == 1 )  CALL timing_stop('dom_msk') 
    355353      ! 
Note: See TracChangeset for help on using the changeset viewer.