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_r12377_KERNEL-06_techene_e3/tests/OVERFLOW/MY_SRC – NEMO

source: NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/tests/OVERFLOW/MY_SRC/usrdef_zgr.F90 @ 12779

Last change on this file since 12779 was 12377, checked in by acc, 4 years ago

The big one. Merging all 2019 developments from the option 1 branch back onto the trunk.

This changeset reproduces 2019/dev_r11943_MERGE_2019 on the trunk using a 2-URL merge
onto a working copy of the trunk. I.e.:

svn merge --ignore-ancestry \

svn+ssh://acc@forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/NEMO/trunk \
svn+ssh://acc@forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/NEMO/branches/2019/dev_r11943_MERGE_2019 ./

The --ignore-ancestry flag avoids problems that may otherwise arise from the fact that
the merge history been trunk and branch may have been applied in a different order but
care has been taken before this step to ensure that all applicable fixes and updates
are present in the merge branch.

The trunk state just before this step has been branched to releases/release-4.0-HEAD
and that branch has been immediately tagged as releases/release-4.0.2. Any fixes
or additions in response to tickets on 4.0, 4.0.1 or 4.0.2 should be done on
releases/release-4.0-HEAD. From now on future 'point' releases (e.g. 4.0.2) will
remain unchanged with periodic releases as needs demand. Note release-4.0-HEAD is a
transitional naming convention. Future full releases, say 4.2, will have a release-4.2
branch which fulfills this role and the first point release (e.g. 4.2.0) will be made
immediately following the release branch creation.

2020 developments can be started from any trunk revision later than this one.

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