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 @ 789

Last change on this file since 789 was 789, checked in by rblod, 16 years ago

Suppress jki routines and associated key_mpp_omp

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