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_merge_2017/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfdrg.F90 @ 9190

Last change on this file since 9190 was 9190, checked in by gm, 6 years ago

dev_merge_2017: OPA_SRC: style only, results unchanged

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