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.
scoord_gen.F90 in branches/2015/dev_r5187_UKMO13_simplification/NEMOGCM/TOOLS/SCOORD_GEN/src – NEMO

source: branches/2015/dev_r5187_UKMO13_simplification/NEMOGCM/TOOLS/SCOORD_GEN/src/scoord_gen.F90 @ 5295

Last change on this file since 5295 was 5295, checked in by timgraham, 9 years ago

Minor bugfixes for sh94 and tanh stretching function cases

File size: 35.7 KB
Line 
1PROGRAM SCOORD_GEN
2      !!----------------------------------------------------------------------
3      !!                  ***  PROGRAM scoord_gen  ***
4      !!                     
5      !! ** Purpose :   define the s-coordinate system
6      !!
7      !! ** Method  :   s-coordinate
8      !!         The depth of model levels is defined as the product of an
9      !!      analytical function by the local bathymetry, while the vertical
10      !!      scale factors are defined as the product of the first derivative
11      !!      of the analytical function by the bathymetry.
12      !!      (this solution save memory as depth and scale factors are not
13      !!      3d fields)
14      !!          - Read bathymetry (in meters) at t-point and compute the
15      !!         bathymetry at u-, v-, and f-points.
16      !!            hbatu = mi( hbatt )
17      !!            hbatv = mj( hbatt )
18      !!            hbatf = mi( mj( hbatt ) )
19      !!          - Compute z_gsigt, z_gsigw, z_esigt, z_esigw from an analytical
20      !!         function and its derivative given as function.
21      !!            z_gsigt(k) = fssig (k    )
22      !!            z_gsigw(k) = fssig (k-0.5)
23      !!            z_esigt(k) = fsdsig(k    )
24      !!            z_esigw(k) = fsdsig(k-0.5)
25      !!      Three options for stretching are give, and they can be modified
26      !!      following the users requirements. Nevertheless, the output as
27      !!      well as the way to compute the model levels and scale factors
28      !!      must be respected in order to insure second order accuracy
29      !!      schemes.
30      !!
31      !!      The three methods for stretching available are:
32      !!
33      !!           s_sh94 (Song and Haidvogel 1994)
34      !!                a sinh/tanh function that allows sigma and stretched sigma
35      !!
36      !!           s_sf12 (Siddorn and Furner 2012?)
37      !!                allows the maintenance of fixed surface and or
38      !!                bottom cell resolutions (cf. geopotential coordinates)
39      !!                within an analytically derived stretched S-coordinate framework.
40      !!
41      !!          s_tanh  (Madec et al 1996)
42      !!                a cosh/tanh function that gives stretched coordinates       
43      !!
44      !!----------------------------------------------------------------------
45      !
46      USE utils
47      IMPLICIT NONE
48
49
50     !!----------------------------------------------------------------------
51      !
52      !
53      OPEN( UNIT=numnam, FILE='namelist', FORM='FORMATTED', STATUS='OLD' )
54      READ( numnam, namzgr_sco )
55      CLOSE( numnam )
56      jpk = INT(rn_jpk)
57
58      WRITE(*,*)
59      WRITE(*,*) 'domzgr_sco : s-coordinate or hybrid z-s-coordinate'
60      WRITE(*,*) '~~~~~~~~~~~'
61      WRITE(*,*) '   Namelist namzgr_sco'
62      WRITE(*,*) '     stretching coeffs '
63      WRITE(*,*) '        Number of levels                             rn_jpk        = ',rn_jpk
64      WRITE(*,*) '        maximum depth of s-bottom surface (>0)       rn_sbot_max   = ',rn_sbot_max
65      WRITE(*,*) '        minimum depth of s-bottom surface (>0)       rn_sbot_min   = ',rn_sbot_min
66      WRITE(*,*) '        Critical depth                               rn_hc         = ',rn_hc
67      WRITE(*,*) '        maximum cut-off r-value allowed              rn_rmax       = ',rn_rmax
68      WRITE(*,*) '     Song and Haidvogel 1994 stretching              ln_s_sh94     = ',ln_s_sh94
69      WRITE(*,*) '        Song and Haidvogel 1994 stretching coefficients'
70      WRITE(*,*) '        surface control parameter (0<=rn_theta<=20)  rn_theta      = ',rn_theta
71      WRITE(*,*) '        bottom  control parameter (0<=rn_thetb<= 1)  rn_thetb      = ',rn_thetb
72      WRITE(*,*) '        stretching parameter (song and haidvogel)    rn_bb         = ',rn_bb
73      WRITE(*,*) '     Siddorn and Furner 2012 stretching              ln_s_sf12     = ',ln_s_sf12
74      WRITE(*,*) '        switching to sigma (T) or Z (F) at H<Hc      ln_sigcrit    = ',ln_sigcrit
75      WRITE(*,*) '        Siddorn and Furner 2012 stretching coefficients'
76      WRITE(*,*) '        stretchin parameter ( >1 surface; <1 bottom) rn_alpha      = ',rn_alpha
77      WRITE(*,*) '        e-fold length scale for transition region    rn_efold      = ',rn_efold
78      WRITE(*,*) '        Surface cell depth (Zs) (m)                  rn_zs         = ',rn_zs
79      WRITE(*,*) '        Bathymetry multiplier for Zb                 rn_zb_a       = ',rn_zb_a
80      WRITE(*,*) '        Offset for Zb                                rn_zb_b       = ',rn_zb_b
81      WRITE(*,*) '        Bottom cell (Zb) (m) = H*rn_zb_a + rn_zb_b'
82
83     
84      ! Read in bathy, jpi and jpj from bathy.nc
85      CALL read_bathy()
86
87      !Allocate all other allocatable variables
88      ios = dom_oce_alloc()
89      IF( ios .NE. 0) THEN
90         WRITE(0,*) 'Unable to allocate all arrays'
91         STOP 1
92      ENDIF
93
94      hift(:,:) = rn_sbot_min                     ! set the minimum depth for the s-coordinate
95      hifu(:,:) = rn_sbot_min
96      hifv(:,:) = rn_sbot_min
97      hiff(:,:) = rn_sbot_min
98
99      !                                        ! set maximum ocean depth
100      bathy(:,:) = MIN( rn_sbot_max, bathy(:,:) )
101
102      DO jj = 1, jpj
103         DO ji = 1, jpi
104           IF( bathy(ji,jj) > 0. )   bathy(ji,jj) = MAX( rn_sbot_min, bathy(ji,jj) )
105         END DO
106      END DO
107      !                                        ! =============================
108      !                                        ! Define the envelop bathymetry   (hbatt)
109      !                                        ! =============================
110      ! use r-value to create hybrid coordinates
111      zenv(:,:) = bathy(:,:)
112      !
113      ! set first land point adjacent to a wet cell to sbot_min as this needs to be included in smoothing
114      DO jj = 1, jpj
115         DO ji = 1, jpi
116           IF( bathy(ji,jj) == 0. ) THEN
117             iip1 = MIN( ji+1, jpi )
118             ijp1 = MIN( jj+1, jpj )
119             iim1 = MAX( ji-1, 1 )
120             ijm1 = MAX( jj-1, 1 )
121             IF( (bathy(iip1,jj) + bathy(iim1,jj) + bathy(ji,ijp1) + bathy(ji,ijm1) +              &
122        &         bathy(iip1,ijp1) + bathy(iim1,ijm1) + bathy(iip1,ijp1) + bathy(iim1,ijm1)) > 0. ) THEN
123               zenv(ji,jj) = rn_sbot_min
124             ENDIF
125           ENDIF
126         END DO
127      END DO
128      WRITE(*,*) 'domzgr_sco print', bathy(196,147)
129      !
130      ! smooth the bathymetry (if required)
131      scosrf(:,:) = 0.             ! ocean surface depth (here zero: no under ice-shelf sea)
132      scobot(:,:) = bathy(:,:)        ! ocean bottom  depth
133      !
134      jl = 0
135      zrmax = 1.
136      !   
137      !     
138      ! set scaling factor used in reducing vertical gradients
139      zrfact = ( 1. - rn_rmax ) / ( 1. + rn_rmax )
140      !
141      ! initialise temporary evelope depth arrays
142      ztmpi1(:,:) = zenv(:,:)
143      ztmpi2(:,:) = zenv(:,:)
144      ztmpj1(:,:) = zenv(:,:)
145      ztmpj2(:,:) = zenv(:,:)
146      !
147      ! initialise temporary r-value arrays
148      zri(:,:) = 1.
149      zrj(:,:) = 1.
150
151      !                                                            ! ================ !
152      DO WHILE( jl <= 10000 .AND. ( zrmax - rn_rmax ) > 1.e-8 ) !  Iterative loop  !
153         !                                                         ! ================ !
154         jl = jl + 1
155         zrmax = 0.
156         ! we set zrmax from previous r-values (zri and zrj) first
157         ! if set after current r-value calculation (as previously)
158         ! we could exit DO WHILE prematurely before checking r-value
159         ! of current zenv
160         DO jj = 1, jpj 
161            DO ji = 1, jpi
162               zrmax = MAX( zrmax, ABS(zri(ji,jj)), ABS(zrj(ji,jj)) )
163            END DO
164         END DO
165         zri(:,:) = 0.
166         zrj(:,:) = 0.
167         DO jj = 1, jpj
168            DO ji = 1, jpi
169               iip1 = MIN(ji+1,jpi) 
170               ijp1 = MIN(jj+1,jpj)     
171               IF( (zenv(ji,jj) > 0.) .AND. (zenv(iip1,jj) > 0.)) THEN
172                  zri(ji,jj) = ( zenv(iip1,jj  ) - zenv(ji,jj) ) / ( zenv(iip1,jj  ) + zenv(ji,jj) )
173               END IF
174               IF( (zenv(ji,jj) > 0.) .AND. (zenv(ji,ijp1) > 0.)) THEN
175                  zrj(ji,jj) = ( zenv(ji  ,ijp1) - zenv(ji,jj) ) / ( zenv(ji  ,ijp1) + zenv(ji,jj) )
176               END IF
177               IF( zri(ji,jj) >  rn_rmax )   ztmpi1(ji  ,jj  ) = zenv(iip1,jj  ) * zrfact
178               IF( zri(ji,jj) < -rn_rmax )   ztmpi2(iip1,jj  ) = zenv(ji  ,jj  ) * zrfact
179               IF( zrj(ji,jj) >  rn_rmax )   ztmpj1(ji  ,jj  ) = zenv(ji  ,ijp1) * zrfact
180               IF( zrj(ji,jj) < -rn_rmax )   ztmpj2(ji  ,ijp1) = zenv(ji  ,jj  ) * zrfact
181            END DO
182         END DO
183         !
184         WRITE(*,*) 'zgr_sco :   iter= ',jl, ' rmax= ', zrmax
185         !
186         DO jj = 1, jpj
187            DO ji = 1, jpi
188               zenv(ji,jj) = MAX(zenv(ji,jj), ztmpi1(ji,jj), ztmpi2(ji,jj), ztmpj1(ji,jj), ztmpj2(ji,jj) )
189            END DO
190         END DO
191         !                                                  ! ================ !
192      END DO                                                !     End loop     !
193      !                                                     ! ================ !
194      DO jj = 1, jpj
195         DO ji = 1, jpi
196            zenv(ji,jj) = MAX( zenv(ji,jj), rn_sbot_min ) ! set all points to avoid undefined scale value warnings
197         END DO
198      END DO
199      !
200      ! Envelope bathymetry saved in hbatt
201      ! TODO - get this section to work
202      hbatt(:,:) = zenv(:,:) 
203!      IF( MINVAL( gphit(:,:) ) * MAXVAL( gphit(:,:) ) <= 0. ) THEN
204 !        CALL ctl_warn( ' s-coordinates are tapered in vicinity of the Equator' )
205 !        DO jj = 1, jpj
206 !           DO ji = 1, jpi
207 !              ztaper = EXP( -(gphit(ji,jj)/8.)**2. )
208 !              hbatt(ji,jj) = rn_sbot_max * ztaper + hbatt(ji,jj) * ( 1. - ztaper )
209 !           END DO
210 !        END DO
211 !     ENDIF
212      !
213      !                                        ! ==============================
214      !                                        !   hbatu, hbatv, hbatf fields
215      !                                        ! ==============================
216      WRITE(*,*)
217      WRITE(*,*) ' zgr_sco: minimum depth of the envelop topography set to : ', rn_sbot_min
218      hbatu(:,:) = rn_sbot_min
219      hbatv(:,:) = rn_sbot_min
220      hbatf(:,:) = rn_sbot_min
221      DO jj = 1, jpj-1
222        DO ji = 1, jpi-1   ! NO vector opt.
223           hbatu(ji,jj) = 0.50 * ( hbatt(ji  ,jj) + hbatt(ji+1,jj  ) )
224           hbatv(ji,jj) = 0.50 * ( hbatt(ji  ,jj) + hbatt(ji  ,jj+1) )
225           hbatf(ji,jj) = 0.25 * ( hbatt(ji  ,jj) + hbatt(ji  ,jj+1)   &
226              &                     + hbatt(ji+1,jj) + hbatt(ji+1,jj+1) )
227        END DO
228      END DO
229      !
230      zhbat(:,:) = hbatu(:,:)   
231      DO jj = 1, jpj
232         DO ji = 1, jpi
233            IF( hbatu(ji,jj) == 0. ) THEN
234               IF( zhbat(ji,jj) == 0. )   hbatu(ji,jj) = rn_sbot_min
235               IF( zhbat(ji,jj) /= 0. )   hbatu(ji,jj) = zhbat(ji,jj)
236            ENDIF
237         END DO
238      END DO
239      zhbat(:,:) = hbatv(:,:) 
240      DO jj = 1, jpj
241         DO ji = 1, jpi
242            IF( hbatv(ji,jj) == 0. ) THEN
243               IF( zhbat(ji,jj) == 0. )   hbatv(ji,jj) = rn_sbot_min
244               IF( zhbat(ji,jj) /= 0. )   hbatv(ji,jj) = zhbat(ji,jj)
245            ENDIF
246         END DO
247      END DO
248      zhbat(:,:) = hbatf(:,:)   
249      DO jj = 1, jpj
250         DO ji = 1, jpi
251            IF( hbatf(ji,jj) == 0. ) THEN
252               IF( zhbat(ji,jj) == 0. )   hbatf(ji,jj) = rn_sbot_min
253               IF( zhbat(ji,jj) /= 0. )   hbatf(ji,jj) = zhbat(ji,jj)
254            ENDIF
255         END DO
256      END DO
257
258!!bug:  key_helsinki a verifer
259      hift(:,:) = MIN( hift(:,:), hbatt(:,:) )
260      hifu(:,:) = MIN( hifu(:,:), hbatu(:,:) )
261      hifv(:,:) = MIN( hifv(:,:), hbatv(:,:) )
262      hiff(:,:) = MIN( hiff(:,:), hbatf(:,:) )
263
264         WRITE(*,*) ' MAX val hif   t ', MAXVAL( hift (:,:) ), ' f ', MAXVAL( hiff (:,:) ),  &
265            &                        ' u ',   MAXVAL( hifu (:,:) ), ' v ', MAXVAL( hifv (:,:) )
266         WRITE(*,*) ' MIN val hif   t ', MINVAL( hift (:,:) ), ' f ', MINVAL( hiff (:,:) ),  &
267            &                        ' u ',   MINVAL( hifu (:,:) ), ' v ', MINVAL( hifv (:,:) )
268         WRITE(*,*) ' MAX val hbat  t ', MAXVAL( hbatt(:,:) ), ' f ', MAXVAL( hbatf(:,:) ),  &
269            &                        ' u ',   MAXVAL( hbatu(:,:) ), ' v ', MAXVAL( hbatv(:,:) )
270         WRITE(*,*) ' MIN val hbat  t ', MINVAL( hbatt(:,:) ), ' f ', MINVAL( hbatf(:,:) ),  &
271            &                        ' u ',   MINVAL( hbatu(:,:) ), ' v ', MINVAL( hbatv(:,:) )
272     
273      ! Create the output file
274      CALL make_coord_file()
275
276      !                                            ! =======================
277      !                                            !   s-coordinate fields     (gdep., e3.)
278      !                                            ! =======================
279      !
280      ! non-dimensional "sigma" for model level depth at w- and t-levels
281
282
283!========================================================================
284! Song and Haidvogel  1994 (ln_s_sh94=T)
285! Siddorn and Furner 2012 (ln_sf12=T)
286! or  tanh function       (both false)                   
287! To reduce memory loop over jpk and write each level to file
288!========================================================================
289      IF      ( ln_s_sh94 ) THEN
290                           CALL s_sh94()
291      ELSE IF ( ln_s_sf12 ) THEN
292                           CALL s_sf12()
293      ELSE                 
294                           CALL s_tanh()
295      ENDIF
296
297      CALL check_nf90( nf90_put_var(ncout, var_ids(11), mbathy) )
298      CALL check_nf90( nf90_close(ncout) )
299
300 
301      !
302END PROGRAM SCOORD_GEN
303
304!!======================================================================
305   SUBROUTINE s_sh94()
306
307      !!----------------------------------------------------------------------
308      !!                  ***  ROUTINE s_sh94  ***
309      !!                     
310      !! ** Purpose :   stretch the s-coordinate system
311      !!
312      !! ** Method  :   s-coordinate stretch using the Song and Haidvogel 1994
313      !!                mixed S/sigma coordinate
314      !!
315      !! Reference : Song and Haidvogel 1994.
316      !!----------------------------------------------------------------------
317      !
318      USE utils
319      REAL(wp) ::   zcoeft, zcoefw   ! temporary scalars
320      !
321      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: z_gsigw3, z_gsigt3, z_gsi3w3
322      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: z_gsigw3p1, z_gsigt3m1, z_gsi3w3m1
323      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: z_esigt3, z_esigw3, z_esigtu3 
324      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: z_esigtv3, z_esigtf3, z_esigwu3, z_esigwv3           
325
326      ALLOCATE( z_gsigw3(jpi,jpj), z_gsigt3(jpi,jpj), z_gsi3w3(jpi,jpj) )
327      ALLOCATE( z_esigt3(jpi,jpj), z_esigw3(jpi,jpj), z_esigtu3(jpi,jpj)) 
328      ALLOCATE( z_esigtv3(jpi,jpj), z_esigtf3(jpi,jpj), z_esigwu3(jpi,jpj)) 
329      ALLOCATE( z_esigwv3(jpi,jpj), z_gsigw3p1(jpi,jpj), z_gsigt3m1(jpi,jpj) )
330      ALLOCATE( z_gsi3w3m1(jpi,jpj) )
331 
332      z_gsigw3  = 0.   ;   z_gsigt3  = 0.   ;   z_gsi3w3  = 0.
333      z_esigt3  = 0.   ;   z_esigw3  = 0. 
334      z_esigt3p1  = 0. ;   z_esigw3p1  = 0. 
335      z_esigtu3 = 0.   ;   z_esigtv3 = 0.   ;   z_esigtf3 = 0.
336      z_esigwu3 = 0.   ;   z_esigwv3 = 0.
337
338      DO jk = 1,jpk
339         DO ji = 1, jpi
340            DO jj = 1, jpj
341
342               IF( hbatt(ji,jj) > rn_hc ) THEN    !deep water, stretched sigma
343                     z_gsigw3(ji,jj) = -fssig1( REAL(jk,wp)-0.5, rn_bb )
344                     z_gsigw3p1(ji,jj) = -fssig1( REAL(jk+1,wp)-0.5, rn_bb )
345                     z_gsigt3(ji,jj) = -fssig1( REAL(jk,wp)    , rn_bb )
346               ELSE ! shallow water, uniform sigma
347                     z_gsigw3(ji,jj) =   REAL(jk-1,wp)            / REAL(jpk-1,wp)
348                     z_gsigw3p1(ji,jj) =   REAL(jk,wp)            / REAL(jpk-1,wp)
349                     z_gsigt3(ji,jj) = ( REAL(jk-1,wp) + 0.5 ) / REAL(jpk-1,wp)
350               ENDIF
351               !
352             !gsi3w3m1 & gsigt3m1 only used if jk /= 1 and is set at the end of the loop over jk
353               IF( jk .EQ. 1) THEN
354                  z_esigw3(ji,jj ) = 2. * ( z_gsigt3(ji,jj ) - z_gsigw3(ji,jj ) )
355                  z_gsi3w3(ji,jj) = 0.5 * z_esigw3(ji,jj)
356               ELSE
357                  z_esigw3(ji,jj) = z_gsigt3(ji,jj) - z_gsigt3m1(ji,jj)
358                  z_gsi3w3(ji,jj) = z_gsi3w3m1(ji,jj) + z_esigw3(ji,jj)
359               ENDIF
360               IF( jk .EQ. jpk) THEN
361                  z_esigt3(ji,jj) = 2. * ( z_gsigt3(ji,jj) - z_gsigw3(ji,jj) )
362               ELSE
363                  z_esigt3(ji,jj ) = z_gsigw3p1(ji,jj) - z_gsigw3(ji,jj)
364               ENDIF
365               !
366
367               zcoeft = ( REAL(jk,wp) - 0.5 ) / REAL(jpk-1,wp)
368               zcoefw = ( REAL(jk,wp) - 1.0 ) / REAL(jpk-1,wp)
369               gdept_0(ji,jj) = ( scosrf(ji,jj) + (hbatt(ji,jj)-rn_hc)*z_gsigt3(ji,jj)+rn_hc*zcoeft )
370               gdepw_0(ji,jj) = ( scosrf(ji,jj) + (hbatt(ji,jj)-rn_hc)*z_gsigw3(ji,jj)+rn_hc*zcoefw )
371               gdep3w_0(ji,jj) = ( scosrf(ji,jj) + (hbatt(ji,jj)-rn_hc)*z_gsi3w3(ji,jj)+rn_hc*zcoeft )
372              !
373            END DO   ! for all jj's
374         END DO    ! for all ji's
375
376         DO ji = 1, jpi-1
377            DO jj = 1, jpj-1
378               z_esigtu3(ji,jj) = ( hbatt(ji,jj)*z_esigt3(ji,jj)+hbatt(ji+1,jj)*z_esigt3(ji+1,jj) )   &
379                  &              / ( hbatt(ji,jj)+hbatt(ji+1,jj) )
380               z_esigtv3(ji,jj) = ( hbatt(ji,jj)*z_esigt3(ji,jj)+hbatt(ji,jj+1)*z_esigt3(ji,jj+1) )   &
381                  &              / ( hbatt(ji,jj)+hbatt(ji,jj+1) )
382               z_esigtf3(ji,jj) = ( hbatt(ji,jj)*z_esigt3(ji,jj)+hbatt(ji+1,jj)*z_esigt3(ji+1,jj)     &
383                  &                + hbatt(ji,jj+1)*z_esigt3(ji,jj+1)+hbatt(ji+1,jj+1)*z_esigt3(ji+1,jj+1) )   &
384                  &              / ( hbatt(ji,jj)+hbatt(ji+1,jj)+hbatt(ji,jj+1)+hbatt(ji+1,jj+1) )
385               z_esigwu3(ji,jj) = ( hbatt(ji,jj)*z_esigw3(ji,jj)+hbatt(ji+1,jj)*z_esigw3(ji+1,jj) )   &
386                  &              / ( hbatt(ji,jj)+hbatt(ji+1,jj) )
387               z_esigwv3(ji,jj) = ( hbatt(ji,jj)*z_esigw3(ji,jj)+hbatt(ji,jj+1)*z_esigw3(ji,jj+1) )   &
388                  &              / ( hbatt(ji,jj)+hbatt(ji,jj+1) )
389               !
390               e3t_0(ji,jj) = ( (hbatt(ji,jj)-rn_hc)*z_esigt3 (ji,jj) + rn_hc/REAL(jpk-1,wp) )
391               e3u_0(ji,jj) = ( (hbatu(ji,jj)-rn_hc)*z_esigtu3(ji,jj) + rn_hc/REAL(jpk-1,wp) )
392               e3v_0(ji,jj) = ( (hbatv(ji,jj)-rn_hc)*z_esigtv3(ji,jj) + rn_hc/REAL(jpk-1,wp) )
393               e3f_0(ji,jj) = ( (hbatf(ji,jj)-rn_hc)*z_esigtf3(ji,jj) + rn_hc/REAL(jpk-1,wp) )
394               !
395               e3w_0 (ji,jj) = ( (hbatt(ji,jj)-rn_hc)*z_esigw3 (ji,jj) + rn_hc/REAL(jpk-1,wp) )
396               e3uw_0(ji,jj) = ( (hbatu(ji,jj)-rn_hc)*z_esigwu3(ji,jj) + rn_hc/REAL(jpk-1,wp) )
397               e3vw_0(ji,jj) = ( (hbatv(ji,jj)-rn_hc)*z_esigwv3(ji,jj) + rn_hc/REAL(jpk-1,wp) )
398            END DO
399         END DO
400
401         z_gsigt3m1 = z_gsigt3
402         z_gsiw3m1 = z_gsiw3
403 
404         where (e3t_0   (:,:).eq.0.0)  e3t_0(:,:) = 1.0
405         where (e3u_0   (:,:).eq.0.0)  e3u_0(:,:) = 1.0
406         where (e3v_0   (:,:).eq.0.0)  e3v_0(:,:) = 1.0
407         where (e3f_0   (:,:).eq.0.0)  e3f_0(:,:) = 1.0
408         where (e3w_0   (:,:).eq.0.0)  e3w_0(:,:) = 1.0
409         where (e3uw_0  (:,:).eq.0.0)  e3uw_0(:,:) = 1.0
410         where (e3vw_0  (:,:).eq.0.0)  e3vw_0(:,:) = 1.0
411       
412         CALL write_netcdf_vars(jk)
413         DO jj = 1, jpj
414            DO ji = 1, jpi
415                  IF( scobot(ji,jj) >= gdept_0(ji,jj) )   mbathy(ji,jj) = MAX( 2, jk )
416                  IF( scobot(ji,jj) == 0.             )   mbathy(ji,jj) = 0
417            END DO
418         END DO
419      END DO !End of loop over jk
420
421      DEALLOCATE( z_gsigw3, z_gsigt3, z_gsi3w3, z_gsigw3p1, z_gsigt3m1, z_gsi3w3m1)
422      DEALLOCATE( z_esigt3, z_esigw3, z_esigtu3, z_esigtv3, z_esigtf3, z_esigwu3, z_esigwv3 )
423
424   END SUBROUTINE s_sh94
425
426   SUBROUTINE s_sf12
427
428      !!----------------------------------------------------------------------
429      !!                  ***  ROUTINE s_sf12 ***
430      !!                     
431      !! ** Purpose :   stretch the s-coordinate system
432      !!
433      !! ** Method  :   s-coordinate stretch using the Siddorn and Furner 2012?
434      !!                mixed S/sigma/Z coordinate
435      !!
436      !!                This method allows the maintenance of fixed surface and or
437      !!                bottom cell resolutions (cf. geopotential coordinates)
438      !!                within an analytically derived stretched S-coordinate framework.
439      !!
440      !!
441      !! Reference : Siddorn and Furner 2012 (submitted Ocean modelling).
442      !!----------------------------------------------------------------------
443      !
444      USE utils
445      REAL(wp) ::   zsmth               ! smoothing around critical depth
446      REAL(wp) ::   zzs, zzb           ! Surface and bottom cell thickness in sigma space
447      !
448      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: z_gsigw3, z_gsigt3, z_gsi3w3
449      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: z_gsigw3p1, z_gsigt3m1, z_gsi3w3m1
450      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: z_esigt3, z_esigw3, z_esigtu3 
451      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: z_esigtv3, z_esigtf3, z_esigwu3, z_esigwv3           
452
453      ALLOCATE( z_gsigw3(jpi,jpj), z_gsigt3(jpi,jpj), z_gsi3w3(jpi,jpj) )
454      ALLOCATE( z_esigt3(jpi,jpj), z_esigw3(jpi,jpj), z_esigtu3(jpi,jpj)) 
455      ALLOCATE( z_esigtv3(jpi,jpj), z_esigtf3(jpi,jpj), z_esigwu3(jpi,jpj)) 
456      ALLOCATE( z_esigwv3(jpi,jpj), z_gsigw3p1(jpi,jpj), z_gsigt3m1(jpi,jpj) )
457      ALLOCATE( z_gsi3w3m1(jpi,jpj) )
458
459      z_gsigw3  = 0.   ;   z_gsigt3  = 0.   ;   z_gsi3w3  = 0.
460      z_esigt3  = 0.   ;   z_esigw3  = 0. 
461      z_esigtu3 = 0.   ;   z_esigtv3 = 0.   ;   z_esigtf3 = 0.
462      z_esigwu3 = 0.   ;   z_esigwv3 = 0.
463     
464
465      DO jk = 1,jpk
466         DO ji = 1, jpi
467            DO jj = 1, jpj
468             IF (hbatt(ji,jj)>rn_hc) THEN !deep water, stretched sigma
469             
470                zzb = hbatt(ji,jj)*rn_zb_a + rn_zb_b   ! this forces a linear bottom cell depth relationship with H,.
471                                                     ! could be changed by users but care must be taken to do so carefully
472                zzb = 1.0-(zzb/hbatt(ji,jj))
473           
474                zzs = rn_zs / hbatt(ji,jj) 
475             
476                IF (rn_efold /= 0.0) THEN
477                   zsmth   = tanh( (hbatt(ji,jj)- rn_hc ) / rn_efold )
478                ELSE
479                   zsmth = 1.0 
480                ENDIF
481                 
482                z_gsigw3(ji,jj) =  REAL(jk-1,wp)        /REAL(jpk-1,wp)
483                z_gsigw3p1(ji,jj) =  REAL(jk,wp)        /REAL(jpk-1,wp)
484                z_gsigt3(ji,jj) = (REAL(jk-1,wp)+0.5)/REAL(jpk-1,wp)
485                z_gsigw3(ji,jj) = fgamma( z_gsigw3(ji,jj), zzb, zzs, zsmth  )
486                z_gsigw3p1(ji,jj) = fgamma( z_gsigw3p1(ji,jj), zzb, zzs, zsmth  )
487                z_gsigt3(ji,jj) = fgamma( z_gsigt3(ji,jj), zzb, zzs, zsmth  )
488 
489             ELSE IF(ln_sigcrit) THEN ! shallow water, uniform sigma
490
491                z_gsigw3(ji,jj) =  REAL(jk-1,wp)     /REAL(jpk-1,wp)
492                z_gsigw3p1(ji,jj) =  REAL(jk,wp)     /REAL(jpk-1,wp)
493                z_gsigt3(ji,jj) = (REAL(jk-1,wp)+0.5)/REAL(jpk-1,wp)
494
495             ELSE  ! shallow water, z coordinates
496
497               z_gsigw3(ji,jj) =  REAL(jk-1,wp)        /REAL(jpk-1,wp)*(rn_hc/hbatt(ji,jj))
498               z_gsigw3p1(ji,jj) =  REAL(jk,wp)        /REAL(jpk-1,wp)*(rn_hc/hbatt(ji,jj))
499               z_gsigt3(ji,jj) = (REAL(jk-1,wp)+0.5)/REAL(jpk-1,wp)*(rn_hc/hbatt(ji,jj))
500
501             ENDIF
502
503             !gsi3w3m1 & z_gsigt3m1 only used if jk /= 1 and is set at the end of the loop over jk
504             IF( jk .EQ. 1) THEN
505                z_esigw3(ji,jj  ) = 2.0 * (z_gsigt3(ji,jj ) - z_gsigw3(ji,jj ))
506                z_gsi3w3(ji,jj) = 0.5 * z_esigw3(ji,jj)
507             ELSE
508                z_esigw3(ji,jj) = z_gsigt3(ji,jj) - z_gsigt3m1(ji,jj)
509                z_gsi3w3(ji,jj) = z_gsi3w3m1(ji,jj) + z_esigw3(ji,jj)
510             ENDIF
511             IF( jk .EQ. jpk) THEN
512                z_esigt3(ji,jj) = 2.0 * (z_gsigt3(ji,jj) - z_gsigw3(ji,jj))
513             ELSE
514                z_esigt3(ji,jj) = z_gsigw3p1(ji,jj) - z_gsigw3(ji,jj)
515             ENDIF
516
517             gdept_0(ji,jj) = (scosrf(ji,jj)+hbatt(ji,jj))*z_gsigt3(ji,jj)
518             gdepw_0(ji,jj) = (scosrf(ji,jj)+hbatt(ji,jj))*z_gsigw3(ji,jj)
519             gdep3w_0(ji,jj) = (scosrf(ji,jj)+hbatt(ji,jj))*z_gsi3w3(ji,jj)
520
521            ENDDO   ! for all jj's
522         ENDDO    ! for all ji's
523
524         DO ji=1,jpi-1
525            DO jj=1,jpj-1
526
527               z_esigtu3(ji,jj) = ( hbatt(ji,jj)*z_esigt3(ji,jj)+hbatt(ji+1,jj)*z_esigt3(ji+1,jj) ) / &
528                                     ( hbatt(ji,jj)+hbatt(ji+1,jj) )
529               z_esigtv3(ji,jj) = ( hbatt(ji,jj)*z_esigt3(ji,jj)+hbatt(ji,jj+1)*z_esigt3(ji,jj+1) ) / &
530                                     ( hbatt(ji,jj)+hbatt(ji,jj+1) )
531               z_esigtf3(ji,jj) = ( hbatt(ji,jj)*z_esigt3(ji,jj)+hbatt(ji+1,jj)*z_esigt3(ji+1,jj) +  &
532                                       hbatt(ji,jj+1)*z_esigt3(ji,jj+1)+hbatt(ji+1,jj+1)*z_esigt3(ji+1,jj+1) ) / &
533                                     ( hbatt(ji,jj)+hbatt(ji+1,jj)+hbatt(ji,jj+1)+hbatt(ji+1,jj+1) )
534               z_esigwu3(ji,jj) = ( hbatt(ji,jj)*z_esigw3(ji,jj)+hbatt(ji+1,jj)*z_esigw3(ji+1,jj) ) / &
535                                     ( hbatt(ji,jj)+hbatt(ji+1,jj) )
536               z_esigwv3(ji,jj) = ( hbatt(ji,jj)*z_esigw3(ji,jj)+hbatt(ji,jj+1)*z_esigw3(ji,jj+1) ) / &
537                                     ( hbatt(ji,jj)+hbatt(ji,jj+1) )
538
539               e3t_0(ji,jj)=(scosrf(ji,jj)+hbatt(ji,jj))*z_esigt3(ji,jj)
540               e3u_0(ji,jj)=(scosrf(ji,jj)+hbatu(ji,jj))*z_esigtu3(ji,jj)
541               e3v_0(ji,jj)=(scosrf(ji,jj)+hbatv(ji,jj))*z_esigtv3(ji,jj)
542               e3f_0(ji,jj)=(scosrf(ji,jj)+hbatf(ji,jj))*z_esigtf3(ji,jj)
543               !
544               e3w_0(ji,jj)=hbatt(ji,jj)*z_esigw3(ji,jj)
545               e3uw_0(ji,jj)=hbatu(ji,jj)*z_esigwu3(ji,jj)
546               e3vw_0(ji,jj)=hbatv(ji,jj)*z_esigwv3(ji,jj)
547            ENDDO
548         ENDDO
549         ! Keep some arrays for next level
550         z_gsigt3m1 = z_gsigt3
551         z_gsiw3m1 = z_gsiw3
552 
553         where (e3t_0   (:,:).eq.0.0)  e3t_0(:,:) = 1.0
554         where (e3u_0   (:,:).eq.0.0)  e3u_0(:,:) = 1.0
555         where (e3v_0   (:,:).eq.0.0)  e3v_0(:,:) = 1.0
556         where (e3f_0   (:,:).eq.0.0)  e3f_0(:,:) = 1.0
557         where (e3w_0   (:,:).eq.0.0)  e3w_0(:,:) = 1.0
558         where (e3uw_0  (:,:).eq.0.0)  e3uw_0(:,:) = 1.0
559         where (e3vw_0  (:,:).eq.0.0)  e3vw_0(:,:) = 1.0
560         
561         CALL write_netcdf_vars(jk)
562
563         DO jj = 1, jpj
564            DO ji = 1, jpi
565                  IF( scobot(ji,jj) >= gdept_0(ji,jj) )   mbathy(ji,jj) = MAX( 2, jk )
566                  IF( scobot(ji,jj) == 0.             )   mbathy(ji,jj) = 0
567            END DO
568         END DO
569
570      ENDDO ! End of loop over jk
571
572      DEALLOCATE( z_gsigw3, z_gsigt3, z_gsi3w3, z_gsigw3p1, z_gsigt3m1, z_gsi3w3m1)
573      DEALLOCATE( z_esigt3, z_esigw3, z_esigtu3, z_esigtv3, z_esigtf3, z_esigwu3, z_esigwv3 )
574
575   END SUBROUTINE s_sf12
576
577   SUBROUTINE s_tanh()
578
579      !!----------------------------------------------------------------------
580      !!                  ***  ROUTINE s_tanh***
581      !!                     
582      !! ** Purpose :   stretch the s-coordinate system
583      !!
584      !! ** Method  :   s-coordinate stretch
585      !!
586      !! Reference : Madec, Lott, Delecluse and Crepon, 1996. JPO, 26, 1393-1408.
587      !!----------------------------------------------------------------------
588
589      USE utils
590      REAL(wp) ::   zcoeft, zcoefw   ! temporary scalars
591
592      REAL(wp), ALLOCATABLE, DIMENSION(:) :: z_gsigw, z_gsigt, z_gsi3w
593      REAL(wp), ALLOCATABLE, DIMENSION(:) :: z_esigt, z_esigw
594
595      ALLOCATE( z_gsigw(jpk), z_gsigt(jpk), z_gsi3w(jpk) )
596      ALLOCATE( z_esigt(jpk), z_esigw(jpk) )
597
598      z_gsigw  = 0.   ;   z_gsigt  = 0.   ;   z_gsi3w  = 0.
599      z_esigt  = 0.   ;   z_esigw  = 0. 
600
601      DO jk = 1, jpk
602        z_gsigw(jk) = -fssig( REAL(jk,wp)-0.5 )
603        z_gsigt(jk) = -fssig( REAL(jk,wp)        )
604      END DO
605      WRITE(*,*) 'z_gsigw 1 jpk    ', z_gsigw(1), z_gsigw(jpk)
606      !
607      ! Coefficients for vertical scale factors at w-, t- levels
608!!gm bug :  define it from analytical function, not like juste bellow....
609!!gm        or betteroffer the 2 possibilities....
610      DO jk = 1, jpk-1
611         z_esigt(jk  ) = z_gsigw(jk+1) - z_gsigw(jk)
612         z_esigw(jk+1) = z_gsigt(jk+1) - z_gsigt(jk)
613      END DO
614      z_esigw( 1 ) = 2. * ( z_gsigt(1  ) - z_gsigw(1  ) ) 
615      z_esigt(jpk) = 2. * ( z_gsigt(jpk) - z_gsigw(jpk) )
616      !
617      ! Coefficients for vertical depth as the sum of e3w scale factors
618      z_gsi3w(1) = 0.5 * z_esigw(1)
619      DO jk = 2, jpk
620         z_gsi3w(jk) = z_gsi3w(jk-1) + z_esigw(jk)
621      END DO
622!!gm: depuw, depvw can be suppressed (modif in ldfslp) and depw=dep3w can be set (save 3 3D arrays)
623      DO jk = 1, jpk
624      END DO
625!!gm: e3uw, e3vw can be suppressed  (modif in dynzdf, dynzdf_iso, zdfbfr) (save 2 3D arrays)
626      DO jk = 1, jpk
627         DO jj = 1, jpj
628            DO ji = 1, jpi
629              zcoeft = ( REAL(jk,wp) - 0.5 ) / REAL(jpk-1,wp)
630              zcoefw = ( REAL(jk,wp) - 1.0 ) / REAL(jpk-1,wp)
631              gdept_0(ji,jj) = ( scosrf(ji,jj) + (hbatt(ji,jj)-hift(ji,jj))*z_gsigt(jk) + hift(ji,jj)*zcoeft )
632              gdepw_0(:,:) = ( scosrf(ji,jj) + (hbatt(ji,jj)-hift(ji,jj))*z_gsigw(jk) + hift(ji,jj)*zcoefw )
633              gdep3w_0(:,:) = ( scosrf(ji,jj) + (hbatt(ji,jj)-hift(ji,jj))*z_gsi3w(jk) + hift(ji,jj)*zcoeft )
634              e3t_0(ji,jj) = ( (hbatt(ji,jj)-hift(ji,jj))*z_esigt(jk) + hift(ji,jj)/REAL(jpk-1,wp) )
635              e3u_0(ji,jj) = ( (hbatu(ji,jj)-hifu(ji,jj))*z_esigt(jk) + hifu(ji,jj)/REAL(jpk-1,wp) )
636              e3v_0(ji,jj) = ( (hbatv(ji,jj)-hifv(ji,jj))*z_esigt(jk) + hifv(ji,jj)/REAL(jpk-1,wp) )
637              e3f_0(ji,jj) = ( (hbatf(ji,jj)-hiff(ji,jj))*z_esigt(jk) + hiff(ji,jj)/REAL(jpk-1,wp) )
638              !
639              e3w_0(ji,jj) = ( (hbatt(ji,jj)-hift(ji,jj))*z_esigw(jk) + hift(ji,jj)/REAL(jpk-1,wp) )
640              e3uw_0(ji,jj) = ( (hbatu(ji,jj)-hifu(ji,jj))*z_esigw(jk) + hifu(ji,jj)/REAL(jpk-1,wp) )
641              e3vw_0(ji,jj) = ( (hbatv(ji,jj)-hifv(ji,jj))*z_esigw(jk) + hifv(ji,jj)/REAL(jpk-1,wp) )
642              IF( scobot(ji,jj) >= gdept_0(ji,jj) )   mbathy(ji,jj) = MAX( 2, jk )
643              IF( scobot(ji,jj) == 0.             )   mbathy(ji,jj) = 0
644            END DO
645         END DO
646         where (e3t_0   (:,:).eq.0.0)  e3t_0(:,:) = 1.0
647         where (e3u_0   (:,:).eq.0.0)  e3u_0(:,:) = 1.0
648         where (e3v_0   (:,:).eq.0.0)  e3v_0(:,:) = 1.0
649         where (e3f_0   (:,:).eq.0.0)  e3f_0(:,:) = 1.0
650         where (e3w_0   (:,:).eq.0.0)  e3w_0(:,:) = 1.0
651         where (e3uw_0  (:,:).eq.0.0)  e3uw_0(:,:) = 1.0
652         where (e3vw_0  (:,:).eq.0.0)  e3vw_0(:,:) = 1.0
653 
654         CALL write_netcdf_vars(jk)
655      ENDDO ! End of loop over jk
656
657
658      DEALLOCATE( z_gsigw, z_gsigt, z_gsi3w )
659      DEALLOCATE( z_esigt, z_esigw )
660
661   END SUBROUTINE s_tanh
662
663   FUNCTION fssig( pk ) RESULT( pf )
664      !!----------------------------------------------------------------------
665      !!                 ***  ROUTINE fssig ***
666      !!       
667      !! ** Purpose :   provide the analytical function in s-coordinate
668      !!         
669      !! ** Method  :   the function provide the non-dimensional position of
670      !!                T and W (i.e. between 0 and 1)
671      !!                T-points at integer values (between 1 and jpk)
672      !!                W-points at integer values - 1/2 (between 0.5 and jpk-0.5)
673      !!----------------------------------------------------------------------
674      USE utils, ONLY : wp,rn_theta,rn_thetb,jpk
675      IMPLICIT NONE
676      REAL(wp), INTENT(in) ::   pk   ! continuous "k" coordinate
677      REAL(wp)             ::   pf   ! sigma value
678      !!----------------------------------------------------------------------
679      !
680      pf =   (   TANH( rn_theta * ( -(pk-0.5) / REAL(jpk-1) + rn_thetb )  )   &
681         &     - TANH( rn_thetb * rn_theta                                )  )   &
682         & * (   COSH( rn_theta                           )                      &
683         &     + COSH( rn_theta * ( 2. * rn_thetb - 1. ) )  )              &
684         & / ( 2. * SINH( rn_theta ) )
685      !
686   END FUNCTION fssig
687
688
689   FUNCTION fssig1( pk1, pbb ) RESULT( pf1 )
690      !!----------------------------------------------------------------------
691      !!                 ***  ROUTINE fssig1 ***
692      !!
693      !! ** Purpose :   provide the Song and Haidvogel version of the analytical function in s-coordinate
694      !!
695      !! ** Method  :   the function provides the non-dimensional position of
696      !!                T and W (i.e. between 0 and 1)
697      !!                T-points at integer values (between 1 and jpk)
698      !!                W-points at integer values - 1/2 (between 0.5 and jpk-0.5)
699      !!----------------------------------------------------------------------
700      USE utils, ONLY : wp, jpk, rn_theta
701      IMPLICIT NONE
702      REAL(wp), INTENT(in) ::   pk1   ! continuous "k" coordinate
703      REAL(wp), INTENT(in) ::   pbb   ! Stretching coefficient
704      REAL(wp)             ::   pf1   ! sigma value
705      !!----------------------------------------------------------------------
706      !
707      IF ( rn_theta == 0 ) then      ! uniform sigma
708         pf1 = - ( pk1 - 0.5 ) / REAL( jpk-1 )
709      ELSE                        ! stretched sigma
710         pf1 =   ( 1. - pbb ) * ( SINH( rn_theta*(-(pk1-0.5)/REAL(jpk-1)) ) ) / SINH( rn_theta )              &
711            &  + pbb * (  (TANH( rn_theta*( (-(pk1-0.5)/REAL(jpk-1)) + 0.5) ) - TANH( 0.5 * rn_theta )  )  &
712            &        / ( 2. * TANH( 0.5 * rn_theta ) )  )
713      ENDIF
714      !
715   END FUNCTION fssig1
716
717
718   FUNCTION fgamma( pk1, pzb, pzs, psmth) RESULT( p_gamma )
719      !!----------------------------------------------------------------------
720      !!                 ***  ROUTINE fgamma  ***
721      !!
722      !! ** Purpose :   provide analytical function for the s-coordinate
723      !!
724      !! ** Method  :   the function provides the non-dimensional position of
725      !!                T and W (i.e. between 0 and 1)
726      !!                T-points at integer values (between 1 and jpk)
727      !!                W-points at integer values - 1/2 (between 0.5 and jpk-0.5)
728      !!
729      !!                This method allows the maintenance of fixed surface and or
730      !!                bottom cell resolutions (cf. geopotential coordinates)
731      !!                within an analytically derived stretched S-coordinate framework.
732      !!
733      !! Reference  :   Siddorn and Furner, in prep
734      !!----------------------------------------------------------------------
735      USE utils, ONLY : jpk,wp,rn_alpha
736      IMPLICIT NONE
737      REAL(wp), INTENT(in   ) ::   pk1           ! continuous "k" coordinate
738      REAL(wp)                ::   p_gamma       ! stretched coordinate
739      REAL(wp), INTENT(in   ) ::   pzb           ! Bottom box depth
740      REAL(wp), INTENT(in   ) ::   pzs           ! surface box depth
741      REAL(wp), INTENT(in   ) ::   psmth         ! Smoothing parameter
742      REAL(wp)                ::   za1,za2,za3   ! local variables
743      REAL(wp)                ::   zn1,zn2       ! local variables
744      REAL(wp)                ::   za,zb,zx      ! local variables
745      !!----------------------------------------------------------------------
746      !
747
748      zn1  =  1./(jpk-1.)
749      zn2  =  1. -  zn1
750
751      za1 = (rn_alpha+2.0)*zn1**(rn_alpha+1.0)-(rn_alpha+1.0)*zn1**(rn_alpha+2.0) 
752      za2 = (rn_alpha+2.0)*zn2**(rn_alpha+1.0)-(rn_alpha+1.0)*zn2**(rn_alpha+2.0)
753      za3 = (zn2**3.0 - za2)/( zn1**3.0 - za1)
754     
755      za = pzb - za3*(pzs-za1)-za2
756      za = za/( zn2-0.5*(za2+zn2**2.0) - za3*(zn1-0.5*(za1+zn1**2.0) ) )
757      zb = (pzs - za1 - za*( zn1-0.5*(za1+zn1**2.0 ) ) ) / (zn1**3.0 - za1)
758      zx = 1.0-za/2.0-zb
759
760      p_gamma = za*(pk1*(1.0-pk1/2.0))+zb*pk1**3.0 +  &
761                  & zx*( (rn_alpha+2.0)*pk1**(rn_alpha+1.0)- &
762                  &      (rn_alpha+1.0)*pk1**(rn_alpha+2.0) )
763      p_gamma = p_gamma*psmth+pk1*(1.0-psmth)
764
765      !
766   END FUNCTION fgamma
767
Note: See TracBrowser for help on using the repository browser.