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.
trabbl.F90 in tags/nemo_v3_2/nemo_v3_2/NEMO/OPA_SRC/TRA – NEMO

source: tags/nemo_v3_2/nemo_v3_2/NEMO/OPA_SRC/TRA/trabbl.F90 @ 1878

Last change on this file since 1878 was 1878, checked in by flavoni, 14 years ago

initial test for nemogcm

File size: 21.6 KB
Line 
1MODULE trabbl
2   !!==============================================================================
3   !!                       ***  MODULE  trabbl  ***
4   !! Ocean physics :  advective and/or diffusive bottom boundary layer scheme
5   !!==============================================================================
6   !! History :  8.0  !  96-06  (L. Mortier)  Original code
7   !!            8.0  !  97-11  (G. Madec)  Optimization
8   !!            8.5  !  02-08  (G. Madec)  free form + modules
9   !!----------------------------------------------------------------------
10#if   defined key_trabbl_dif   ||   defined key_trabbl_adv   || defined key_esopa
11   !!----------------------------------------------------------------------
12   !!   'key_trabbl_dif'   or            diffusive bottom boundary layer
13   !!   'key_trabbl_adv'                 advective bottom boundary layer
14   !!----------------------------------------------------------------------
15   !!----------------------------------------------------------------------
16   !!   tra_bbl_dif  : update the active tracer trends due to the bottom
17   !!                  boundary layer (diffusive only)
18   !!   tra_bbl_adv  : update the active tracer trends due to the bottom
19   !!                  boundary layer (advective and/or diffusive)
20   !!   tra_bbl_init : initialization, namlist read, parameters control
21   !!----------------------------------------------------------------------
22   USE oce                ! ocean dynamics and active tracers
23   USE dom_oce            ! ocean space and time domain
24   USE trdmod             ! ocean active tracers trends
25   USE trdmod_oce         ! ocean variables trends
26   USE in_out_manager     ! I/O manager
27   USE lbclnk             ! ocean lateral boundary conditions
28   USE prtctl             ! Print control
29
30   IMPLICIT NONE
31   PRIVATE
32
33   PUBLIC tra_bbl_dif    ! routine called by step.F90
34   PUBLIC tra_bbl_adv    ! routine called by step.F90
35
36   !!* Namelist nambbl: bottom boundary layer
37   REAL(wp), PUBLIC ::   rn_ahtbbl = 1.e+3   !: lateral coeff. for bottom boundary layer scheme (m2/s)
38
39# if defined key_trabbl_dif
40   LOGICAL, PUBLIC, PARAMETER ::   lk_trabbl_dif = .TRUE.          !: diffusive bottom boundary layer flag
41# else
42   LOGICAL, PUBLIC, PARAMETER ::   lk_trabbl_dif = .FALSE.         !: diffusive bottom boundary layer flag
43# endif
44
45# if defined key_trabbl_adv
46   LOGICAL, PUBLIC, PARAMETER ::   lk_trabbl_adv = .TRUE.   !: advective bottom boundary layer flag
47   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   u_bbl      !: 3 components of the velocity
48   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   v_bbl      !: associated with advective BBL
49   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   w_bbl      !: (only affect tracer)
50# else
51   LOGICAL, PUBLIC, PARAMETER ::   lk_trabbl_adv = .FALSE.  !: advective bottom boundary layer flag
52# endif
53
54   INTEGER, DIMENSION(jpi,jpj) ::   mbkt          ! vertical index of the bottom ocean T-level
55   INTEGER, DIMENSION(jpi,jpj) ::   mbku, mbkv    ! vertical index of the bottom ocean U/V-level
56
57   !! * Substitutions
58#  include "domzgr_substitute.h90"
59#  include "vectopt_loop_substitute.h90"
60   !!----------------------------------------------------------------------
61   !!   OPA 9.0 , LOCEAN-IPSL (2006)
62   !! $Id: trabbl.F90 1601 2009-08-11 10:09:19Z ctlod $
63   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
64   !!----------------------------------------------------------------------
65
66CONTAINS
67
68   SUBROUTINE tra_bbl_dif( kt )
69      !!----------------------------------------------------------------------
70      !!                   ***  ROUTINE tra_bbl_dif  ***
71      !!
72      !! ** Purpose :   Compute the before tracer (t & s) trend associated
73      !!      with the bottom boundary layer and add it to the general trend
74      !!      of tracer equations. The bottom boundary layer is supposed to be
75      !!      a purely diffusive bottom boundary layer.
76      !!
77      !! ** Method  :   When the product grad( rho) * grad(h) < 0 (where grad
78      !!      is an along bottom slope gradient) an additional lateral diffu-
79      !!      sive trend along the bottom slope is added to the general tracer
80      !!      trend, otherwise nothing is done.
81      !!      Second order operator (laplacian type) with variable coefficient
82      !!      computed as follow for temperature (idem on s):
83      !!         difft = 1/(e1t*e2t*e3t) { di-1[ ahbt e2u*e3u/e1u di[ztb] ]
84      !!                                 + dj-1[ ahbt e1v*e3v/e2v dj[ztb] ] }
85      !!      where ztb is a 2D array: the bottom ocean temperature and ahtb
86      !!      is a time and space varying diffusive coefficient defined by:
87      !!         ahbt = zahbp    if grad(rho).grad(h) < 0
88      !!              = 0.       otherwise.
89      !!      Note that grad(.) is the along bottom slope gradient. grad(rho)
90      !!      is evaluated using the local density (i.e. referenced at the
91      !!      local depth). Typical value of ahbt is 2000 m2/s (equivalent to
92      !!      a downslope velocity of 20 cm/s if the condition for slope
93      !!      convection is satified)
94      !!      Add this before trend to the general trend (ta,sa) of the
95      !!      botton ocean tracer point:
96      !!         ta = ta + difft
97      !!
98      !! ** Action  : - update (ta,sa) at the bottom level with the bottom
99      !!                boundary layer trend
100      !!              - save the trends in ztrdt/ztrds ('key_trdtra')
101      !!
102      !! References : Beckmann, A., and R. Doscher, 1997, J. Phys.Oceanogr., 581-591.
103      !!----------------------------------------------------------------------
104      USE oce, ONLY :   ztrdt => ua   ! use ua as 3D workspace   
105      USE oce, ONLY :   ztrds => va   ! use va as 3D workspace   
106      USE eosbn2                      ! equation of state
107      !!
108      INTEGER, INTENT( in ) ::   kt   ! ocean time-step
109      !!
110      INTEGER  ::   ji, jj                  ! dummy loop indices
111      INTEGER  ::   ik
112      INTEGER  ::   ii0, ii1, ij0, ij1      ! temporary integers
113      INTEGER  ::   iku1, iku2, ikv1,ikv2   ! temporary intergers
114      REAL(wp) ::   ze3u, ze3v              ! temporary scalars
115      INTEGER  ::   iku, ikv
116      REAL(wp) ::   zsign, zt, zs, zh, zalbet   ! temporary scalars
117      REAL(wp) ::   zgdrho, zbtr, zta, zsa
118      REAL(wp), DIMENSION(jpi,jpj) ::   zki, zkj, zkw, zkx, zky, zkz             ! 2D workspace
119      REAL(wp), DIMENSION(jpi,jpj) ::   ztnb, zsnb, zdep, ztbb, zsbb, zahu, zahv
120      !!
121      REAL(wp) ::    fsalbt, pft, pfs, pfh   ! statement function
122      !!----------------------------------------------------------------------
123      ! ratio alpha/beta
124      ! ================
125      !  fsalbt: ratio of thermal over saline expension coefficients
126      !       pft :  potential temperature in degrees celcius
127      !       pfs :  salinity anomaly (s-35) in psu
128      !       pfh :  depth in meters
129
130      fsalbt( pft, pfs, pfh ) =                                              &
131         ( ( ( -0.255019e-07 * pft + 0.298357e-05 ) * pft                    &
132         &                         - 0.203814e-03 ) * pft                    &
133         &                         + 0.170907e-01 ) * pft                    &
134         &                         + 0.665157e-01                            &
135         +(-0.678662e-05 * pfs - 0.846960e-04 * pft + 0.378110e-02 ) * pfs   &
136         +  ( ( - 0.302285e-13 * pfh                                         &
137         &      - 0.251520e-11 * pfs                                         &
138         &      + 0.512857e-12 * pft * pft          ) * pfh                  &
139         &                           - 0.164759e-06   * pfs                  &
140         &   +(   0.791325e-08 * pft - 0.933746e-06 ) * pft                  &
141         &                           + 0.380374e-04 ) * pfh   
142      !!----------------------------------------------------------------------
143
144      IF( kt == nit000 )   CALL tra_bbl_init
145
146      IF( l_trdtra )   THEN         ! Save ta and sa trends
147         ztrdt(:,:,:) = ta(:,:,:) 
148         ztrds(:,:,:) = sa(:,:,:) 
149      ENDIF
150
151      ! 0. 2D fields of bottom temperature and salinity, and bottom slope
152      ! -----------------------------------------------------------------
153      ! mbathy= number of w-level, minimum value=1 (cf dommsk.F)
154#  if defined key_vectopt_loop
155      DO jj = 1, 1
156         DO ji = 1, jpij   ! vector opt. (forced unrolling)
157#  else
158      DO jj = 1, jpj
159         DO ji = 1, jpi
160#  endif
161            ik = mbkt(ji,jj)                              ! index of the bottom ocean T-level
162            ztnb(ji,jj) = tn(ji,jj,ik) * tmask(ji,jj,1)   ! masked now T and S at ocean bottom
163            zsnb(ji,jj) = sn(ji,jj,ik) * tmask(ji,jj,1)
164            ztbb(ji,jj) = tb(ji,jj,ik) * tmask(ji,jj,1)   ! masked before T and S at ocean bottom
165            zsbb(ji,jj) = sb(ji,jj,ik) * tmask(ji,jj,1)
166            zdep(ji,jj) = fsdept(ji,jj,ik)                ! depth of the ocean bottom T-level
167         END DO
168      END DO
169
170      IF( ln_zps ) THEN      ! partial steps correction
171# if defined key_vectopt_loop
172         DO jj = 1, 1
173            DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
174# else
175         DO jj = 1, jpjm1
176            DO ji = 1, jpim1
177# endif
178               iku1 = MAX( mbathy(ji+1,jj  )-1, 1 )
179               iku2 = MAX( mbathy(ji  ,jj  )-1, 1 )
180               ikv1 = MAX( mbathy(ji  ,jj+1)-1, 1 )
181               ikv2 = MAX( mbathy(ji  ,jj  )-1, 1 )
182               ze3u = MIN( fse3u(ji,jj,iku1), fse3u(ji,jj,iku2) ) 
183               ze3v = MIN( fse3v(ji,jj,ikv1), fse3v(ji,jj,ikv2) ) 
184               zahu(ji,jj) = rn_ahtbbl * e2u(ji,jj) * ze3u / e1u(ji,jj) * umask(ji,jj,1)
185               zahv(ji,jj) = rn_ahtbbl * e1v(ji,jj) * ze3v / e2v(ji,jj) * vmask(ji,jj,1)
186            END DO
187         END DO
188      ELSE                    ! z-coordinate - full steps or s-coordinate
189#   if defined key_vectopt_loop
190         DO jj = 1, 1
191            DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
192#   else
193         DO jj = 1, jpjm1
194            DO ji = 1, jpim1
195#   endif
196               iku = mbku(ji,jj)
197               ikv = mbkv(ji,jj)
198               zahu(ji,jj) = rn_ahtbbl * e2u(ji,jj) * fse3u(ji,jj,iku) / e1u(ji,jj) * umask(ji,jj,1)
199               zahv(ji,jj) = rn_ahtbbl * e1v(ji,jj) * fse3v(ji,jj,ikv) / e2v(ji,jj) * vmask(ji,jj,1)
200            END DO
201         END DO
202      ENDIF
203
204      ! 1. Criteria of additional bottom diffusivity: grad(rho).grad(h)<0
205      ! --------------------------------------------
206      ! Sign of the local density gradient along the i- and j-slopes
207      ! multiplied by the slope of the ocean bottom
208
209      SELECT CASE ( nn_eos )
210      !
211      CASE ( 0 )                 !==  Jackett and McDougall (1994) formulation  ==!
212#  if defined key_vectopt_loop
213         DO jj = 1, 1
214            DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
215#  else
216         DO jj = 1, jpjm1
217            DO ji = 1, jpim1
218#  endif
219               ! temperature, salinity anomalie and depth
220               zt = 0.5 * ( ztnb(ji,jj) + ztnb(ji+1,jj) )
221               zs = 0.5 * ( zsnb(ji,jj) + zsnb(ji+1,jj) ) - 35.0
222               zh = 0.5 * ( zdep(ji,jj) + zdep(ji+1,jj) )
223               ! masked ratio alpha/beta
224               zalbet = fsalbt( zt, zs, zh )*umask(ji,jj,1)
225               ! local density gradient along i-bathymetric slope
226               zgdrho = zalbet * ( ztnb(ji+1,jj) - ztnb(ji,jj) )   &
227                      -          ( zsnb(ji+1,jj) - zsnb(ji,jj) )
228               ! sign of local i-gradient of density multiplied by the i-slope
229               zsign = SIGN( 0.5, - zgdrho * ( zdep(ji+1,jj) - zdep(ji,jj) ) )
230               zki(ji,jj) = ( 0.5 - zsign ) * zahu(ji,jj)
231               !
232               ! temperature, salinity anomalie and depth
233               zt = 0.5 * ( ztnb(ji,jj+1) + ztnb(ji,jj) )
234               zs = 0.5 * ( zsnb(ji,jj+1) + zsnb(ji,jj) ) - 35.0
235               zh = 0.5 * ( zdep(ji,jj+1) + zdep(ji,jj) )
236               ! masked ratio alpha/beta
237               zalbet = fsalbt( zt, zs, zh )*vmask(ji,jj,1)
238               ! local density gradient along j-bathymetric slope
239               zgdrho = zalbet * ( ztnb(ji,jj+1) - ztnb(ji,jj) )   &
240                      -          ( zsnb(ji,jj+1) - zsnb(ji,jj) )
241               ! sign of local j-gradient of density multiplied by the j-slope
242               zsign = sign( 0.5, -zgdrho * ( zdep(ji,jj+1) - zdep(ji,jj) ) )
243               zkj(ji,jj) = ( 0.5 - zsign ) * zahv(ji,jj)
244            END DO
245         END DO
246         !
247      CASE ( 1 )               !==  Linear formulation function of temperature only  ==!
248#  if defined key_vectopt_loop
249         DO jj = 1, 1
250            DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
251#  else
252         DO jj = 1, jpjm1
253            DO ji = 1, jpim1
254#  endif
255               ! local 'density/temperature' gradient along i-bathymetric slope
256               zgdrho =  ztnb(ji+1,jj) - ztnb(ji,jj) 
257               ! sign of local i-gradient of density multiplied by the i-slope
258               zsign = SIGN( 0.5, - zgdrho * ( zdep(ji+1,jj) - zdep(ji,jj) ) )
259               zki(ji,jj) = ( 0.5 - zsign ) * zahu(ji,jj)
260               !
261               ! local density gradient along j-bathymetric slope
262               zgdrho =  ztnb(ji,jj+1) - ztnb(ji,jj) 
263               ! sign of local j-gradient of density multiplied by the j-slope
264               zsign = sign( 0.5, -zgdrho * ( zdep(ji,jj+1) - zdep(ji,jj) ) )
265               zkj(ji,jj) = ( 0.5 - zsign ) * zahv(ji,jj)
266            END DO
267         END DO
268         !
269      CASE ( 2 )               !==  Linear formulation function of temperature and salinity  ==!
270#  if defined key_vectopt_loop
271         DO jj = 1, 1
272            DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
273#  else
274         DO jj = 1, jpjm1
275            DO ji = 1, jpim1
276#  endif     
277               ! local density gradient along i-bathymetric slope
278               zgdrho = - ( rn_beta *( zsnb(ji+1,jj) - zsnb(ji,jj) )   &
279                  &       - rn_alpha*( ztnb(ji+1,jj) - ztnb(ji,jj) ) )
280               ! sign of local i-gradient of density multiplied by the i-slope
281               zsign = SIGN( 0.5, - zgdrho * ( zdep(ji+1,jj) - zdep(ji,jj) ) )
282               zki(ji,jj) = ( 0.5 - zsign ) * zahu(ji,jj)
283               !
284               ! local density gradient along j-bathymetric slope
285               zgdrho = - ( rn_beta *( zsnb(ji,jj+1) - zsnb(ji,jj) )   &
286                  &       - rn_alpha*( ztnb(ji,jj+1) - ztnb(ji,jj) ) )   
287               ! sign of local j-gradient of density multiplied by the j-slope
288               zsign = sign( 0.5, -zgdrho * ( zdep(ji,jj+1) - zdep(ji,jj) ) )
289               zkj(ji,jj) = ( 0.5 - zsign ) * zahv(ji,jj)
290            END DO
291         END DO
292         !
293      END SELECT
294
295      ! 2. Additional second order diffusive trends
296      ! -------------------------------------------
297
298      ! first derivative (gradient)
299#  if defined key_vectopt_loop
300      jj = 1
301      DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
302#  else
303      DO jj = 1, jpjm1
304         DO ji = 1, jpim1
305#  endif
306            zkx(ji,jj) = zki(ji,jj) * ( ztbb(ji+1,jj) - ztbb(ji,jj) )
307            zkz(ji,jj) = zki(ji,jj) * ( zsbb(ji+1,jj) - zsbb(ji,jj) )
308
309            zky(ji,jj) = zkj(ji,jj) * ( ztbb(ji,jj+1) - ztbb(ji,jj) )
310            zkw(ji,jj) = zkj(ji,jj) * ( zsbb(ji,jj+1) - zsbb(ji,jj) )
311#  if ! defined key_vectopt_loop
312         END DO
313#  endif
314      END DO
315
316      IF( cp_cfg == "orca" ) THEN
317         !
318         SELECT CASE ( jp_cfg )
319         !                                           ! =======================
320         CASE ( 2 )                                  !  ORCA_R2 configuration
321            !                                        ! =======================
322            ! Gibraltar enhancement of BBL
323            ij0 = 102   ;   ij1 = 102
324            ii0 = 139   ;   ii1 = 140 
325            zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 4.e0 * zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
326            zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 4.e0 * zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
327            !
328            ! Red Sea enhancement of BBL
329            ij0 =  88   ;   ij1 =  88
330            ii0 = 161   ;   ii1 = 162
331            zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 10.e0 * zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
332            zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 10.e0 * zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
333            !
334            !                                        ! =======================
335         CASE ( 4 )                                  !  ORCA_R4 configuration
336            !                                        ! =======================
337            ! Gibraltar enhancement of BBL
338            ij0 =  52   ;   ij1 =  52
339            ii0 =  70   ;   ii1 =  71 
340            zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 4.e0 * zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
341            zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 4.e0 * zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
342            !
343         END SELECT
344      !
345      ENDIF
346
347
348      ! second derivative (divergence) and add to the general tracer trend
349#  if defined key_vectopt_loop
350      DO jj = 1, 1
351         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
352#  else
353      DO jj = 2, jpjm1
354         DO ji = 2, jpim1
355#  endif
356            ik = max( mbathy(ji,jj)-1, 1 )
357            zbtr = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,ik) )
358            zta = (  zkx(ji,jj) - zkx(ji-1,jj  )    &
359                   + zky(ji,jj) - zky(ji  ,jj-1)  ) * zbtr
360            zsa = (  zkz(ji,jj) - zkz(ji-1,jj  )    &
361                   + zkw(ji,jj) - zkw(ji  ,jj-1)  ) * zbtr
362            ta(ji,jj,ik) = ta(ji,jj,ik) + zta
363            sa(ji,jj,ik) = sa(ji,jj,ik) + zsa
364         END DO
365      END DO
366
367      IF( l_trdtra ) THEN      ! save the BBL lateral diffusion trends for diagnostic
368         ztrdt(:,:,:) = ta(:,:,:) - ztrdt(:,:,:)
369         ztrds(:,:,:) = sa(:,:,:) - ztrds(:,:,:)
370         CALL trd_mod(ztrdt, ztrds, jptra_trd_bbl, 'TRA', kt)
371      ENDIF
372
373      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ta, clinfo1=' bbl  - Ta: ', mask1=tmask,   &
374         &                       tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
375      !
376   END SUBROUTINE tra_bbl_dif
377
378# if defined key_trabbl_adv
379   !!----------------------------------------------------------------------
380   !!   'key_trabbl_adv'                    advective bottom boundary layer
381   !!----------------------------------------------------------------------
382#  include "trabbl_adv.h90"
383# else
384   !!----------------------------------------------------------------------
385   !!   Default option :                 NO advective bottom boundary layer
386   !!----------------------------------------------------------------------
387   SUBROUTINE tra_bbl_adv (kt )              ! Empty routine
388      INTEGER, INTENT(in) :: kt
389      WRITE(*,*) 'tra_bbl_adv: You should not have seen this print! error?', kt
390   END SUBROUTINE tra_bbl_adv
391# endif
392
393   SUBROUTINE tra_bbl_init
394      !!----------------------------------------------------------------------
395      !!                  ***  ROUTINE tra_bbl_init  ***
396      !!
397      !! ** Purpose :   Initialization for the bottom boundary layer scheme.
398      !!
399      !! ** Method  :   Read the nambbl namelist and check the parameters
400      !!      called by tra_bbl at the first timestep (nit000)
401      !!----------------------------------------------------------------------
402      INTEGER ::   ji, jj      ! dummy loop indices
403      REAL(wp),  DIMENSION(jpi,jpj) :: zmbk 
404
405      NAMELIST/nambbl/ rn_ahtbbl
406      !!----------------------------------------------------------------------
407
408      REWIND ( numnam )              ! Read Namelist nambbl : bottom boundary layer scheme
409      READ   ( numnam, nambbl )
410
411      IF(lwp) THEN                   ! Parameter control and print
412         WRITE(numout,*)
413         WRITE(numout,*) 'tra_bbl_init : '
414         WRITE(numout,*) '~~~~~~~~~~~~'
415         IF( lk_trabbl_dif )   WRITE(numout,*) '               * Diffusive Bottom Boundary Layer'
416         IF( lk_trabbl_adv )   WRITE(numout,*) '               * Advective Bottom Boundary Layer'
417         WRITE(numout,*) '       Namelist nambbl : set bbl parameters'
418         WRITE(numout,*) '          bottom boundary layer coef.    rn_ahtbbl = ', rn_ahtbbl
419      ENDIF
420 
421      DO jj = 1, jpj
422         DO ji = 1, jpi
423            mbkt(ji,jj) = MAX( mbathy(ji,jj) - 1, 1 )   ! vertical index of the bottom ocean T-level
424         END DO
425      END DO
426      DO jj = 1, jpjm1
427         DO ji = 1, jpim1
428            mbku(ji,jj) = MAX( MIN( mbathy(ji+1,jj  ), mbathy(ji,jj) ) - 1, 1 )
429            mbkv(ji,jj) = MAX( MIN( mbathy(ji  ,jj+1), mbathy(ji,jj) ) - 1, 1 )
430         END DO
431      END DO
432
433      zmbk(:,:) = FLOAT( mbku (:,:) )   
434      CALL lbc_lnk(zmbk,'U',1.)
435      mbku(:,:) = MAX( INT( zmbk(:,:) ), 1 ) 
436
437      zmbk(:,:) = FLOAT( mbkv (:,:) )   
438      CALL lbc_lnk(zmbk,'V',1.)
439      mbkv(:,:) = MAX( INT( zmbk(:,:) ), 1 ) 
440
441# if defined key_trabbl_adv
442      w_bbl(:,:,:) = 0.e0          ! initialisation of w_bbl to zero
443# endif
444      !
445   END SUBROUTINE tra_bbl_init
446
447#else
448   !!----------------------------------------------------------------------
449   !!   Dummy module :                      No bottom boundary layer scheme
450   !!----------------------------------------------------------------------
451   LOGICAL, PUBLIC, PARAMETER ::   lk_trabbl_dif = .FALSE.   !: diff bbl flag
452   LOGICAL, PUBLIC, PARAMETER ::   lk_trabbl_adv = .FALSE.   !: adv  bbl flag
453CONTAINS
454   SUBROUTINE tra_bbl_dif( kt )              ! Empty routine
455      WRITE(*,*) 'tra_bbl_dif: You should not have seen this print! error?', kt
456   END SUBROUTINE tra_bbl_dif
457   SUBROUTINE tra_bbl_adv( kt )              ! Empty routine
458      WRITE(*,*) 'tra_bbl_adv: You should not have seen this print! error?', kt
459   END SUBROUTINE tra_bbl_adv
460#endif
461
462   !!======================================================================
463END MODULE trabbl
Note: See TracBrowser for help on using the repository browser.