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

source: trunk/NEMO/OPA_SRC/TRA/trabbl.F90 @ 527

Last change on this file since 527 was 503, checked in by opalod, 18 years ago

nemo_v1_update_064 : CT : general trends update including the addition of mean windows analysis possibility in the mixed layer

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