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.
zdfdrg.F90 in branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfdrg.F90 @ 8882

Last change on this file since 8882 was 8882, checked in by flavoni, 6 years ago

dev_CNRS_2017 branch: merged dev_r7881_ENHANCE09_RK3 with trunk r8864

File size: 21.0 KB
Line 
1MODULE zdfdrg
2   !!======================================================================
3   !!                       ***  MODULE  zdfdrg  ***
4   !! Ocean physics: top and/or Bottom friction
5   !!======================================================================
6   !! History :  OPA  ! 1997-06  (G. Madec, A.-M. Treguier)  Original code
7   !!   NEMO     1.0  ! 2002-06  (G. Madec)  F90: Free form and module
8   !!            3.2  ! 2009-09  (A.C.Coward)  Correction to include barotropic contribution
9   !!            3.3  ! 2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase
10   !!            3.4  ! 2011-11  (H. Liu) implementation of semi-implicit bottom friction option
11   !!                 ! 2012-06  (H. Liu) implementation of Log Layer bottom friction option
12   !!            4.0  ! 2017-05  (G. Madec) zdfbfr becomes zdfdrg + variable names change
13   !!                                     + drag defined at t-point + new user interface + top drag (ocean cavities)
14   !!----------------------------------------------------------------------
15
16   !!----------------------------------------------------------------------
17   !!   zdf_drg       : update bottom friction coefficient (non-linear bottom friction only)
18   !!   zdf_drg_init  : read in namdrg namelist and control the bottom friction parameters.
19   !!       drg_init  :
20   !!----------------------------------------------------------------------
21   USE oce            ! ocean dynamics and tracers variables
22   USE phycst  , ONLY : vkarmn
23   USE dom_oce        ! ocean space and time domain variables
24   USE zdf_oce        ! ocean vertical physics variables
25   !
26   USE in_out_manager ! I/O manager
27   USE iom            ! I/O module
28   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
29   USE lib_mpp        ! distributed memory computing
30   USE prtctl         ! Print control
31   USE timing         ! Timing
32
33   IMPLICIT NONE
34   PRIVATE
35
36   PUBLIC   zdf_drg         ! called by zdf_phy
37   PUBLIC   zdf_drg_init    ! called by zdf_phy_init
38
39   !                                 !!* Namelist namdrg: nature of drag coefficient namelist *
40   LOGICAL          ::   ln_NONE      ! free-slip       : Cd = 0
41   LOGICAL          ::   ln_lin       !     linear  drag: Cd = Cd0_lin
42   LOGICAL          ::   ln_non_lin   ! non-linear  drag: Cd = Cd0_nl |U|
43   LOGICAL          ::   ln_loglayer  ! logarithmic drag: Cd = vkarmn/log(z/z0)
44   LOGICAL , PUBLIC ::   ln_drgimp    ! implicit top/bottom friction flag
45
46   !                                 !!* Namelist namdrg_top & _bot: TOP or BOTTOM coefficient namelist *
47   REAL(wp)         ::   rn_Cd0       !: drag coefficient                                           [ - ]
48   REAL(wp)         ::   rn_Uc0       !: characteristic velocity (linear case: tau=rho*Cd0*Uc0*u)   [m/s]
49   REAL(wp)         ::   rn_Cdmax     !: drag value maximum (ln_loglayer=T)                         [ - ]
50   REAL(wp)         ::   rn_z0        !: roughness          (ln_loglayer=T)                         [ m ]
51   REAL(wp)         ::   rn_ke0       !: background kinetic energy (non-linear case)                [m2/s2]
52   LOGICAL          ::   ln_boost     !: =T regional boost of Cd0 ; =F Cd0 horizontally uniform
53   REAL(wp)         ::     rn_boost      !: local boost factor                                       [ - ]
54
55   REAL(wp), PUBLIC ::   r_Cdmin_top, r_Cdmax_top, r_z0_top, r_ke0_top   ! set from namdrg_top namelist values
56   REAL(wp), PUBLIC ::   r_Cdmin_bot, r_Cdmax_bot, r_z0_bot, r_ke0_bot   !  -    -  namdrg_bot    -       -
57
58   INTEGER ::              ndrg       ! choice of the type of drag coefficient
59   !                                  ! associated indices:
60   INTEGER, PARAMETER ::   np_NONE     = 0   ! free-slip: drag set to zero
61   INTEGER, PARAMETER ::   np_lin      = 1   !     linear drag: Cd = Cd0_lin
62   INTEGER, PARAMETER ::   np_non_lin  = 2   ! non-linear drag: Cd = Cd0_nl |U|
63   INTEGER, PARAMETER ::   np_loglayer = 3   ! non linear drag (logarithmic formulation): Cd = vkarmn/log(z/z0)
64
65   LOGICAL , PUBLIC ::   l_zdfdrg           !: flag to update at each time step the top/bottom Cd
66   LOGICAL          ::   l_log_not_linssh   !: flag to update at each time step the position ot the velocity point
67   !
68   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:), PUBLIC ::   rCd0_top, rCd0_bot   !: precomputed top/bottom drag coeff. at t-point (>0)
69   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:), PUBLIC ::   rCdU_top, rCdU_bot   !: top/bottom drag coeff. at t-point (<0)  [m/s]
70
71   !! * Substitutions
72#  include "vectopt_loop_substitute.h90"
73   !!----------------------------------------------------------------------
74   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
75   !! $Id: zdfdrg.F90 7753 2017-03-03 11:46:59Z gm $
76   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
77   !!----------------------------------------------------------------------
78CONTAINS
79
80   SUBROUTINE zdf_drg( kt, k_mk, pCdmin, pCdmax, pz0, pke0, pCd0, pCdU )
81      !!----------------------------------------------------------------------
82      !!                   ***  ROUTINE zdf_drg  ***
83      !!
84      !! ** Purpose :   update the top/bottom drag coefficient (non-linear case only)
85      !!
86      !! ** Method  :   In non linear friction case, the drag coeficient is
87      !!              a function of the velocity:
88      !!                          Cd = cd0 * |U+Ut|   
89      !!              where U is the top or bottom velocity and
90      !!                    Ut a tidal velocity (Ut^2 = Tidal kinetic energy
91      !!                       assumed here here to be constant)
92      !!              Depending on the input variable, the top- or bottom drag is compted
93      !!
94      !! ** Action  :   p_Cd   drag coefficient at t-point
95      !!----------------------------------------------------------------------
96      INTEGER                 , INTENT(in   ) ::   kt       ! ocean time-step index
97      !                       !               !!         !==  top or bottom variables  ==!
98      INTEGER , DIMENSION(:,:), INTENT(in   ) ::   k_mk     ! wet level (1st or last)
99      REAL(wp)                , INTENT(in   ) ::   pCdmin   ! min drag value
100      REAL(wp)                , INTENT(in   ) ::   pCdmax   ! max drag value
101      REAL(wp)                , INTENT(in   ) ::   pz0      ! roughness
102      REAL(wp)                , INTENT(in   ) ::   pke0     ! background tidal KE
103      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   pCd0     ! masked precomputed part of Cd0
104      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   pCdU     ! = - Cd*|U|   (t-points) [m/s]
105      !!
106      INTEGER ::   ji, jj   ! dummy loop indices
107      INTEGER ::   imk      ! local integers
108      REAL(wp)::   zzz, zut, zvt, zcd   ! local scalars
109      !!----------------------------------------------------------------------
110      !
111      IF( ln_timing )   CALL timing_start('zdf_drg')
112      !
113      !
114      IF( l_log_not_linssh ) THEN     !==  "log layer"  ==!   compute Cd and -Cd*|U|
115         DO jj = 2, jpjm1
116            DO ji = 2, jpim1
117               imk = k_mk(ji,jj)          ! ocean bottom level at t-points
118               zut = un(ji,jj,imk) + un(ji-1,jj,imk)     ! 2 x velocity at t-point
119               zvt = vn(ji,jj,imk) + vn(ji,jj-1,imk)
120               zzz = 0.5_wp * e3t_n(ji,jj,imk)           ! altitude below/above (top/bottom) the boundary
121               !
122!!JC: possible WAD implementation should modify line below if layers vanish
123               zcd = (  vkarmn / LOG( zzz / pz0 )  )**2
124               zcd = pCd0(ji,jj) * MIN(  MAX( pCdmin , zcd ) , pCdmax  )   ! here pCd0 = mask*boost
125               pCdU(ji,jj) = - zcd * SQRT(  0.25 * ( zut*zvt + zvt*zvt ) + pke0  )
126            END DO
127         END DO
128      ELSE                                            !==  standard Cd  ==!
129         DO jj = 2, jpjm1
130            DO ji = 2, jpim1
131               imk = k_mk(ji,jj)    ! ocean bottom level at t-points
132               zut = un(ji,jj,imk) + un(ji-1,jj,imk)     ! 2 x velocity at t-point
133               zvt = vn(ji,jj,imk) + vn(ji,jj-1,imk)
134               !                                                           ! here pCd0 = mask*boost * drag
135               pCdU(ji,jj) = - pCd0(ji,jj) * SQRT(  0.25 * ( zut*zvt + zvt*zvt ) + pke0  )
136            END DO
137         END DO
138      ENDIF
139      !
140      IF(ln_ctl)   CALL prt_ctl( tab2d_1=pCdU, clinfo1=' Cd*U ')
141      !
142      IF( ln_timing )   CALL timing_stop('zdf_drg')
143      !
144   END SUBROUTINE zdf_drg
145
146
147   SUBROUTINE zdf_drg_init
148      !!----------------------------------------------------------------------
149      !!                  ***  ROUTINE zdf_brg_init  ***
150      !!
151      !! ** Purpose :   Initialization of the bottom friction
152      !!
153      !! ** Method  :   Read the namdrg namelist and check their consistency
154      !!                called at the first timestep (nit000)
155      !!----------------------------------------------------------------------
156      INTEGER   ::   ji, jj      ! dummy loop indexes
157      INTEGER   ::   ios, ioptio   ! local integers
158      !!
159      NAMELIST/namdrg/ ln_NONE, ln_lin, ln_non_lin, ln_loglayer, ln_drgimp
160      !!----------------------------------------------------------------------
161      !
162      !                     !==  drag nature  ==!
163      !
164      REWIND( numnam_ref )                   ! Namelist namdrg in reference namelist
165      READ  ( numnam_ref, namdrg, IOSTAT = ios, ERR = 901)
166901   IF( ios /= 0 )   CALL ctl_nam( ios , 'namdrg in reference namelist', lwp )
167      REWIND( numnam_cfg )                   ! Namelist namdrg in configuration namelist
168      READ  ( numnam_cfg, namdrg, IOSTAT = ios, ERR = 902 )
169902   IF( ios /= 0 )   CALL ctl_nam( ios , 'namdrg in configuration namelist', lwp )
170      IF(lwm) WRITE ( numond, namdrg )
171      !
172      IF(lwp) THEN
173         WRITE(numout,*)
174         WRITE(numout,*) 'zdf_drg_init : top and/or bottom drag setting'
175         WRITE(numout,*) '~~~~~~~~~~~~'
176         WRITE(numout,*) '   Namelist namdrg : top/bottom friction choices'
177         WRITE(numout,*) '      free-slip       : Cd = 0                  ln_NONE     = ', ln_NONE
178         WRITE(numout,*) '      linear  drag    : Cd = Cd0                ln_lin      = ', ln_lin
179         WRITE(numout,*) '      non-linear  drag: Cd = Cd0_nl |U|         ln_non_lin  = ', ln_non_lin
180         WRITE(numout,*) '      logarithmic drag: Cd = vkarmn/log(z/z0)   ln_loglayer = ', ln_loglayer
181         WRITE(numout,*) '      implicit friction                         ln_drgimp   = ', ln_drgimp
182      ENDIF
183      !
184      ioptio = 0                       ! set ndrg and control check
185      IF( ln_NONE     ) THEN   ;   ndrg = np_NONE       ;   ioptio = ioptio + 1   ;   ENDIF
186      IF( ln_lin      ) THEN   ;   ndrg = np_lin        ;   ioptio = ioptio + 1   ;   ENDIF
187      IF( ln_non_lin  ) THEN   ;   ndrg = np_non_lin    ;   ioptio = ioptio + 1   ;   ENDIF
188      IF( ln_loglayer ) THEN   ;   ndrg = np_loglayer   ;   ioptio = ioptio + 1   ;   ENDIF
189      !
190      IF( ioptio /= 1 )   CALL ctl_stop( 'zdf_drg_init: Choose ONE type of drag coef in namdrg' )
191      !
192      !
193      !                     !==  BOTTOM drag setting  ==!   (applied at seafloor)
194      !
195      ALLOCATE( rCd0_bot(jpi,jpj), rCdU_bot(jpi,jpj) )
196      CALL drg_init( 'BOTTOM'   , mbkt       ,                                         &   ! <== in
197         &           r_Cdmin_bot, r_Cdmax_bot, r_z0_bot, r_ke0_bot, rCd0_bot, rCdU_bot )   ! ==> out
198
199      !
200      !                     !==  TOP drag setting  ==!   (applied at the top of ocean cavities)
201      !
202      IF ( ln_isfcav ) THEN             ! Ocean cavities: top friction setting
203         ALLOCATE( rCd0_top(jpi,jpj), rCdU_top(jpi,jpj) )
204         CALL drg_init( 'TOP   '   , mikt       ,                                         &   ! <== in
205            &           r_Cdmin_top, r_Cdmax_top, r_z0_top, r_ke0_top, rCd0_bot, rCdU_bot )   ! ==> out
206      ENDIF
207      !
208   END SUBROUTINE zdf_drg_init
209
210
211   SUBROUTINE drg_init( cd_topbot, k_mk,  &
212      &                 pCdmin, pCdmax, pz0, pke0, pCd0, pCdU ) 
213      !!----------------------------------------------------------------------
214      !!                  ***  ROUTINE drg_init  ***
215      !!
216      !! ** Purpose :   Initialization of the top/bottom friction CdO and Cd
217      !!              from namelist parameters
218      !!----------------------------------------------------------------------
219      CHARACTER(len=6)        , INTENT(in   ) ::   cd_topbot       ! top/ bot indicator
220      INTEGER , DIMENSION(:,:), INTENT(in   ) ::   k_mk            ! 1st/last  wet level
221      REAL(wp)                , INTENT(  out) ::   pCdmin, pCdmax  ! min and max drag coef. [-]
222      REAL(wp)                , INTENT(  out) ::   pz0             ! roughness              [m]
223      REAL(wp)                , INTENT(  out) ::   pke0            ! background KE          [m2/s2]
224      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   pCd0            ! masked precomputed part of the non-linear drag coefficient
225      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   pCdU            ! linear drag*|U| at t-points  [m/s]
226      !!
227      CHARACTER(len=40) ::   cl_namdrg, cl_file, cl_varname, cl_namref, cl_namcfg  ! local names
228      INTEGER ::   ji, jj              ! dummy loop indexes
229      LOGICAL ::   ll_top, ll_bot      ! local logical
230      INTEGER ::   ios, inum, imk      ! local integers
231      REAL(wp)::   zmsk, zzz, zcd      ! local scalars
232      REAL(wp), DIMENSION(jpi,jpj) ::   zmsk_boost   ! 2D workspace
233      !!
234      NAMELIST/namdrg_top/ rn_Cd0, rn_Uc0, rn_Cdmax, rn_ke0, rn_z0, ln_boost, rn_boost
235      NAMELIST/namdrg_bot/ rn_Cd0, rn_Uc0, rn_Cdmax, rn_ke0, rn_z0, ln_boost, rn_boost
236      !!----------------------------------------------------------------------
237      !
238      !                          !==  set TOP / BOTTOM specificities  ==!
239      ll_top = .FALSE.
240      ll_bot = .FALSE.
241      !
242      SELECT CASE (cd_topbot)
243      CASE( 'TOP   ' )
244         ll_top = .TRUE.
245         cl_namdrg  = 'namdrg_top'
246         cl_namref  = 'namdrg_top in reference     namelist'
247         cl_namcfg  = 'namdrg_top in configuration namelist'
248         cl_file    = 'tfr_coef.nc'
249         cl_varname = 'tfr_coef'
250      CASE( 'BOTTOM' )
251         ll_bot = .TRUE.
252         cl_namdrg  = 'namdrg_bot'
253         cl_namref  = 'namdrg_bot  in reference     namelist'
254         cl_namcfg  = 'namdrg_bot  in configuration namelist'
255         cl_file    = 'bfr_coef.nc'
256         cl_varname = 'tfr_coef'
257      CASE DEFAULT
258         CALL ctl_stop( 'drg_init: bad value for cd_topbot ' )
259      END SELECT
260      !
261      !                          !==  read namlist  ==!
262      !
263      REWIND( numnam_ref )                   ! Namelist cl_namdrg in reference namelist
264      IF(ll_top)   READ  ( numnam_ref, namdrg_top, IOSTAT = ios, ERR = 901)
265      IF(ll_bot)   READ  ( numnam_ref, namdrg_bot, IOSTAT = ios, ERR = 901)
266901   IF( ios /= 0 )   CALL ctl_nam( ios , TRIM(cl_namref), lwp )
267      REWIND( numnam_cfg )                   ! Namelist cd_namdrg in configuration namelist
268      IF(ll_top)   READ  ( numnam_cfg, namdrg_top, IOSTAT = ios, ERR = 902 )
269      IF(ll_bot)   READ  ( numnam_cfg, namdrg_bot, IOSTAT = ios, ERR = 902 )
270902   IF( ios /= 0 )   CALL ctl_nam( ios , TRIM(cl_namcfg), lwp )
271      IF(lwm .AND. ll_top)   WRITE ( numond, namdrg_top )
272      IF(lwm .AND. ll_bot)   WRITE ( numond, namdrg_bot )
273      !
274      IF(lwp) THEN
275         WRITE(numout,*)
276         WRITE(numout,*) '   Namelist ',TRIM(cl_namdrg),' : set ',TRIM(cd_topbot),' friction parameters'
277         WRITE(numout,*) '      drag coefficient                        rn_Cd0   = ', rn_Cd0
278         WRITE(numout,*) '      characteristic velocity (linear case)   rn_Uc0   = ', rn_Uc0, ' m/s'
279         WRITE(numout,*) '      non-linear drag maximum                 rn_Cdmax = ', rn_Cdmax
280         WRITE(numout,*) '      background kinetic energy  (n-l case)   rn_ke0   = ', rn_ke0
281         WRITE(numout,*) '      bottom roughness           (n-l case)   rn_z0    = ', rn_z0
282         WRITE(numout,*) '      set a regional boost of Cd0             ln_boost = ', ln_boost
283         WRITE(numout,*) '         associated boost factor              rn_boost = ', rn_boost
284      ENDIF
285      !
286      !                          !==  return some namelist parametres  ==!   (used in non_lin and loglayer cases)
287      pCdmin = rn_Cd0
288      pCdmax = rn_Cdmax
289      pz0    = rn_z0
290      pke0   = rn_ke0
291      !
292      !                          !==  mask * boost factor  ==!
293      !
294      IF( ln_boost ) THEN           !* regional boost:   boost factor = 1 + regional boost
295         IF(lwp) WRITE(numout,*)
296         IF(lwp) WRITE(numout,*) '   ==>>   use a regional boost read in ', TRIM(cl_file), ' file'
297         IF(lwp) WRITE(numout,*) '          using enhancement factor of ', rn_boost
298         ! cl_varname is a coefficient in [0,1] giving where to apply the regional boost
299         CALL iom_open ( TRIM(cl_file), inum )
300         CALL iom_get  ( inum, jpdom_data, TRIM(cl_varname), zmsk_boost, 1 )
301         CALL iom_close( inum)
302         zmsk_boost(:,:) = 1._wp + rn_boost * zmsk_boost(:,:)
303         !
304      ELSE                          !* no boost:   boost factor = 1
305         zmsk_boost(:,:) = 1._wp
306      ENDIF
307      !                             !* mask outside ocean cavities area (top) or land area (bot)
308      IF(ll_top)   zmsk_boost(:,:) = zmsk_boost(:,:) * ssmask(:,:) * (1. - tmask(:,:,1) )  ! none zero in ocean cavities only
309      IF(ll_bot)   zmsk_boost(:,:) = zmsk_boost(:,:) * ssmask(:,:)                         ! x seafloor mask
310      !
311      !
312      SELECT CASE( ndrg )
313      !
314      CASE( np_NONE )            !==  No top/bottom friction  ==!   (pCdU = 0)
315         IF(lwp) WRITE(numout,*)
316         IF(lwp) WRITE(numout,*) '   ==>>   ',TRIM(cd_topbot),' free-slip, friction set to zero'
317         !
318         l_zdfdrg = .FALSE.         ! no time variation of the drag: set it one for all
319         !
320         pCdU(:,:) = 0._wp          ! pCd0 never used
321         !
322      CASE( np_lin )             !==  linear friction  ==!   (pCdU = Cd0 * Uc0)
323         IF(lwp) WRITE(numout,*)
324         IF(lwp) WRITE(numout,*) '   ==>>   linear ',TRIM(cd_topbot),' friction (constant coef = Cd0*Uc0 = ', rn_Cd0*rn_Uc0, ')'
325         !
326         l_zdfdrg = .FALSE.         ! no time variation of the Cd*|U| : set it one for all
327         !                     
328         pCdU(:,:) = - rn_Cd0 * rn_Uc0 * zmsk_boost(:,:)   ! pCd0 never used: remain undefined
329         !
330      CASE( np_non_lin )         !== non-linear friction  ==!   (pCd0 = Cd0 )
331         IF(lwp) WRITE(numout,*)
332         IF(lwp) WRITE(numout,*) '   ==>>   quadratic ',TRIM(cd_topbot),' friction (propotional to module of the velocity)'
333         IF(lwp) WRITE(numout,*) '   with   Cd0 = ', rn_Cd0, ', and',   &
334            &                             ' a background velocity module of (rn_ke0)^1/2 = ', SQRT(rn_ke0), 'm/s)'
335         !
336         l_zdfdrg = .TRUE.          !* Cd*|U| updated at each time-step (it depends on ocean velocity)
337         !
338         pCd0(:,:) = rn_Cd0 * zmsk_boost(:,:)  !* constant in time proportionality coefficient (= mask (and boost) Cd0)
339         !
340      CASE( np_loglayer )       !== logarithmic layer formulation of friction  ==!   (CdU = (vkarman log(z/z0))^2 |U| )
341         IF(lwp) WRITE(numout,*)
342         IF(lwp) WRITE(numout,*) '   ==>>   quadratic ',TRIM(cd_topbot),' drag (propotional to module of the velocity)'
343         IF(lwp) WRITE(numout,*) '   with   a logarithmic Cd0 formulation Cd0 = ( vkarman log(z/z0) )^2 ,'
344         IF(lwp) WRITE(numout,*) '          a background velocity module of (rn_ke0)^1/2 = ', SQRT(pke0), 'm/s), '
345         IF(lwp) WRITE(numout,*) '          a logarithmic formulation: a roughness of ', pz0, ' meters,   and '
346         IF(lwp) WRITE(numout,*) '          a proportionality factor bounded by min/max values of ', pCdmin, pCdmax
347         !
348         l_zdfdrg = .TRUE.          !* Cd*|U| updated at each time-step (it depends on ocean velocity)
349         !
350         IF( ln_linssh ) THEN       !* pCd0 = (v log(z/z0))^2   as velocity points have a fixed z position
351         IF(lwp) WRITE(numout,*)
352         IF(lwp) WRITE(numout,*) '   N.B.   linear free surface case, Cd0 computed one for all'
353            !
354            l_log_not_linssh = .FALSE.    !- don't update Cd at each time step
355            !
356            DO jj = 1, jpj                   ! pCd0 = mask (and boosted) logarithmic drag coef.
357               DO ji = 1, jpi
358                  zzz =  0.5_wp * e3t_0(ji,jj,k_mk(ji,jj))
359                  zcd = (  vkarmn / LOG( zzz / rn_z0 )  )**2
360                  pCd0(ji,jj) = zmsk_boost(ji,jj) * MIN(  MAX( rn_Cd0 , zcd ) , rn_Cdmax  )  ! rn_Cd0 < Cd0 < rn_Cdmax
361               END DO
362            END DO
363         ELSE                       !* Cd updated at each time-step ==> pCd0 = mask * boost
364            IF(lwp) WRITE(numout,*)
365            IF(lwp) WRITE(numout,*) '   N.B.   non-linear free surface case, Cd0 updated at each time-step '
366            !
367            l_log_not_linssh = .TRUE.     ! compute the drag coef. at each time-step
368            !
369            pCd0(:,:) = zmsk_boost(:,:)
370         ENDIF
371         !
372      CASE DEFAULT
373         CALL ctl_stop( 'drg_init: bad flag value for ndrg ' )
374      END SELECT
375      !
376   END SUBROUTINE drg_init
377
378   !!======================================================================
379END MODULE zdfdrg
Note: See TracBrowser for help on using the repository browser.