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

Last change on this file since 484 was 481, checked in by opalod, 18 years ago

nemo_v1_bugfix_047 : CT : correction in order to take into account the partial steps for advective bbl

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