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 branches/dev_001_GM/NEMO/OPA_SRC/TRA – NEMO

source: branches/dev_001_GM/NEMO/OPA_SRC/TRA/trabbl.F90 @ 790

Last change on this file since 790 was 790, checked in by gm, 16 years ago

dev_001_GM - complete theprevious comit with omitted routines

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 22.5 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 ::   atrbbl = 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   !! $Header$
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 arrays
119      REAL(wp), DIMENSION(jpi,jpj) ::   ztnb, zsnb, zdep
120      REAL(wp), DIMENSION(jpi,jpj) ::   ztbb, zsbb, zahu, zahv
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         ! vector opt.
156         DO ji = 1, jpij   ! forced the loop collapse
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             ! vector opt.
173         DO ji = 1, jpij-jpi   ! forced the loop collapse
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) = atrbbl * e2u(ji,jj) * ze3u / e1u(ji,jj) * umask(ji,jj,1)
185               zahv(ji,jj) = atrbbl * 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   &&   ! defined key_mpp_omp
190      DO jj = 1, 1             ! vector opt.
191         DO ji = 1, jpij-jpi   ! forced the loop collapse
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) = atrbbl * e2u(ji,jj) * fse3u(ji,jj,iku) / e1u(ji,jj) * umask(ji,jj,1)
199               zahv(ji,jj) = atrbbl * 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 ( neos )
210
211      CASE ( 0 )                 ! Jackett and McDougall (1994) formulation
212
213#  if defined key_vectopt_loop   &&   ! defined key_mpp_omp
214      DO jj = 1, 1             ! vector opt.
215         DO ji = 1, jpij-jpi   ! forced the loop collapse
216#  else
217      DO jj = 1, jpjm1
218         DO ji = 1, jpim1
219#  endif
220            ! temperature, salinity anomalie and depth
221            zt = 0.5 * ( ztnb(ji,jj) + ztnb(ji+1,jj) )
222            zs = 0.5 * ( zsnb(ji,jj) + zsnb(ji+1,jj) ) - 35.0
223            zh = 0.5 * ( zdep(ji,jj) + zdep(ji+1,jj) )
224            ! masked ratio alpha/beta
225            zalbet = fsalbt( zt, zs, zh )*umask(ji,jj,1)
226            ! local density gradient along i-bathymetric slope
227            zgdrho = zalbet * ( ztnb(ji+1,jj) - ztnb(ji,jj) )   &
228                   -          ( zsnb(ji+1,jj) - zsnb(ji,jj) )
229            ! sign of local i-gradient of density multiplied by the i-slope
230            zsign = SIGN( 0.5, - zgdrho * ( zdep(ji+1,jj) - zdep(ji,jj) ) )
231            zki(ji,jj) = ( 0.5 - zsign ) * zahu(ji,jj)
232         END DO
233      END DO
234
235#  if defined key_vectopt_loop   &&   ! defined key_mpp_omp
236      DO jj = 1, 1             ! vector opt.
237         DO ji = 1, jpij-jpi   ! forced the loop collapse
238#  else
239      DO jj = 1, jpjm1
240         DO ji = 1, jpim1
241#  endif
242            ! temperature, salinity anomalie and depth
243            zt = 0.5 * ( ztnb(ji,jj+1) + ztnb(ji,jj) )
244            zs = 0.5 * ( zsnb(ji,jj+1) + zsnb(ji,jj) ) - 35.0
245            zh = 0.5 * ( zdep(ji,jj+1) + zdep(ji,jj) )
246            ! masked ratio alpha/beta
247            zalbet = fsalbt( zt, zs, zh )*vmask(ji,jj,1)
248            ! local density gradient along j-bathymetric slope
249            zgdrho = zalbet * ( ztnb(ji,jj+1) - ztnb(ji,jj) )   &
250               &   -          ( zsnb(ji,jj+1) - zsnb(ji,jj) )
251            ! sign of local j-gradient of density multiplied by the j-slope
252            zsign = sign( 0.5, -zgdrho * ( zdep(ji,jj+1) - zdep(ji,jj) ) )
253            zkj(ji,jj) = ( 0.5 - zsign ) * zahv(ji,jj)
254         END DO
255      END DO
256
257      CASE ( 1 )               ! Linear formulation function of temperature only
258                               !
259#  if defined key_vectopt_loop   &&   ! defined key_mpp_omp
260      DO jj = 1, 1             ! vector opt.
261         DO ji = 1, jpij-jpi   ! forced the loop collapse
262#  else
263      DO jj = 1, jpjm1
264         DO ji = 1, jpim1
265#  endif
266            ! local 'density/temperature' gradient along i-bathymetric slope
267            zgdrho =  ztnb(ji+1,jj) - ztnb(ji,jj) 
268            ! sign of local i-gradient of density multiplied by the i-slope
269            zsign = SIGN( 0.5, - zgdrho * ( zdep(ji+1,jj) - zdep(ji,jj) ) )
270            zki(ji,jj) = ( 0.5 - zsign ) * zahu(ji,jj)
271         END DO
272      END DO
273
274#  if defined key_vectopt_loop   &&   ! defined key_mpp_omp
275      DO jj = 1, 1             ! vector opt.
276         DO ji = 1, jpij-jpi   ! forced the loop collapse
277#  else
278      DO jj = 1, jpjm1
279         DO ji = 1, jpim1
280#  endif
281            ! local density gradient along j-bathymetric slope
282            zgdrho =  ztnb(ji,jj+1) - ztnb(ji,jj) 
283            ! sign of local j-gradient of density multiplied by the j-slope
284            zsign = sign( 0.5, -zgdrho * ( zdep(ji,jj+1) - zdep(ji,jj) ) )
285            zkj(ji,jj) = ( 0.5 - zsign ) * zahv(ji,jj)
286         END DO
287      END DO
288
289      CASE ( 2 )               ! Linear formulation function of temperature and salinity
290
291#  if defined key_vectopt_loop   &&   ! defined key_mpp_omp
292      DO jj = 1, 1             ! vector opt.
293         DO ji = 1, jpij-jpi   ! forced the loop collapse
294#  else
295      DO jj = 1, jpjm1
296         DO ji = 1, jpim1
297#  endif     
298            ! local density gradient along i-bathymetric slope
299            zgdrho = - ( rbeta*( zsnb(ji+1,jj) - zsnb(ji,jj) )   &
300               &     -  ralpha*( ztnb(ji+1,jj) - ztnb(ji,jj) ) )
301            ! sign of local i-gradient of density multiplied by the i-slope
302            zsign = SIGN( 0.5, - zgdrho * ( zdep(ji+1,jj) - zdep(ji,jj) ) )
303            zki(ji,jj) = ( 0.5 - zsign ) * zahu(ji,jj)
304         END DO
305      END DO
306
307#  if defined key_vectopt_loop   &&   ! defined key_mpp_omp
308      DO jj = 1, 1             ! vector opt.
309         DO ji = 1, jpij-jpi   ! forced the loop collapse
310#  else
311      DO jj = 1, jpjm1
312         DO ji = 1, jpim1
313#  endif     
314            ! local density gradient along j-bathymetric slope
315            zgdrho = - ( rbeta*( zsnb(ji,jj+1) - zsnb(ji,jj) )   &
316               &     -  ralpha*( ztnb(ji,jj+1) - ztnb(ji,jj) ) )   
317            ! sign of local j-gradient of density multiplied by the j-slope
318            zsign = SIGN( 0.5, -zgdrho * ( zdep(ji,jj+1) - zdep(ji,jj) ) )
319            zkj(ji,jj) = ( 0.5 - zsign ) * zahv(ji,jj)
320         END DO
321      END DO
322     
323      CASE DEFAULT
324
325         WRITE(ctmp1,*) '          bad flag value for neos = ', neos
326         CALL ctl_stop(ctmp1)
327
328      END SELECT
329
330      ! 2. Additional second order diffusive trends
331      ! -------------------------------------------
332
333      ! first derivative (gradient)
334#  if defined key_vectopt_loop   &&   ! defined key_mpp_omp
335      DO jj = 1, 1             ! vector opt.
336         DO ji = 1, jpij-jpi   ! forced the loop collapse
337#  else
338      DO jj = 1, jpjm1
339         DO ji = 1, jpim1
340#  endif
341            zkx(ji,jj) = zki(ji,jj) * ( ztbb(ji+1,jj) - ztbb(ji,jj) )
342            zkz(ji,jj) = zki(ji,jj) * ( zsbb(ji+1,jj) - zsbb(ji,jj) )
343
344            zky(ji,jj) = zkj(ji,jj) * ( ztbb(ji,jj+1) - ztbb(ji,jj) )
345            zkw(ji,jj) = zkj(ji,jj) * ( zsbb(ji,jj+1) - zsbb(ji,jj) )
346         END DO
347      END DO
348
349      IF( cp_cfg == "orca" ) THEN
350         !
351         SELECT CASE ( jp_cfg )
352         !                                           ! =======================
353         CASE ( 2 )                                  !  ORCA_R2 configuration
354            !                                        ! =======================
355            ! Gibraltar enhancement of BBL
356            ij0 = 102   ;   ij1 = 102
357            ii0 = 139   ;   ii1 = 140 
358            zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 4.e0 * zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
359            zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 4.e0 * zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
360            !
361            ! Red Sea enhancement of BBL
362            ij0 =  88   ;   ij1 =  88
363            ii0 = 161   ;   ii1 = 162
364            zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 10.e0 * zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
365            zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 10.e0 * zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
366            !
367            !                                        ! =======================
368         CASE ( 4 )                                  !  ORCA_R4 configuration
369            !                                        ! =======================
370            ! Gibraltar enhancement of BBL
371            ij0 =  52   ;   ij1 =  52
372            ii0 =  70   ;   ii1 =  71 
373            zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 4.e0 * zkx( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
374            zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) = 4.e0 * zky( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) )
375            !
376         END SELECT
377      !
378      ENDIF
379
380
381      ! second derivative (divergence) and add to the general tracer trend
382#  if defined key_vectopt_loop   &&   ! defined key_mpp_omp
383      DO jj = 1, 1                   ! vector opt.
384         DO ji = jpi+2, jpij-jpi-1   ! forced the loop collapse
385#  else
386      DO jj = 2, jpjm1
387         DO ji = 2, jpim1
388#  endif
389            ik = max( mbathy(ji,jj)-1, 1 )
390            zbtr = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,ik) )
391            zta = (  zkx(ji,jj) - zkx(ji-1,jj  )    &
392               &   + zky(ji,jj) - zky(ji  ,jj-1)  ) * zbtr
393            zsa = (  zkz(ji,jj) - zkz(ji-1,jj  )    &
394               &   + zkw(ji,jj) - zkw(ji  ,jj-1)  ) * zbtr
395            ta(ji,jj,ik) = ta(ji,jj,ik) + zta
396            sa(ji,jj,ik) = sa(ji,jj,ik) + zsa
397         END DO
398      END DO
399
400      IF( l_trdtra ) THEN      ! save the BBL lateral diffusion trends for diagnostic
401         ztrdt(:,:,:) = ta(:,:,:) - ztrdt(:,:,:)
402         ztrds(:,:,:) = sa(:,:,:) - ztrds(:,:,:)
403         CALL trd_tra( kt, jp_tem, jpt_trd_bbl, 'TRA', ptrd3d=ztrdt)
404         CALL trd_tra( kt, jp_sal, jpt_trd_bbl, 'TRA', ptrd3d=ztrds)
405      ENDIF
406
407      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ta, clinfo1=' bbl  - Ta: ', mask1=tmask,   &
408         &                       tab3d_2=sa, clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
409      !
410   END SUBROUTINE tra_bbl_dif
411
412# if defined key_trabbl_adv
413   !!----------------------------------------------------------------------
414   !!   'key_trabbl_adv'                    advective bottom boundary layer
415   !!----------------------------------------------------------------------
416#  include "trabbl_adv.h90"
417# else
418   !!----------------------------------------------------------------------
419   !!   Default option :                 NO advective bottom boundary layer
420   !!----------------------------------------------------------------------
421   SUBROUTINE tra_bbl_adv( kt )              ! Empty routine
422      INTEGER :: kt
423      WRITE(*,*) 'tra_bbl_adv: You should not have seen this print! error?', kt
424   END SUBROUTINE tra_bbl_adv
425# endif
426
427   SUBROUTINE tra_bbl_init
428      !!----------------------------------------------------------------------
429      !!                  ***  ROUTINE tra_bbl_init  ***
430      !!
431      !! ** Purpose :   Initialization for the bottom boundary layer scheme.
432      !!
433      !! ** Method  :   Read the nambbl namelist and check the parameters
434      !!      called by tra_bbl at the first timestep (nit000)
435      !!----------------------------------------------------------------------
436      INTEGER ::   ji, jj      ! dummy loop indices
437      REAL(wp),  DIMENSION(jpi,jpj) :: zmbk 
438
439      NAMELIST/nambbl/ atrbbl
440      !!----------------------------------------------------------------------
441
442      REWIND ( numnam )              ! Read Namelist nambbl : bottom boundary layer scheme
443      READ   ( numnam, nambbl )
444
445      IF(lwp) THEN                   ! Parameter control and print
446         WRITE(numout,*)
447         WRITE(numout,*) 'tra_bbl_init : '
448         WRITE(numout,*) '~~~~~~~~~~~~'
449         IF (lk_trabbl_dif )   WRITE(numout,*) '               * Diffusive Bottom Boundary Layer'
450         IF( lk_trabbl_adv )   WRITE(numout,*) '               * Advective Bottom Boundary Layer'
451         WRITE(numout,*) '       Namelist nambbl : set bbl parameters'
452         WRITE(numout,*) '          bottom boundary layer coef.    atrbbl = ', atrbbl
453      ENDIF
454 
455      DO jj = 1, jpj
456         DO ji = 1, jpi
457            mbkt(ji,jj) = MAX( mbathy(ji,jj) - 1, 1 )   ! vertical index of the bottom ocean T-level
458         END DO
459      END DO
460      DO jj = 1, jpjm1
461         DO ji = 1, jpim1
462            mbku(ji,jj) = MAX( MIN( mbathy(ji+1,jj  ), mbathy(ji,jj) ) - 1, 1 )
463            mbkv(ji,jj) = MAX( MIN( mbathy(ji  ,jj+1), mbathy(ji,jj) ) - 1, 1 )
464         END DO
465      END DO
466
467      zmbk(:,:) = FLOAT( mbku (:,:) )   
468      CALL lbc_lnk(zmbk,'U',1.)
469      mbku(:,:) = MAX( INT( zmbk(:,:) ), 1 ) 
470
471      zmbk(:,:) = FLOAT( mbkv (:,:) )   
472      CALL lbc_lnk(zmbk,'V',1.)
473      mbkv(:,:) = MAX( INT( zmbk(:,:) ), 1 ) 
474
475# if defined key_trabbl_adv
476      w_bbl(:,:,:) = 0.e0          ! initialisation of w_bbl to zero
477# endif
478      !
479   END SUBROUTINE tra_bbl_init
480
481#else
482   !!----------------------------------------------------------------------
483   !!   Dummy module :                      No bottom boundary layer scheme
484   !!----------------------------------------------------------------------
485   LOGICAL, PUBLIC, PARAMETER ::   lk_trabbl_dif = .FALSE.   !: diff bbl flag
486   LOGICAL, PUBLIC, PARAMETER ::   lk_trabbl_adv = .FALSE.   !: adv  bbl flag
487CONTAINS
488   SUBROUTINE tra_bbl_dif( kt )              ! Empty routine
489      WRITE(*,*) 'tra_bbl_dif: You should not have seen this print! error?', kt
490   END SUBROUTINE tra_bbl_dif
491   SUBROUTINE tra_bbl_adv( kt )              ! Empty routine
492      WRITE(*,*) 'tra_bbl_adv: You should not have seen this print! error?', kt
493   END SUBROUTINE tra_bbl_adv
494#endif
495
496   !!======================================================================
497END MODULE trabbl
Note: See TracBrowser for help on using the repository browser.