Changeset 2443 for branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/DOM
- Timestamp:
- 2010-11-29T08:52:36+01:00 (14 years ago)
- Location:
- branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/DOM
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/DOM/dom_oce.F90
r2392 r2443 6 6 !!====================================================================== 7 7 !! History : 1.0 ! 2005-10 (A. Beckmann, G. Madec) reactivate s-coordinate 8 !!---------------------------------------------------------------------- 9 !! NEMO/OPA 3.3 , NEMO Consortium (2010) 10 !! $Id$ 11 !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) 8 !! 3.3 ! 2010-11 (G. Madec) add mbk. arrays associated to the deepest ocean level 12 9 !!---------------------------------------------------------------------- 13 10 USE par_oce ! ocean parameters … … 172 169 !! masks, bathymetry 173 170 !! --------------------------------------------------------------------- 174 INTEGER , PUBLIC, DIMENSION(jpi,jpj) :: mbathy !: number of ocean level (=0, 1, ... , jpk-1) 175 REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: bathy !: ocean depth (meters) 176 REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: tmask_i !: interior domain T-point mask 177 REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: bmask !: land/ocean mask of barotropic stream function 178 179 REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) :: tmask, umask, vmask, fmask !: land/ocean mask at T-, U-, V- and F-points 180 181 REAL(wp), PUBLIC, DIMENSION(jpiglo) :: tpol, fpol !: north fold mask (nperio= 3 or 4) 171 INTEGER , PUBLIC, DIMENSION(jpi,jpj) :: mbathy !: number of ocean level (=0, 1, ... , jpk-1) 172 INTEGER , PUBLIC, DIMENSION(jpi,jpj) :: mbkt !: vertical index of the bottom last T- ocean level 173 INTEGER , PUBLIC, DIMENSION(jpi,jpj) :: mbku, mbkv !: vertical index of the bottom last U- and W- ocean level 174 REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: bathy !: ocean depth (meters) 175 REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: tmask_i !: interior domain T-point mask 176 REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: bmask !: land/ocean mask of barotropic stream function 177 178 REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) :: tmask, umask, vmask, fmask !: land/ocean mask at T-, U-, V- and F-pts 179 180 REAL(wp), PUBLIC, DIMENSION(jpiglo) :: tpol, fpol !: north fold mask (jperio= 3 or 4) 182 181 183 182 #if defined key_noslip_accurate … … 233 232 END FUNCTION Agrif_CFixed 234 233 #endif 234 !!---------------------------------------------------------------------- 235 !! NEMO/OPA 3.3 , NEMO Consortium (2010) 236 !! $Id$ 237 !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) 235 238 !!====================================================================== 236 239 END MODULE dom_oce -
branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/DOM/domzgr.F90
r2436 r2443 14 14 !! 3.0 ! 2008-06 (G. Madec) insertion of domzgr_zps.h90 & conding style 15 15 !! 3.2 ! 2009-07 (R. Benshila) Suppression of rigid-lid option 16 !! 3.3 ! 2010-11 (G. Madec) add mbk. arrays associated to the deepest ocean level 16 17 !!---------------------------------------------------------------------- 17 18 … … 497 498 IF( .NOT.lk_c1d ) CALL zgr_bat_ctl ! Bathymetry check ! 498 499 ! ! =================== ! 500 ! 501 ! ! ========================= ! 502 CALL zgr_bot_level ! level of ocean bottom ! 503 ! ! ========================= ! 499 504 END SUBROUTINE zgr_bat 500 505 … … 682 687 ! 683 688 END SUBROUTINE zgr_bat_ctl 689 690 691 SUBROUTINE zgr_bot_level 692 !!---------------------------------------------------------------------- 693 !! *** ROUTINE zgr_bot_level *** 694 !! 695 !! ** Purpose : defines the vertical index of ocean bottom (mbk. arrays) 696 !! 697 !! ** Method : computes from mbathy with a minimum value of 1 over land 698 !! 699 !! ** Action : mbkt, mbku, mbkv : vertical indices of the deeptest 700 !! ocean level at t-, u- & v-points 701 !! (min value = 1 over land) 702 !!---------------------------------------------------------------------- 703 INTEGER :: ji, jj ! dummy loop indices 704 REAL(wp), DIMENSION(jpi,jpj) :: zmbk ! 2D workspace 705 !!---------------------------------------------------------------------- 706 ! 707 IF(lwp) WRITE(numout,*) 708 IF(lwp) WRITE(numout,*) ' zgr_bot_level : ocean bottom k-index of T-, U-, V- and W-levels ' 709 IF(lwp) WRITE(numout,*) ' ~~~~~~~~~~~~~' 710 ! 711 mbkt(:,:) = MAX( mbathy(:,:) , 1 ) ! bottom k-index of T-level (=1 over land) 712 ! ! bottom k-index of W-level = mbkt+1 713 DO jj = 1, jpjm1 ! bottom k-index of u- (v-) level 714 DO ji = 1, jpim1 715 mbku(ji,jj) = MIN( mbkt(ji+1,jj ) , mbkt(ji,jj) ) 716 mbkv(ji,jj) = MIN( mbkt(ji ,jj+1) , mbkt(ji,jj) ) 717 END DO 718 END DO 719 ! converte into REAL to use lbc_lnk ; impose a min value of 1 as a zero can be set in lbclnk 720 zmbk(:,:) = REAL( mbku(:,:), wp ) ; CALL lbc_lnk(zmbk,'U',1.) ; mbku (:,:) = MAX( INT( zmbk(:,:) ), 1 ) 721 zmbk(:,:) = REAL( mbkv(:,:), wp ) ; CALL lbc_lnk(zmbk,'V',1.) ; mbkv (:,:) = MAX( INT( zmbk(:,:) ), 1 ) 722 ! 723 END SUBROUTINE zgr_bot_level 684 724 685 725
Note: See TracChangeset
for help on using the changeset viewer.