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.
vgrid.f90 in branches/2015/nemo_v3_6_STABLE/NEMOGCM/TOOLS/SIREN/src – NEMO

source: branches/2015/nemo_v3_6_STABLE/NEMOGCM/TOOLS/SIREN/src/vgrid.f90 @ 5608

Last change on this file since 5608 was 5608, checked in by jpaul, 9 years ago

commit changes/bugfix/... for SIREN; see ticket #1580

File size: 32.0 KB
Line 
1!----------------------------------------------------------------------
2! NEMO system team, System and Interface for oceanic RElocable Nesting
3!----------------------------------------------------------------------
4!
5! MODULE: vgrid
6!
7! DESCRIPTION:
8!> @brief This module manage vertical grid.
9!>
10!> @details
11!>    to set the depth of model levels and the resulting vertical scale
12!> factors:<br/>
13!> @code
14!>    CALL vgrid_zgr_z(dd_gdepw(:), dd_gdept(:), dd_e3w(:), dd_e3t(:),
15!>                     dd_ppkth, dd_ppkth2, dd_ppacr, dd_ppacr2,
16!>                     dd_ppdzmin, dd_pphmax, dd_pp_to_be_computed,
17!>                     dd_ppa0, dd_ppa1, dd_ppa2, dd_ppsur)
18!> @endcode
19!>       - dd_gdepw is array of depth value on W point
20!>       - dd_gdept is array of depth value on T point
21!>       - dd_e3w   is array of vertical mesh size on W point
22!>       - dd_e3t   is array of vertical mesh size on T point
23!>       - dd_ppkth              see NEMO documentation
24!>       - dd_ppkth2             see NEMO documentation
25!>       - dd_ppacr              see NEMO documentation
26!>       - dd_ppdzmin            see NEMO documentation
27!>       - dd_pphmax             see NEMO documentation
28!>       - dd_pp_to_be_computed  see NEMO documentation
29!>       - dd_ppa1               see NEMO documentation
30!>       - dd_ppa2               see NEMO documentation
31!>       - dd_ppa0               see NEMO documentation
32!>       - dd_ppsur              see NEMO documentation
33!>   
34!>
35!>    to set the depth and vertical scale factor in partial step z-coordinate
36!>  case:<br/>
37!> @code
38!>    CALL vgrid_zgr_zps(id_mbathy(:,:), dd_bathy(:,:), id_jpkmax, dd_gdepw(:),
39!>                       dd_e3t(:), dd_e3zps_min, dd_e3zps_rat)
40!> @endcode
41!>       - id_mbathy is array of bathymetry level
42!>       - dd_bathy  is array of bathymetry
43!>       - id_jpkmax is the maximum number of level to be used
44!>       - dd_gdepw  is array of vertical mesh size on W point
45!>       - dd_e3t    is array of vertical mesh size on T point
46!>       - dd_e3zps_min    see NEMO documentation
47!>       - dd_e3zps_rat    see NEMO documentation
48!>
49!>    to check the bathymetry in levels:<br/>
50!> @code
51!>    CALL vgrid_zgr_bat_ctl(id_mbathy, id_jpkmax, id_jpk)
52!> @endcode
53!>       - id_mbathy is array of bathymetry level
54!>       - id_jpkmax is the maximum number of level to be used
55!>       - id_jpk    is the number of level
56!>   
57!>    to compute bathy level in T,U,V,F point from  Bathymetry file:<br/>
58!> @code
59!>    tl_level(:)=vgrid_get_level(td_bathy, [cd_namelist,] [td_dom,] [id_nlevel])
60!> @endcode
61!>       - td_bathy is Bathymetry file structure
62!>       - cd_namelist is namelist [optional]
63!>       - td_dom is domain structure [optional]
64!>       - id_nlevel is number of lelvel to be used [optional]
65!>   
66!> @author
67!> J.Paul
68! REVISION HISTORY:
69!> @date November, 2013 - Initial Version
70!> @date Spetember, 2014
71!> - add header
72!> @date June, 2015 - update subroutine with NEMO 3.6
73!>
74!> @note Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
75!----------------------------------------------------------------------
76MODULE vgrid
77   USE netcdf                          ! nf90 library
78   USE kind                            ! F90 kind parameter
79   USE fct                             ! basic usefull function
80   USE global                          ! global parameter
81   USE phycst                          ! physical constant
82   USE logger                          ! log file manager
83   USE file                            ! file manager
84   USE var                             ! variable manager
85   USE dim                             ! dimension manager
86   USE dom                             ! domain manager
87   USE grid                            ! grid manager
88   USE iom                             ! I/O manager
89   USE mpp                             ! MPP manager
90   USE iom_mpp                         ! I/O MPP manager
91   IMPLICIT NONE
92   ! NOTE_avoid_public_variables_if_possible
93
94   ! type and variable
95
96   ! function and subroutine
97   PUBLIC :: vgrid_zgr_z 
98   PUBLIC :: vgrid_zgr_zps
99   PUBLIC :: vgrid_zgr_bat_ctl
100   PUBLIC :: vgrid_get_level
101
102CONTAINS
103   !-------------------------------------------------------------------
104   !> @brief This subroutine set the depth of model levels and the resulting
105   !>      vertical scale factors.
106   !
107   !> @details
108   !> ** Method  :   z-coordinate system (use in all type of coordinate)
109   !>        The depth of model levels is defined from an analytical
110   !>      function the derivative of which gives the scale factors.
111   !>        both depth and scale factors only depend on k (1d arrays). <\br>
112   !>              w-level: gdepw  = fsdep(k)                           <\br>
113   !>                       e3w(k) = dk(fsdep)(k)     = fse3(k)         <\br>
114   !>              t-level: gdept  = fsdep(k+0.5)                       <\br>
115   !>                       e3t(k) = dk(fsdep)(k+0.5) = fse3(k+0.5)     <\br>
116   !>
117   !> ** Action  : - gdept, gdepw : depth of T- and W-point (m)          <\br>
118   !>              -  e3t, e3w    : scale factors at T- and W-levels (m) <\br>
119   !>
120   !> @author G. Madec
121   !> - 03,08-  G. Madec: F90: Free form and module
122   !
123   !> @note Reference : Marti, Madec & Delecluse, 1992, JGR, 97, No8, 12,763-12,766.
124   !>
125   !> @param[inout] dd_gdepw
126   !> @param[inout] dd_gedpt
127   !> @param[inout] dd_e3w
128   !> @param[inout] dd_e2t
129   !> @param[in] dd_ppkth
130   !> @param[in] dd_ppkth2
131   !> @param[in] dd_ppacr
132   !> @param[in] dd_ppacr2
133   !> @param[in] dd_ppdzmin
134   !> @param[in] dd_pphmax
135   !> @param[in] dd_pp_to_be_computed
136   !> @param[in] dd_ppa1
137   !> @param[in] dd_ppa2
138   !> @param[in] dd_ppa0
139   !> @param[in] dd_ppsur
140   !-------------------------------------------------------------------
141   SUBROUTINE vgrid_zgr_z( dd_gdepw, dd_gdept, dd_e3w, dd_e3t,          &
142   &                       dd_e3w_1d, dd_e3t_1d, &
143   &                       dd_ppkth, dd_ppkth2, dd_ppacr, dd_ppacr2,    &
144   &                       dd_ppdzmin, dd_pphmax, dd_pp_to_be_computed, &
145   &                       dd_ppa0, dd_ppa1, dd_ppa2, dd_ppsur )
146      IMPLICIT NONE
147      ! Argument     
148      REAL(dp), DIMENSION(:), INTENT(INOUT) :: dd_gdepw
149      REAL(dp), DIMENSION(:), INTENT(INOUT) :: dd_gdept
150      REAL(dp), DIMENSION(:), INTENT(INOUT) :: dd_e3w
151      REAL(dp), DIMENSION(:), INTENT(INOUT) :: dd_e3t
152      REAL(dp), DIMENSION(:), INTENT(INOUT) :: dd_e3w_1d
153      REAL(dp), DIMENSION(:), INTENT(INOUT) :: dd_e3t_1d
154
155      REAL(dp)              , INTENT(IN   ) :: dd_ppkth
156      REAL(dp)              , INTENT(IN   ) :: dd_ppkth2
157      REAL(dp)              , INTENT(IN   ) :: dd_ppacr
158      REAL(dp)              , INTENT(IN   ) :: dd_ppacr2
159
160      REAL(dp)              , INTENT(IN   ) :: dd_ppdzmin
161      REAL(dp)              , INTENT(IN   ) :: dd_pphmax
162      REAL(dp)              , INTENT(IN   ) :: dd_pp_to_be_computed
163
164      REAL(dp)              , INTENT(IN   ) :: dd_ppa0
165      REAL(dp)              , INTENT(IN   ) :: dd_ppa1
166      REAL(dp)              , INTENT(IN   ) :: dd_ppa2
167      REAL(dp)              , INTENT(IN   ) :: dd_ppsur
168
169      ! local variable
170      REAL(dp) :: dl_zkth
171      REAL(dp) :: dl_zkth2
172      REAL(dp) :: dl_zdzmin
173      REAL(dp) :: dl_zhmax
174      REAL(dp) :: dl_zacr
175      REAL(dp) :: dl_zacr2
176
177      REAL(dp) :: dl_ppacr
178      REAL(dp) :: dl_ppacr2
179
180      REAL(dp) :: dl_za0
181      REAL(dp) :: dl_za1
182      REAL(dp) :: dl_za2
183      REAL(dp) :: dl_zsur
184      REAL(dp) :: dl_zw
185      REAL(dp) :: dl_zt
186
187      INTEGER(i4) :: il_jpk
188
189      ! loop indices
190      INTEGER(i4) :: jk
191      !----------------------------------------------------------------
192
193       dl_ppacr = dd_ppacr
194       IF( dd_ppa1 == 0._dp ) dl_ppacr =1.0
195       dl_ppacr2= dd_ppacr2
196       IF( dd_ppa2 == 0._dp ) dl_ppacr2=1.0
197
198      ! Set variables from parameters
199      ! ------------------------------
200       dl_zkth   = dd_ppkth       ;   dl_zacr  = dl_ppacr
201       dl_zdzmin = dd_ppdzmin     ;   dl_zhmax = dd_pphmax
202       dl_zkth2  = dd_ppkth2      ;   dl_zacr2 = dl_ppacr2
203
204       il_jpk = SIZE(dd_gdepw(:))
205
206      ! If ppa1 and ppa0 and ppsur are et to pp_to_be_computed
207      !  za0, za1, zsur are computed from ppdzmin , pphmax, ppkth, ppacr
208      !
209       IF(  dd_ppa1  == dd_pp_to_be_computed  .AND.  &
210         &  dd_ppa0  == dd_pp_to_be_computed  .AND.  &
211         &  dd_ppsur == dd_pp_to_be_computed           ) THEN
212         dl_za1 = ( dl_zdzmin - dl_zhmax / REAL((il_jpk-1),dp) ) &
213             &     / ( TANH((1-dl_zkth)/dl_zacr) - dl_zacr/REAL((il_jpk-1),dp) &
214             &     * (  LOG( COSH( (REAL(il_jpk,dp) - dl_zkth) / dl_zacr) )    &
215             &        - LOG( COSH(( 1.0  - dl_zkth) / dl_zacr) ) )  )
216
217         dl_za0  = dl_zdzmin - dl_za1 * TANH( (1.0-dl_zkth) / dl_zacr )
218         dl_zsur = - dl_za0 - dl_za1 * dl_zacr * LOG( COSH( (1-dl_zkth) / dl_zacr )  )
219
220       ELSE
221         dl_za1 = dd_ppa1 ;       dl_za0 = dd_ppa0 ;  dl_zsur = dd_ppsur
222         dl_za2 = dd_ppa2
223       ENDIF
224
225      ! Reference z-coordinate (depth - scale factor at T- and W-points)
226      ! ======================
227      IF(  dd_ppkth == 0. )THEN            !  uniform vertical grid       
228
229         dl_za1 = dl_zhmax/REAL((il_jpk-1),dp) 
230         DO jk = 1, il_jpk
231            dl_zw = REAL(jk,dp)
232            dl_zt = REAL(jk,dp) + 0.5_dp
233            dd_gdepw(jk) = ( dl_zw - 1.0 ) * dl_za1
234            dd_gdept(jk) = ( dl_zt - 1.0 ) * dl_za1
235            dd_e3w  (jk) =  dl_za1
236            dd_e3t  (jk) =  dl_za1
237         END DO
238
239      ELSE
240
241         DO jk = 1, il_jpk
242            dl_zw = REAL( jk,dp)
243            dl_zt = REAL( jk,dp) + 0.5_dp
244            dd_gdepw(jk) = ( dl_zsur + dl_za0 * dl_zw + &
245            &                dl_za1 * dl_zacr * LOG( COSH( (dl_zw-dl_zkth)/dl_zacr ) ) + &
246            &                dl_za2 * dl_zacr2* LOG( COSH( (dl_zw-dl_zkth2)/dl_zacr2 ) ) )
247            dd_gdept(jk) = ( dl_zsur + dl_za0 * dl_zt + &
248            &                dl_za1 * dl_zacr * LOG( COSH( (dl_zt-dl_zkth)/dl_zacr ) ) + &
249            &                dl_za2 * dl_zacr2* LOG( COSH( (dl_zt-dl_zkth2)/dl_zacr2 ) ) )
250            dd_e3w  (jk) =             dl_za0 + &
251            &                          dl_za1 * TANH(      (dl_zw-dl_zkth)/dl_zacr   ) + &
252            &                          dl_za2 * TANH(      (dl_zw-dl_zkth2)/dl_zacr2 )
253            dd_e3t  (jk) =             dl_za0 + &
254            &                          dl_za1 * TANH(      (dl_zt-dl_zkth)/dl_zacr   ) + &
255            &                          dl_za2 * TANH(      (dl_zt-dl_zkth2)/dl_zacr2 )
256         END DO
257         dd_gdepw(1) = 0.e0   ! force first w-level to be exactly at zero
258
259      ENDIF
260
261   ! need to be like this to compute the pressure gradient with ISF.
262   ! If not, level beneath the ISF are not aligned (sum(e3t) /= depth)
263   ! define e3t_0 and e3w_0 as the differences between gdept and gdepw respectively
264      DO jk = 1, il_jpk-1
265         dd_e3t_1d(jk) = dd_gdepw(jk+1)-dd_gdepw(jk) 
266      END DO
267      dd_e3t_1d(il_jpk) = dd_e3t_1d(il_jpk-1) ! we don't care because this level is masked in NEMO
268
269      DO jk = 2, il_jpk
270         dd_e3w_1d(jk) = dd_gdept(jk) - dd_gdept(jk-1) 
271      END DO
272      dd_e3w_1d(1  ) = 2._dp * (dd_gdept(1) - dd_gdepw(1))
273
274      ! Control and  print
275      ! ==================
276
277      DO jk = 1, il_jpk
278         IF( dd_e3w(jk)  <= 0. .OR. dd_e3t(jk)  <= 0. )then
279            CALL logger_debug("VGRID ZGR Z: e3w or e3t <= 0 ")
280         ENDIF   
281
282         IF( dd_e3w_1d(jk)  <= 0. .OR. dd_e3t_1d(jk)  <= 0. )then
283            CALL logger_debug("VGRID ZGR Z: e3w_1d or e3t_1d <= 0 ")
284         ENDIF
285
286         IF( dd_gdepw(jk) < 0. .OR. dd_gdept(jk) < 0. )then
287            CALL logger_debug("VGRID ZGR Z: gdepw or gdept < 0 ")
288         ENDIF
289      END DO
290
291   END SUBROUTINE vgrid_zgr_z
292   !-------------------------------------------------------------------
293   !-------------------------------------------------------------------
294   SUBROUTINE vgrid_zgr_bat( dd_bathy, dd_gdepw, dd_hmin, dd_fill )
295      IMPLICIT NONE
296      ! Argument
297      REAL(dp), DIMENSION(:,:), INTENT(INOUT) :: dd_bathy 
298      REAL(dp), DIMENSION(:)  , INTENT(IN   ) :: dd_gdepw 
299      REAL(dp)                , INTENT(IN   ) :: dd_hmin
300      REAL(dp)                , INTENT(IN   ), OPTIONAL :: dd_fill
301
302      ! local
303      INTEGER(i4) :: il_jpk
304     
305      REAL(dp)    :: dl_hmin
306      REAL(dp)    :: dl_fill
307
308      ! loop indices
309      INTEGER(i4) :: jk
310      !----------------------------------------------------------------
311      il_jpk = SIZE(dd_gdepw(:))
312
313      dl_fill=0._dp
314      IF( PRESENT(dd_fill) ) dl_fill=dd_fill
315
316      IF( dd_hmin < 0._dp ) THEN
317         jk = - INT( dd_hmin )     ! from a nb of level
318      ELSE
319         jk = MINLOC( dd_gdepw, mask = dd_gdepw > dd_hmin, dim = 1 )  ! from a depth
320      ENDIF
321     
322      dl_hmin = dd_gdepw(jk+1) ! minimum depth = ik+1 w-levels
323      WHERE( dd_bathy(:,:) <= 0._wp .OR. dd_bathy(:,:) == dl_fill )
324         dd_bathy(:,:) = dl_fill                         ! min=0     over the lands
325      ELSE WHERE
326         dd_bathy(:,:) = MAX(  dl_hmin , dd_bathy(:,:)  )   ! min=dl_hmin over the oceans
327      END WHERE
328      WRITE(*,*) 'Minimum ocean depth: ', dl_hmin, ' minimum number of ocean levels : ', jk     
329
330   END SUBROUTINE vgrid_zgr_bat
331   !-------------------------------------------------------------------
332   !> @brief This subroutine set the depth and vertical scale factor in partial step
333   !>      z-coordinate case
334   !
335   !> @details
336   !> ** Method  :   Partial steps : computes the 3D vertical scale factors
337   !>      of T-, U-, V-, W-, UW-, VW and F-points that are associated with
338   !>      a partial step representation of bottom topography.
339   !>
340   !>        The reference depth of model levels is defined from an analytical
341   !>      function the derivative of which gives the reference vertical
342   !>      scale factors.
343   !>      From  depth and scale factors reference, we compute there new value
344   !>      with partial steps  on 3d arrays ( i, j, k ).
345   !>
346   !>      w-level:
347   !>          - gdepw_ps(i,j,k)  = fsdep(k)
348   !>          - e3w_ps(i,j,k) = dk(fsdep)(k)     = fse3(i,j,k)
349   !>      t-level:
350   !>          - gdept_ps(i,j,k)  = fsdep(k+0.5)
351   !>          - e3t_ps(i,j,k) = dk(fsdep)(k+0.5) = fse3(i,j,k+0.5)
352   !>
353   !>      With the help of the bathymetric file ( bathymetry_depth_ORCA_R2.nc),
354   !>      we find the mbathy index of the depth at each grid point.
355   !>      This leads us to three cases:
356   !>          - bathy = 0 => mbathy = 0
357   !>          - 1 < mbathy < jpkm1   
358   !>          - bathy > gdepw(jpk) => mbathy = jpkm1 
359   !>
360   !>      Then, for each case, we find the new depth at t- and w- levels
361   !>      and the new vertical scale factors at t-, u-, v-, w-, uw-, vw-
362   !>      and f-points.
363   !>
364   !>        This routine is given as an example, it must be modified
365   !>      following the user s desiderata. nevertheless, the output as
366   !>      well as the way to compute the model levels and scale factors
367   !>      must be respected in order to insure second order accuracy
368   !>      schemes.
369   !>
370   !>  @warning
371   !>         - gdept, gdepw and e3 are positives
372   !>         - gdept_ps, gdepw_ps and e3_ps are positives
373   !
374   !> @author A. Bozec, G. Madec
375   !> - 02-09 (A. Bozec, G. Madec) F90: Free form and module
376   !> - 02-09 (A. de Miranda)  rigid-lid + islands
377   !>
378   !> @note Reference : Pacanowsky & Gnanadesikan 1997, Mon. Wea. Rev., 126, 3248-3270.
379   !>
380   !> @param[inout] id_mbathy
381   !> @param[inout] dd_bathy
382   !> @param[inout] id_jpkmax
383   !> @param[in] dd_gdepw
384   !> @param[in] dd_e3t
385   !> @param[in] dd_e3zps_min
386   !> @param[in] dd_e3zps_rat
387   !-------------------------------------------------------------------
388   SUBROUTINE vgrid_zgr_zps( id_mbathy, dd_bathy, id_jpkmax, &
389   &                          dd_gdepw, dd_e3t,               &
390   &                          dd_e3zps_min, dd_e3zps_rat,     &
391   &                          dd_fill )
392      IMPLICIT NONE
393      ! Argument     
394      INTEGER(i4), DIMENSION(:,:), INTENT(  OUT) :: id_mbathy
395      REAL(dp)   , DIMENSION(:,:), INTENT(INOUT) :: dd_bathy
396      INTEGER(i4)                , INTENT(INOUT) :: id_jpkmax
397      REAL(dp)   , DIMENSION(:)  , INTENT(IN   ) :: dd_gdepw
398      REAL(dp)   , DIMENSION(:)  , INTENT(IN   ) :: dd_e3t
399      REAL(dp)                   , INTENT(IN   ) :: dd_e3zps_min
400      REAL(dp)                   , INTENT(IN   ) :: dd_e3zps_rat
401      REAL(dp)                   , INTENT(IN   ), OPTIONAL :: dd_fill
402
403      ! local variable
404      REAL(dp) :: dl_zmax     ! Maximum depth
405      !REAL(dp) :: dl_zmin     ! Minimum depth
406      REAL(dp) :: dl_zdepth   ! Ajusted ocean depth to avoid too small e3t
407      REAL(dp) :: dl_fill     
408
409      INTEGER(i4) :: il_jpk
410      INTEGER(i4) :: il_jpkm1
411      INTEGER(i4) :: il_jpiglo
412      INTEGER(i4) :: il_jpjglo
413
414      ! loop indices
415      INTEGER(i4) :: ji
416      INTEGER(i4) :: jj
417      INTEGER(i4) :: jk
418      !----------------------------------------------------------------
419
420      il_jpk=SIZE(dd_gdepw(:))
421      il_jpiglo=SIZE(id_mbathy(:,:),DIM=1)
422      il_jpjglo=SIZE(id_mbathy(:,:),DIM=2)
423
424      dl_fill=0._dp
425      IF( PRESENT(dd_fill) ) dl_fill=dd_fill
426
427      ! Initialization of constant
428      dl_zmax = dd_gdepw(il_jpk) + dd_e3t(il_jpk) ! maximum depth (i.e. the last ocean level thickness <= 2*e3t_1d(jpkm1) )
429
430      ! bounded value of bathy (min already set at the end of zgr_bat)
431      WHERE( dd_bathy(:,:) /= dl_fill )
432         dd_bathy(:,:) = MIN( dl_zmax ,  dd_bathy(:,:) )
433      END WHERE
434
435      ! bathymetry in level (from bathy_meter)
436      ! ===================
437      il_jpkm1=il_jpk-1
438      ! initialize mbathy to the maximum ocean level available
439      id_mbathy(:,:) = il_jpkm1
440
441      ! storage of land and island's number (zera and negative values) in mbathy
442      DO jj = 1, il_jpjglo
443         DO ji= 1, il_jpiglo
444            IF( dd_bathy(ji,jj) <= 0._dp )THEN
445               id_mbathy(ji,jj) = INT(dd_bathy(ji,jj),i4)
446            ELSEIF( dd_bathy(ji,jj) == dl_fill )THEN
447               id_mbathy(ji,jj) = 0_i4
448            ENDIF
449         END DO
450      END DO
451
452      ! Compute mbathy for ocean points (i.e. the number of ocean levels)
453      ! find the number of ocean levels such that the last level thickness
454      ! is larger than the minimum of e3zps_min and e3zps_rat * e3t (where
455      ! e3t is the reference level thickness
456
457      DO jk = il_jpkm1, 1, -1
458         dl_zdepth = dd_gdepw(jk) + MIN( dd_e3zps_min, dd_e3t(jk)*dd_e3zps_rat )
459
460         DO jj = 1, il_jpjglo
461            DO ji = 1, il_jpiglo
462               IF( dd_bathy(ji,jj) /= dl_fill )THEN
463                  IF( 0. < dd_bathy(ji,jj) .AND. &
464                  &       dd_bathy(ji,jj) <= dl_zdepth ) id_mbathy(ji,jj) = jk-1
465               ENDIF
466            END DO
467         END DO
468      END DO
469
470      ! ================
471      ! Bathymetry check
472      ! ================
473
474      CALL vgrid_zgr_bat_ctl( id_mbathy, id_jpkmax, il_jpk)
475
476   END SUBROUTINE vgrid_zgr_zps
477   !-------------------------------------------------------------------
478   !> @brief This subroutine check the bathymetry in levels
479   !
480   !> @details
481   !> ** Method  :   The array mbathy is checked to verified its consistency
482   !>      with the model options. in particular:
483   !>            mbathy must have at least 1 land grid-points (mbathy<=0)
484   !>                  along closed boundary.
485   !>            mbathy must be cyclic IF jperio=1.
486   !>            mbathy must be lower or equal to jpk-1.
487   !>            isolated ocean grid points are suppressed from mbathy
488   !>                  since they are only connected to remaining
489   !>                  ocean through vertical diffusion.
490   !>      C A U T I O N : mbathy will be modified during the initializa-
491   !>      tion phase to become the number of non-zero w-levels of a water
492   !>      column, with a minimum value of 1.
493   !>
494   !> ** Action  : - update mbathy: level bathymetry (in level index)
495   !>              - update bathy : meter bathymetry (in meters)
496
497   !> @author G.Madec
498   !> - 03-08 Original code
499   !
500   !> @param[in] id_mbathy
501   !> @param[in] id_jpkmax
502   !> @param[in] id_jpk
503   !-------------------------------------------------------------------
504   SUBROUTINE vgrid_zgr_bat_ctl( id_mbathy, id_jpkmax, id_jpk)
505      IMPLICIT NONE
506      ! Argument     
507      INTEGER(i4), DIMENSION(:,:), INTENT(INOUT) :: id_mbathy
508      INTEGER(i4)                , INTENT(INOUT) :: id_jpkmax
509      INTEGER(i4)                , INTENT(INOUT) :: id_jpk
510
511      ! local variable
512      INTEGER(i4) :: il_jpiglo
513      INTEGER(i4) :: il_jpjglo
514
515      INTEGER(i4) :: il_icompt
516      INTEGER(i4) :: il_ibtest
517      INTEGER(i4) :: il_ikmax
518      INTEGER(i4) :: il_jpkm1
519
520      INTEGER(i4) :: il_jim
521      INTEGER(i4) :: il_jip
522      INTEGER(i4) :: il_jjm
523      INTEGER(i4) :: il_jjp
524
525      ! loop indices
526      INTEGER(i4) :: jl
527      INTEGER(i4) :: ji
528      INTEGER(i4) :: jj
529      !----------------------------------------------------------------
530
531      il_jpiglo=SIZE(id_mbathy(:,:),DIM=1)
532      il_jpjglo=SIZE(id_mbathy(:,:),DIM=2)
533
534      ! ================
535      ! Bathymetry check
536      ! ================
537
538      ! suppress isolated ocean grid points'
539      il_icompt = 0
540
541      DO jl = 1, 2
542         DO jj = 1, il_jpjglo
543            DO ji = 1, il_jpiglo
544               il_jim=max(ji-1,1) ; il_jip=min(ji+1,il_jpiglo)
545               il_jjm=max(jj-1,1) ; il_jjp=min(jj+1,il_jpjglo)
546
547               if(il_jim==ji) il_jim=il_jip ; if(il_jip==ji) il_jip=il_jim
548               if(il_jjm==jj) il_jjm=il_jjp ; if(il_jjp==jj) il_jjp=il_jjm
549
550               il_ibtest = MAX( id_mbathy(il_jim,jj), id_mbathy(il_jip,jj),  &
551                                id_mbathy(ji,il_jjm),id_mbathy(ji,il_jjp) )
552
553               IF( il_ibtest < id_mbathy(ji,jj) ) THEN
554                  id_mbathy(ji,jj) = il_ibtest
555                  il_icompt = il_icompt + 1
556               ENDIF
557            END DO
558         END DO
559
560      END DO
561      IF( il_icompt == 0 ) THEN
562         CALL logger_info("VGRID ZGR BAT CTL:  no isolated ocean grid points")
563      ELSE
564         CALL logger_info("VGRID ZGR BAT CTL:"//TRIM(fct_str(il_icompt))//&
565         &              " ocean grid points suppressed")
566      ENDIF
567
568      id_mbathy(:,:) = MAX( 0, id_mbathy(:,:))
569
570      ! Number of ocean level inferior or equal to jpkm1
571
572      il_ikmax = 0
573      DO jj = 1, il_jpjglo
574         DO ji = 1, il_jpiglo
575            il_ikmax = MAX( il_ikmax, id_mbathy(ji,jj) )
576         END DO
577      END DO
578
579      id_jpkmax=id_jpk
580
581      il_jpkm1=id_jpk-1
582      IF( il_ikmax > il_jpkm1 ) THEN
583         CALL logger_error("VGRID ZGR BAT CTL: maximum number of ocean level = "//&
584         &              TRIM(fct_str(il_ikmax))//" >  jpk-1."//&
585         &              " Change jpk to "//TRIM(fct_str(il_ikmax+1))//&
586         &              " to use the exact ead bathymetry" )
587      ELSE IF( il_ikmax < il_jpkm1 ) THEN
588         id_jpkmax=il_ikmax+1
589      ENDIF
590
591   END SUBROUTINE vgrid_zgr_bat_ctl
592   !-------------------------------------------------------------------
593   !> @brief This function compute bathy level in T,U,V,F point, and return
594   !> them as array of variable structure
595   !
596   !> @details
597   !> Bathymetry is read on Bathymetry file, then bathy level is computed
598   !> on T point, and finally fit to U,V,F point.
599   !>
600   !> you could specify :<br/>
601   !> - namelist where find parameter to set the depth of model levels
602   !> (default use GLORYS 75 levels parameters)
603   !> - domain structure to specify on e area to work on
604   !> - number of level to be used
605   !>
606   !> @author J.Paul
607   !> - November, 2013- Initial Version
608   !
609   !> @param[in] td_bathy     Bathymetry file structure
610   !> @param[in] cd_namelist  namelist
611   !> @param[in] td_dom       domain structure
612   !> @param[in] id_nlevel    number of lelvel to be used
613   !> @return array of level on T,U,V,F point (variable structure)
614   !-------------------------------------------------------------------
615   FUNCTION vgrid_get_level(td_bathy, cd_namelist, td_dom, id_nlevel)
616      IMPLICIT NONE
617      ! Argument
618      TYPE(TMPP)      , INTENT(IN) :: td_bathy
619      CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: cd_namelist
620      TYPE(TDOM)      , INTENT(IN), OPTIONAL :: td_dom
621      INTEGER(i4)     , INTENT(IN), OPTIONAL :: id_nlevel
622
623      ! function
624      TYPE(TVAR), DIMENSION(ip_npoint) :: vgrid_get_level
625
626      ! local variable
627      REAL(dp)   , DIMENSION(:)      , ALLOCATABLE :: dl_gdepw 
628      REAL(dp)   , DIMENSION(:)      , ALLOCATABLE :: dl_gdept 
629      REAL(dp)   , DIMENSION(:)      , ALLOCATABLE :: dl_e3w 
630      REAL(dp)   , DIMENSION(:)      , ALLOCATABLE :: dl_e3t
631      REAL(dp)   , DIMENSION(:)      , ALLOCATABLE :: dl_e3w_1d 
632      REAL(dp)   , DIMENSION(:)      , ALLOCATABLE :: dl_e3t_1d
633
634      INTEGER(i4)                                  :: il_status
635      INTEGER(i4)                                  :: il_fileid
636      INTEGER(i4)                                  :: il_jpkmax
637      INTEGER(i4), DIMENSION(2,2)                  :: il_xghost
638      INTEGER(i4), DIMENSION(:,:)    , ALLOCATABLE :: il_mbathy
639      INTEGER(i4), DIMENSION(:,:,:,:), ALLOCATABLE :: il_level
640     
641      LOGICAL                                      :: ll_exist
642
643      TYPE(TDIM) , DIMENSION(ip_maxdim)            :: tl_dim
644
645      TYPE(TDOM)                                   :: tl_dom
646
647      TYPE(TVAR)                                   :: tl_var
648
649      TYPE(TMPP)                                   :: tl_bathy
650
651      ! loop indices
652      INTEGER(i4) :: ji
653      INTEGER(i4) :: jj
654
655      INTEGER(i4) :: jip
656      INTEGER(i4) :: jjp
657
658      !namelist (intialise with GLORYS 75 levels parameters)
659      REAL(dp)                                :: dn_pp_to_be_computed = 0._dp
660      REAL(dp)                                :: dn_ppsur     = -3958.951371276829_dp
661      REAL(dp)                                :: dn_ppa0      =   103.9530096000000_dp
662      REAL(dp)                                :: dn_ppa1      =     2.4159512690000_dp
663      REAL(dp)                                :: dn_ppa2      =   100.7609285000000_dp
664      REAL(dp)                                :: dn_ppkth     =    15.3510137000000_dp
665      REAL(dp)                                :: dn_ppkth2    =    48.0298937200000_dp
666      REAL(dp)                                :: dn_ppacr     =     7.0000000000000_dp
667      REAL(dp)                                :: dn_ppacr2    =    13.000000000000_dp
668      REAL(dp)                                :: dn_ppdzmin   = 6._dp
669      REAL(dp)                                :: dn_pphmax    = 5750._dp
670      INTEGER(i4)                             :: in_nlevel    = 75
671
672      REAL(dp)                                :: dn_e3zps_min = 25._dp
673      REAL(dp)                                :: dn_e3zps_rat = 0.2_dp
674      !----------------------------------------------------------------
675      NAMELIST /namzgr/ &
676      &  dn_pp_to_be_computed, &
677      &  dn_ppsur,     &
678      &  dn_ppa0,      &
679      &  dn_ppa1,      &
680      &  dn_ppa2,      &
681      &  dn_ppkth,     &
682      &  dn_ppkth2,    &
683      &  dn_ppacr,     &
684      &  dn_ppacr2,    &
685      &  dn_ppdzmin,   &
686      &  dn_pphmax,    &
687      &  in_nlevel
688
689      NAMELIST /namzps/ &
690      &  dn_e3zps_min, &
691      &  dn_e3zps_rat
692      !----------------------------------------------------------------
693
694      IF( PRESENT(cd_namelist) )THEN
695         !1- read namelist
696         INQUIRE(FILE=TRIM(cd_namelist), EXIST=ll_exist)
697         IF( ll_exist )THEN
698 
699            il_fileid=fct_getunit()
700
701            OPEN( il_fileid, FILE=TRIM(cd_namelist), &
702            &                FORM='FORMATTED',       &
703            &                ACCESS='SEQUENTIAL',    &
704            &                STATUS='OLD',           &
705            &                ACTION='READ',          &
706            &                IOSTAT=il_status)
707            CALL fct_err(il_status)
708            IF( il_status /= 0 )THEN
709               CALL logger_fatal("VGRID GET LEVEL: ERROR opening "//&
710               &  TRIM(cd_namelist))
711            ENDIF
712
713            READ( il_fileid, NML = namzgr )
714            READ( il_fileid, NML = namzps )
715
716            CLOSE( il_fileid, IOSTAT=il_status )
717            CALL fct_err(il_status)
718            IF( il_status /= 0 )THEN
719               CALL logger_error("VGRID GET LEVELL: ERROR closing "//&
720               &  TRIM(cd_namelist))
721            ENDIF
722
723         ELSE
724
725            CALL logger_fatal("VGRID GET LEVEL: ERROR. can not find "//&
726            &  TRIM(cd_namelist))
727
728         ENDIF
729      ENDIF
730
731      ! copy structure
732      tl_bathy=mpp_copy(td_bathy)
733
734      ! get domain
735      IF( PRESENT(td_dom) )THEN
736         tl_dom=dom_copy(td_dom)
737      ELSE
738         CALL logger_debug("VGRID GET LEVEL: get dom from "//&
739         &  TRIM(tl_bathy%c_name))
740         tl_dom=dom_init(tl_bathy)
741      ENDIF
742
743      ! get ghoste cell
744      il_xghost(:,:)=grid_get_ghost(tl_bathy)
745
746      ! open mpp files
747      CALL iom_dom_open(tl_bathy, tl_dom)
748
749      ! check namelist
750      IF( PRESENT(id_nlevel) ) in_nlevel=id_nlevel
751      IF( in_nlevel == 0 )THEN
752         CALL logger_fatal("VGRID GET LEVEL: number of level to be used "//&
753         &                 "is not specify. check namelist.")
754      ENDIF
755
756      ! read bathymetry
757      tl_var=iom_dom_read_var(tl_bathy,'bathymetry',tl_dom)
758      ! clean
759      CALL dom_clean(tl_dom)
760
761      ! remove ghost cell
762      CALL grid_del_ghost(tl_var, il_xghost(:,:))
763
764      ! force _FillValue (land) to be 0
765      WHERE( tl_var%d_value(:,:,1,1) == tl_var%d_fill )
766         tl_var%d_value(:,:,1,1)=0
767      END WHERE
768
769      ! clean
770      CALL iom_dom_close(tl_bathy)
771      CALL mpp_clean(tl_bathy)
772
773      ! compute vertical grid
774      ALLOCATE( dl_gdepw(in_nlevel), dl_gdept(in_nlevel) ) 
775      ALLOCATE(   dl_e3w(in_nlevel),   dl_e3t(in_nlevel) ) 
776      ALLOCATE(   dl_e3w_1d(in_nlevel),   dl_e3t_1d(in_nlevel) ) 
777      CALL vgrid_zgr_z( dl_gdepw(:), dl_gdept(:), dl_e3w(:), dl_e3t(:), &
778      &                 dl_e3w_1d, dl_e3t_1d, &
779      &                 dn_ppkth, dn_ppkth2, dn_ppacr, dn_ppacr2,       &
780      &                 dn_ppdzmin, dn_pphmax, dn_pp_to_be_computed,    &
781      &                 dn_ppa0, dn_ppa1, dn_ppa2, dn_ppsur )
782
783      ! compute bathy level on T point
784      ALLOCATE( il_mbathy(tl_var%t_dim(1)%i_len, &
785      &                   tl_var%t_dim(2)%i_len ) )
786      CALL vgrid_zgr_zps( il_mbathy(:,:), tl_var%d_value(:,:,1,1), il_jpkmax, &
787      &                   dl_gdepw(:), dl_e3t(:),               &
788      &                   dn_e3zps_min, dn_e3zps_rat )
789
790      DEALLOCATE( dl_gdepw, dl_gdept ) 
791      DEALLOCATE(   dl_e3w,   dl_e3t )
792
793      ! compute bathy level in T,U,V,F point
794      ALLOCATE( il_level(tl_var%t_dim(1)%i_len, &
795      &                  tl_var%t_dim(2)%i_len, &
796      &                  ip_npoint,1) )
797
798      DO jj=1,tl_var%t_dim(2)%i_len
799         DO ji= 1,tl_var%t_dim(1)%i_len
800
801         jip=MIN(ji+1,tl_var%t_dim(1)%i_len)
802         jjp=MIN(jj+1,tl_var%t_dim(2)%i_len)
803
804         ! T point
805         il_level(ji,jj,jp_T,1)=il_mbathy(ji,jj)
806         ! U point
807         il_level(ji,jj,jp_U,1)=MIN( il_mbathy(ji, jj ), il_mbathy(jip, jj ))
808         ! V point
809         il_level(ji,jj,jp_V,1)=MIN( il_mbathy(ji, jj ), il_mbathy(ji , jjp))
810         ! F point
811         il_level(ji,jj,jp_F,1)=MIN( il_mbathy(ji, jj ), il_mbathy(jip, jj ), &
812         &                           il_mbathy(ji, jjp), il_mbathy(jip, jjp))
813
814         ENDDO
815      ENDDO
816
817      DEALLOCATE( il_mbathy )
818
819      tl_dim(:)=dim_copy(tl_var%t_dim(:))
820      ! clean
821      CALL var_clean(tl_var)
822
823      ! only 2 first dimension to be used
824      tl_dim(3:4)%l_use=.FALSE.
825
826      vgrid_get_level(jp_T)=var_init( 'tlevel', il_level(:,:,jp_T:jp_T,:), &
827      &                                td_dim=tl_dim(:) )
828      vgrid_get_level(jp_U)=var_init( 'ulevel', il_level(:,:,jp_U:jp_U,:), &
829      &                                td_dim=tl_dim(:))
830      vgrid_get_level(jp_V)=var_init( 'vlevel', il_level(:,:,jp_V:jp_V,:), &
831      &                                td_dim=tl_dim(:))
832      vgrid_get_level(jp_F)=var_init( 'flevel', il_level(:,:,jp_F:jp_F,:), &
833      &                                td_dim=tl_dim(:))
834
835      DEALLOCATE( il_level )
836
837      CALL grid_add_ghost( vgrid_get_level(jp_T), il_xghost(:,:) )
838      CALL grid_add_ghost( vgrid_get_level(jp_U), il_xghost(:,:) )
839      CALL grid_add_ghost( vgrid_get_level(jp_V), il_xghost(:,:) )
840      CALL grid_add_ghost( vgrid_get_level(jp_F), il_xghost(:,:) )
841
842      ! clean
843      CALL dim_clean(tl_dim(:))
844
845   END FUNCTION vgrid_get_level
846END MODULE vgrid
847
Note: See TracBrowser for help on using the repository browser.