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.
usrdef_zgr.F90 in NEMO/branches/2020/dev_r13898_Tiling_Cleanup_MPI3/tests/OVERFLOW/MY_SRC – NEMO

source: NEMO/branches/2020/dev_r13898_Tiling_Cleanup_MPI3/tests/OVERFLOW/MY_SRC/usrdef_zgr.F90 @ 13906

Last change on this file since 13906 was 13906, checked in by mocavero, 4 years ago

Merge with dev_r13296_HPC-07_mocavero_mpi3

  • Property svn:keywords set to Id
File size: 13.1 KB
Line 
1MODULE usrdef_zgr
2   !!======================================================================
3   !!                     ***  MODULE usrdef_zgr  ***
4   !!
5   !!                       ===  OVERFLOW case  ===
6   !!
7   !! user defined :  vertical coordinate system of a user configuration
8   !!======================================================================
9   !! History :  4.0  ! 2016-08  (G. Madec, S. Flavoni)  Original code
10   !!----------------------------------------------------------------------
11
12   !!----------------------------------------------------------------------
13   !!   usr_def_zgr   : user defined vertical coordinate system (required)
14   !!       zgr_z1d   : reference 1D z-coordinate
15   !!---------------------------------------------------------------------
16   USE oce            ! ocean variables
17   USE dom_oce ,  ONLY: mi0, mi1   ! ocean space and time domain
18   USE dom_oce ,  ONLY: glamt      ! ocean space and time domain
19   USE usrdef_nam     ! User defined : namelist variables
20   !
21   USE in_out_manager ! I/O manager
22   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
23   USE lib_mpp        ! distributed memory computing library
24   USE timing         ! Timing
25
26   IMPLICIT NONE
27   PRIVATE
28
29   PUBLIC   usr_def_zgr   ! called by domzgr.F90
30
31   !! * Substitutions
32#  include "do_loop_substitute.h90"
33   !!----------------------------------------------------------------------
34   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
35   !! $Id$
36   !! Software governed by the CeCILL license (see ./LICENSE)
37   !!----------------------------------------------------------------------
38CONTAINS             
39
40   SUBROUTINE usr_def_zgr( ld_zco  , ld_zps  , ld_sco  , ld_isfcav,    &   ! type of vertical coordinate
41      &                    pdept_1d, pdepw_1d, pe3t_1d , pe3w_1d  ,    &   ! 1D reference vertical coordinate
42      &                    pdept , pdepw ,                             &   ! 3D t & w-points depth
43      &                    pe3t  , pe3u  , pe3v , pe3f ,               &   ! vertical scale factors
44      &                    pe3w  , pe3uw , pe3vw,                      &   !     -      -      -
45      &                    k_top  , k_bot    )                             ! top & bottom ocean level
46      !!---------------------------------------------------------------------
47      !!              ***  ROUTINE usr_def_zgr  ***
48      !!
49      !! ** Purpose :   User defined the vertical coordinates
50      !!
51      !!----------------------------------------------------------------------
52      LOGICAL                   , INTENT(in   ) ::   ld_zco, ld_zps, ld_sco      ! vertical coordinate flags ( read in namusr_def )
53      LOGICAL                   , INTENT(  out) ::   ld_isfcav                   ! under iceshelf cavity flag
54      REAL(wp), DIMENSION(:)    , INTENT(  out) ::   pdept_1d, pdepw_1d          ! 1D grid-point depth     [m]
55      REAL(wp), DIMENSION(:)    , INTENT(  out) ::   pe3t_1d , pe3w_1d           ! 1D grid-point depth     [m]
56      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pdept, pdepw                ! grid-point depth        [m]
57      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pe3t , pe3u , pe3v , pe3f   ! vertical scale factors  [m]
58      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pe3w , pe3uw, pe3vw         ! i-scale factors
59      INTEGER , DIMENSION(:,:)  , INTENT(  out) ::   k_top, k_bot                ! first & last ocean level
60      !
61      INTEGER  ::   ji, jj, jk        ! dummy indices
62      INTEGER  ::   ik                ! local integers
63      REAL(wp) ::   zfact, z1_jpkm1   ! local scalar
64      REAL(wp) ::   ze3min            ! local scalar
65      REAL(wp), DIMENSION(jpi,jpj) ::   zht, zhu, z2d   ! 2D workspace
66      !!----------------------------------------------------------------------
67      !
68      IF(lwp) WRITE(numout,*)
69      IF(lwp) WRITE(numout,*) 'usr_def_zgr : OVERFLOW configuration (z(ps)- or s-coordinate closed box ocean without cavities)'
70      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
71      !
72      !
73      ! type of vertical coordinate
74      ! ---------------------------
75      ! already set in usrdef_nam.F90 by reading the namusr_def namelist except for ISF
76      ld_isfcav = .FALSE.
77      !
78      !
79      ! Build the vertical coordinate system
80      ! ------------------------------------
81      !
82      !                       !==  UNmasked meter bathymetry  ==!
83      !
84      ! western continental shelf (500m deep) and eastern deep ocean (2000m deep)
85      ! (set through the jpk and jpi (see userdef_nam.F90))
86      ! with a hyperbolic tangent transition centered at 40km
87      ! NB: here glamt is in kilometers
88      !
89      zht(:,:) = + (  500. + 0.5 * 1500. * ( 1.0 + tanh( (glamt(:,:) - 40.) / 7. ) )  )
90      !
91      ! at u-point: averaging zht
92      DO ji = 1, jpim1
93         zhu(ji,:) = 0.5_wp * ( zht(ji,:) + zht(ji+1,:) )
94      END DO
95#if defined key_mpi3
96      CALL lbc_lnk_nc_multi( 'usrdef_zgr', zhu, 'U', 1. )     ! boundary condition: this mask the surrouding grid-points
97#else
98      CALL lbc_lnk( 'usrdef_zgr', zhu, 'U', 1. )     ! boundary condition: this mask the surrouding grid-points
99#endif
100      !                                ! ==>>>  set by hand non-zero value on first/last columns & rows
101      DO ji = mi0(1), mi1(1)              ! first row of global domain only
102         zhu(ji,2) = zht(ji,2)
103      END DO
104       DO ji = mi0(jpiglo), mi1(jpiglo)   ! last  row of global domain only
105         zhu(ji,2) = zht(ji,2)
106      END DO
107      zhu(:,1) = zhu(:,2)
108      zhu(:,3) = zhu(:,2)
109      !     
110      CALL zgr_z1d( pdept_1d, pdepw_1d, pe3t_1d , pe3w_1d )   ! Reference z-coordinate system
111      !
112      !
113      !                       !==  top masked level bathymetry  ==!  (all coordinates)
114      !
115      ! no ocean cavities : top ocean level is ONE, except over land
116      ! the ocean basin surrounded by land (1 grid-point) set through lbc_lnk call as jperio=0
117      z2d(:,:) = 1._wp                    ! surface ocean is the 1st level
118#if defined key_mpi3
119      CALL lbc_lnk_nc_multi( 'usrdef_zgr', z2d, 'T', 1. )        ! closed basin since jperio = 0 (see userdef_nam.F90)
120#else
121      CALL lbc_lnk( 'usrdef_zgr', z2d, 'T', 1. )        ! closed basin since jperio = 0 (see userdef_nam.F90)
122#endif
123      k_top(:,:) = NINT( z2d(:,:) )
124      !
125      !                             
126      !
127      IF ( ld_sco ) THEN      !==  s-coordinate  ==!   (terrain-following coordinate)
128         !
129         k_bot(:,:) = jpkm1 * k_top(:,:)  !* bottom ocean = jpk-1 (here use k_top as a land mask)
130         !
131         !                                !* terrain-following coordinate with e3.(k)=cst)
132         !                                !  OVERFLOW case : identical with j-index (T=V, U=F)
133         z1_jpkm1 = 1._wp / REAL( jpkm1 , wp)
134         DO jk = 1, jpk
135            pdept(:,:,jk) = zht(:,:) * z1_jpkm1 * ( REAL( jk   , wp ) - 0.5_wp )
136            pdepw(:,:,jk) = zht(:,:) * z1_jpkm1 * ( REAL( jk-1 , wp )          )
137            pe3t (:,:,jk) = zht(:,:) * z1_jpkm1
138            pe3u (:,:,jk) = zhu(:,:) * z1_jpkm1
139            pe3v (:,:,jk) = zht(:,:) * z1_jpkm1
140            pe3f (:,:,jk) = zhu(:,:) * z1_jpkm1
141            pe3w (:,:,jk) = zht(:,:) * z1_jpkm1
142            pe3uw(:,:,jk) = zhu(:,:) * z1_jpkm1
143            pe3vw(:,:,jk) = zht(:,:) * z1_jpkm1
144         END DO     
145      ENDIF
146      !
147      !
148      IF ( ld_zco ) THEN      !==  z-coordinate  ==!   (step-like topography)
149         !
150         !                                !* bottom ocean compute from the depth of grid-points
151         k_bot(:,:) = jpkm1 * k_top(:,:)     ! here use k_top as a land mask
152         DO jk = 1, jpkm1
153            WHERE( pdept_1d(jk) < zht(:,:) .AND. zht(:,:) <= pdept_1d(jk+1) )   k_bot(:,:) = jk * k_top(:,:)
154         END DO
155         !                                !* horizontally uniform coordinate (reference z-co everywhere)
156         DO jk = 1, jpk
157            pdept(:,:,jk) = pdept_1d(jk)
158            pdepw(:,:,jk) = pdepw_1d(jk)
159            pe3t (:,:,jk) = pe3t_1d (jk)
160            pe3u (:,:,jk) = pe3t_1d (jk)
161            pe3v (:,:,jk) = pe3t_1d (jk)
162            pe3f (:,:,jk) = pe3t_1d (jk)
163            pe3w (:,:,jk) = pe3w_1d (jk)
164            pe3uw(:,:,jk) = pe3w_1d (jk)
165            pe3vw(:,:,jk) = pe3w_1d (jk)
166         END DO
167      ENDIF
168      !
169      !
170      IF ( ld_zps ) THEN      !==  zps-coordinate  ==!   (partial bottom-steps)
171         !
172         ze3min = 0.1_wp * rn_dz
173         IF(lwp) WRITE(numout,*) '   minimum thickness of the partial cells = 10 % of e3 = ', ze3min
174         !
175         !
176         !                                !* bottom ocean compute from the depth of grid-points
177         k_bot(:,:) = jpkm1
178         DO jk = jpkm1, 1, -1
179            WHERE( zht(:,:) < pdepw_1d(jk) + ze3min )   k_bot(:,:) = jk-1
180         END DO
181         !
182         !                                !* vertical coordinate system
183         DO jk = 1, jpk                      ! initialization to the reference z-coordinate
184            pdept(:,:,jk) = pdept_1d(jk)
185            pdepw(:,:,jk) = pdepw_1d(jk)
186            pe3t (:,:,jk) = pe3t_1d (jk)
187            pe3u (:,:,jk) = pe3t_1d (jk)
188            pe3v (:,:,jk) = pe3t_1d (jk)
189            pe3f (:,:,jk) = pe3t_1d (jk)
190            pe3w (:,:,jk) = pe3w_1d (jk)
191            pe3uw(:,:,jk) = pe3w_1d (jk)
192            pe3vw(:,:,jk) = pe3w_1d (jk)
193         END DO
194         DO_2D( 1, 1, 1, 1 )
195            ik = k_bot(ji,jj)
196            pdepw(ji,jj,ik+1) = MIN( zht(ji,jj) , pdepw_1d(ik+1) )
197            pe3t (ji,jj,ik  ) = pdepw(ji,jj,ik+1) - pdepw(ji,jj,ik)
198            pe3t (ji,jj,ik+1) = pe3t (ji,jj,ik  ) 
199            !
200            pdept(ji,jj,ik  ) = pdepw(ji,jj,ik  ) + pe3t (ji,jj,ik  ) * 0.5_wp
201            pdept(ji,jj,ik+1) = pdepw(ji,jj,ik+1) + pe3t (ji,jj,ik+1) * 0.5_wp
202            pe3w (ji,jj,ik+1) = pdept(ji,jj,ik+1) - pdept(ji,jj,ik)              ! = pe3t (ji,jj,ik  )
203         END_2D         
204         !                                   ! bottom scale factors and depth at  U-, V-, UW and VW-points
205         !                                   ! usually Computed as the minimum of neighbooring scale factors
206         pe3u (:,:,:) = pe3t(:,:,:)          ! HERE OVERFLOW configuration :
207         pe3v (:,:,:) = pe3t(:,:,:)          !    e3 increases with i-index and identical with j-index
208         pe3f (:,:,:) = pe3t(:,:,:)          !    so e3 minimum of (i,i+1) points is (i) point
209         pe3uw(:,:,:) = pe3w(:,:,:)          !    in j-direction e3v=e3t and e3f=e3v
210         pe3vw(:,:,:) = pe3w(:,:,:)          !    ==>>  no need of lbc_lnk calls
211         !     
212      ENDIF
213      !
214   END SUBROUTINE usr_def_zgr
215
216
217   SUBROUTINE zgr_z1d( pdept_1d, pdepw_1d, pe3t_1d , pe3w_1d )   ! 1D reference vertical coordinate
218      !!----------------------------------------------------------------------
219      !!                   ***  ROUTINE zgr_z1d  ***
220      !!
221      !! ** Purpose :   set the depth of model levels and the resulting
222      !!      vertical scale factors.
223      !!
224      !! ** Method  :   1D z-coordinate system (use in all type of coordinate)
225      !!       The depth of model levels is set from dep(k), an analytical function:
226      !!                   w-level: depw_1d  = dep(k)
227      !!                   t-level: dept_1d  = dep(k+0.5)
228      !!       The scale factors are the discrete derivative of the depth:
229      !!                   e3w_1d(jk) = dk[ dept_1d ]
230      !!                   e3t_1d(jk) = dk[ depw_1d ]
231      !!
232      !!            ===    Here constant vertical resolution   ===
233      !!
234      !! ** Action  : - pdept_1d, pdepw_1d : depth of T- and W-point (m)
235      !!              - pe3t_1d , pe3w_1d  : scale factors at T- and W-levels (m)
236      !!----------------------------------------------------------------------
237      REAL(wp), DIMENSION(:), INTENT(out) ::   pdept_1d, pdepw_1d   ! 1D grid-point depth        [m]
238      REAL(wp), DIMENSION(:), INTENT(out) ::   pe3t_1d , pe3w_1d    ! 1D vertical scale factors  [m]
239      !
240      INTEGER  ::   jk       ! dummy loop indices
241      REAL(wp) ::   zt, zw   ! local scalar
242      !!----------------------------------------------------------------------
243      !
244      IF(lwp) THEN                         ! Parameter print
245         WRITE(numout,*)
246         WRITE(numout,*) '    zgr_z1d : Reference vertical z-coordinates: uniform dz = ', rn_dz
247         WRITE(numout,*) '    ~~~~~~~'
248      ENDIF
249      !
250      ! Reference z-coordinate (depth - scale factor at T- and W-points)   ! Madec & Imbard 1996 function
251      ! ----------------------
252      DO jk = 1, jpk
253         zw = REAL( jk , wp )
254         zt = REAL( jk , wp ) + 0.5_wp
255         pdepw_1d(jk) =    rn_dz *   REAL( jk-1 , wp )
256         pdept_1d(jk) =    rn_dz * ( REAL( jk-1 , wp ) + 0.5_wp )
257         pe3w_1d (jk) =    rn_dz
258         pe3t_1d (jk) =    rn_dz
259      END DO
260      !
261      IF(lwp) THEN                        ! control print
262         WRITE(numout,*)
263         WRITE(numout,*) '              Reference 1D z-coordinate depth and scale factors:'
264         WRITE(numout, "(9x,' level  gdept_1d  gdepw_1d  e3t_1d   e3w_1d  ')" )
265         WRITE(numout, "(10x, i4, 4f9.2)" ) ( jk, pdept_1d(jk), pdepw_1d(jk), pe3t_1d(jk), pe3w_1d(jk), jk = 1, jpk )
266      ENDIF
267      !
268   END SUBROUTINE zgr_z1d
269   
270   !!======================================================================
271END MODULE usrdef_zgr
Note: See TracBrowser for help on using the repository browser.