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 NEMO/branches/2019/dev_r11351_fldread_with_XIOS/src/OCE/ZDF – NEMO

source: NEMO/branches/2019/dev_r11351_fldread_with_XIOS/src/OCE/ZDF/zdfdrg.F90 @ 13463

Last change on this file since 13463 was 13463, checked in by andmirek, 4 years ago

Ticket #2195:update to trunk 13461

  • Property svn:keywords set to Id
File size: 24.5 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_exp   : compute the top & bottom friction in explicit case
19   !!   zdf_drg_init  : read in namdrg namelist and control the bottom friction parameters.
20   !!       drg_init  :
21   !!----------------------------------------------------------------------
22   USE oce            ! ocean dynamics and tracers variables
23   USE phycst  , ONLY : vkarmn
24   USE dom_oce        ! ocean space and time domain variables
25   USE zdf_oce        ! ocean vertical physics variables
26   USE trd_oce        ! trends: ocean variables
27   USE trddyn         ! trend manager: dynamics
28   !
29   USE in_out_manager ! I/O manager
30   USE iom            ! I/O module
31   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
32   USE lib_mpp        ! distributed memory computing
33   USE prtctl         ! Print control
34
35   IMPLICIT NONE
36   PRIVATE
37
38   PUBLIC   zdf_drg         ! called by zdf_phy
39   PUBLIC   zdf_drg_exp     ! called by dyn_zdf
40   PUBLIC   zdf_drg_init    ! called by zdf_phy_init
41
42   !                                 !!* Namelist namdrg: nature of drag coefficient namelist *
43   LOGICAL , PUBLIC ::   ln_drg_OFF   ! free-slip       : Cd = 0
44   LOGICAL          ::   ln_lin       !     linear  drag: Cd = Cd0_lin
45   LOGICAL          ::   ln_non_lin   ! non-linear  drag: Cd = Cd0_nl |U|
46   LOGICAL          ::   ln_loglayer  ! logarithmic drag: Cd = vkarmn/log(z/z0)
47   LOGICAL , PUBLIC ::   ln_drgimp    ! implicit top/bottom friction flag
48
49   !                                 !!* Namelist namdrg_top & _bot: TOP or BOTTOM coefficient namelist *
50   REAL(wp)         ::   rn_Cd0       !: drag coefficient                                           [ - ]
51   REAL(wp)         ::   rn_Uc0       !: characteristic velocity (linear case: tau=rho*Cd0*Uc0*u)   [m/s]
52   REAL(wp)         ::   rn_Cdmax     !: drag value maximum (ln_loglayer=T)                         [ - ]
53   REAL(wp)         ::   rn_z0        !: roughness          (ln_loglayer=T)                         [ m ]
54   REAL(wp)         ::   rn_ke0       !: background kinetic energy (non-linear case)                [m2/s2]
55   LOGICAL          ::   ln_boost     !: =T regional boost of Cd0 ; =F Cd0 horizontally uniform
56   REAL(wp)         ::   rn_boost     !: local boost factor                                         [ - ]
57
58   REAL(wp), PUBLIC ::   r_Cdmin_top, r_Cdmax_top, r_z0_top, r_ke0_top   ! set from namdrg_top namelist values
59   REAL(wp), PUBLIC ::   r_Cdmin_bot, r_Cdmax_bot, r_z0_bot, r_ke0_bot   !  -    -  namdrg_bot    -       -
60
61   INTEGER ::              ndrg       ! choice of the type of drag coefficient
62   !                                  ! associated indices:
63   INTEGER, PARAMETER ::   np_OFF      = 0   ! free-slip: drag set to zero
64   INTEGER, PARAMETER ::   np_lin      = 1   !     linear drag: Cd = Cd0_lin
65   INTEGER, PARAMETER ::   np_non_lin  = 2   ! non-linear drag: Cd = Cd0_nl |U|
66   INTEGER, PARAMETER ::   np_loglayer = 3   ! non linear drag (logarithmic formulation): Cd = vkarmn/log(z/z0)
67
68   LOGICAL , PUBLIC ::   l_zdfdrg           !: flag to update at each time step the top/bottom Cd
69   LOGICAL          ::   l_log_not_linssh   !: flag to update at each time step the position ot the velocity point
70   !
71   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:), PUBLIC ::   rCd0_top, rCd0_bot   !: precomputed top/bottom drag coeff. at t-point (>0)
72   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:), PUBLIC ::   rCdU_top, rCdU_bot   !: top/bottom drag coeff. at t-point (<0)  [m/s]
73
74   !! * Substitutions
75#  include "do_loop_substitute.h90"
76#  include "domzgr_substitute.h90"
77   !!----------------------------------------------------------------------
78   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
79   !! $Id$
80   !! Software governed by the CeCILL license (see ./LICENSE)
81   !!----------------------------------------------------------------------
82CONTAINS
83
84   SUBROUTINE zdf_drg( kt, Kmm, k_mk, pCdmin, pCdmax, pz0, pke0, pCd0,   &   ! <<== in
85      &                                                     pCdU )      ! ==>> out : bottom drag [m/s]
86      !!----------------------------------------------------------------------
87      !!                   ***  ROUTINE zdf_drg  ***
88      !!
89      !! ** Purpose :   update the top/bottom drag coefficient (non-linear case only)
90      !!
91      !! ** Method  :   In non linear friction case, the drag coeficient is
92      !!              a function of the velocity:
93      !!                          Cd = cd0 * |U+Ut|   
94      !!              where U is the top or bottom velocity and
95      !!                    Ut a tidal velocity (Ut^2 = Tidal kinetic energy
96      !!                       assumed here here to be constant)
97      !!              Depending on the input variable, the top- or bottom drag is compted
98      !!
99      !! ** Action  :   p_Cd   drag coefficient at t-point
100      !!----------------------------------------------------------------------
101      INTEGER                 , INTENT(in   ) ::   kt       ! ocean time-step index
102      INTEGER                 , INTENT(in   ) ::   Kmm      ! ocean time level index
103      !                       !               !!         !==  top or bottom variables  ==!
104      INTEGER , DIMENSION(:,:), INTENT(in   ) ::   k_mk     ! wet level (1st or last)
105      REAL(wp)                , INTENT(in   ) ::   pCdmin   ! min drag value
106      REAL(wp)                , INTENT(in   ) ::   pCdmax   ! max drag value
107      REAL(wp)                , INTENT(in   ) ::   pz0      ! roughness
108      REAL(wp)                , INTENT(in   ) ::   pke0     ! background tidal KE
109      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   pCd0     ! masked precomputed part of Cd0
110      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   pCdU     ! = - Cd*|U|   (t-points) [m/s]
111      !!
112      INTEGER ::   ji, jj   ! dummy loop indices
113      INTEGER ::   imk      ! local integers
114      REAL(wp)::   zzz, zut, zvt, zcd   ! local scalars
115      !!----------------------------------------------------------------------
116      !
117      IF( l_log_not_linssh ) THEN     !==  "log layer"  ==!   compute Cd and -Cd*|U|
118         DO_2D( 0, 0, 0, 0 )
119            imk = k_mk(ji,jj)          ! ocean bottom level at t-points
120            zut = uu(ji,jj,imk,Kmm) + uu(ji-1,jj,imk,Kmm)     ! 2 x velocity at t-point
121            zvt = vv(ji,jj,imk,Kmm) + vv(ji,jj-1,imk,Kmm)
122            zzz = 0.5_wp * e3t(ji,jj,imk,Kmm)           ! altitude below/above (top/bottom) the boundary
123            !
124!!JC: possible WAD implementation should modify line below if layers vanish
125            zcd = (  vkarmn / LOG( zzz / pz0 )  )**2
126            zcd = pCd0(ji,jj) * MIN(  MAX( pCdmin , zcd ) , pCdmax  )   ! here pCd0 = mask*boost
127            pCdU(ji,jj) = - zcd * SQRT(  0.25 * ( zut*zut + zvt*zvt ) + pke0  )
128         END_2D
129      ELSE                                            !==  standard Cd  ==!
130         DO_2D( 0, 0, 0, 0 )
131            imk = k_mk(ji,jj)    ! ocean bottom level at t-points
132            zut = uu(ji,jj,imk,Kmm) + uu(ji-1,jj,imk,Kmm)     ! 2 x velocity at t-point
133            zvt = vv(ji,jj,imk,Kmm) + vv(ji,jj-1,imk,Kmm)
134            !                                                           ! here pCd0 = mask*boost * drag
135            pCdU(ji,jj) = - pCd0(ji,jj) * SQRT(  0.25 * ( zut*zut + zvt*zvt ) + pke0  )
136         END_2D
137      ENDIF
138      !
139      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=pCdU, clinfo1=' Cd*U ')
140      !
141   END SUBROUTINE zdf_drg
142
143
144   SUBROUTINE zdf_drg_exp( kt, Kmm, pub, pvb, pua, pva )
145      !!----------------------------------------------------------------------
146      !!                  ***  ROUTINE zdf_drg_exp  ***
147      !!
148      !! ** Purpose :   compute and add the explicit top and bottom frictions.
149      !!
150      !! ** Method  :   in explicit case,
151      !!
152      !!              NB: in implicit case the calculation is performed in dynzdf.F90
153      !!
154      !! ** Action  :   (pua,pva)   momentum trend increased by top & bottom friction trend
155      !!---------------------------------------------------------------------
156      INTEGER                         , INTENT(in   ) ::   kt         ! ocean time-step index
157      INTEGER                         , INTENT(in   ) ::   Kmm        ! time level indices
158      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pub, pvb   ! the two components of the before velocity
159      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pua, pva   ! the two components of the velocity tendency
160      !!
161      INTEGER  ::   ji, jj       ! dummy loop indexes
162      INTEGER  ::   ikbu, ikbv   ! local integers
163      REAL(wp) ::   zm1_2dt      ! local scalar
164      REAL(wp) ::   zCdu, zCdv   !   -      -
165      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ztrdu, ztrdv
166      !!---------------------------------------------------------------------
167      !
168!!gm bug : time step is only rn_Dt (not 2 rn_Dt if euler start !)
169      zm1_2dt = - 1._wp / ( 2._wp * rn_Dt )
170
171      IF( l_trddyn ) THEN      ! trends: store the input trends
172         ALLOCATE( ztrdu(jpi,jpj,jpk) , ztrdv(jpi,jpj,jpk) )
173         ztrdu(:,:,:) = pua(:,:,:)
174         ztrdv(:,:,:) = pva(:,:,:)
175      ENDIF
176
177      DO_2D( 0, 0, 0, 0 )
178         ikbu = mbku(ji,jj)          ! deepest wet ocean u- & v-levels
179         ikbv = mbkv(ji,jj)
180         !
181         ! Apply stability criteria on absolute value  : abs(bfr/e3) < 1/(2dt) => bfr/e3 > -1/(2dt)
182         zCdu = 0.5*( rCdU_bot(ji+1,jj)+rCdU_bot(ji,jj) ) / e3u(ji,jj,ikbu,Kmm)
183         zCdv = 0.5*( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj) ) / e3v(ji,jj,ikbv,Kmm)
184         !
185         pua(ji,jj,ikbu) = pua(ji,jj,ikbu) + MAX(  zCdu , zm1_2dt  ) * pub(ji,jj,ikbu)
186         pva(ji,jj,ikbv) = pva(ji,jj,ikbv) + MAX(  zCdv , zm1_2dt  ) * pvb(ji,jj,ikbv)
187      END_2D
188      !
189      IF( ln_isfcav ) THEN        ! ocean cavities
190         DO_2D( 0, 0, 0, 0 )
191            ikbu = miku(ji,jj)          ! first wet ocean u- & v-levels
192            ikbv = mikv(ji,jj)
193            !
194            ! Apply stability criteria on absolute value  : abs(bfr/e3) < 1/(2dt) => bfr/e3 > -1/(2dt)
195            zCdu = 0.5*( rCdU_top(ji+1,jj)+rCdU_top(ji,jj) ) / e3u(ji,jj,ikbu,Kmm)    ! NB: Cdtop masked
196            zCdv = 0.5*( rCdU_top(ji,jj+1)+rCdU_top(ji,jj) ) / e3v(ji,jj,ikbv,Kmm)
197            !
198            pua(ji,jj,ikbu) = pua(ji,jj,ikbu) + MAX(  zCdu , zm1_2dt  ) * pub(ji,jj,ikbu)
199            pva(ji,jj,ikbv) = pva(ji,jj,ikbv) + MAX(  zCdv , zm1_2dt  ) * pvb(ji,jj,ikbv)
200         END_2D
201      ENDIF
202      !
203      IF( l_trddyn ) THEN      ! trends: send trends to trddyn for further diagnostics
204         ztrdu(:,:,:) = pua(:,:,:) - ztrdu(:,:,:)
205         ztrdv(:,:,:) = pva(:,:,:) - ztrdv(:,:,:)
206         CALL trd_dyn( ztrdu(:,:,:), ztrdv(:,:,:), jpdyn_bfr, kt, Kmm )
207         DEALLOCATE( ztrdu, ztrdv )
208      ENDIF
209      !                                          ! print mean trends (used for debugging)
210      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab3d_1=pua, clinfo1=' bfr  - Ua: ', mask1=umask,               &
211         &                                  tab3d_2=pva, clinfo2=       ' Va: ', mask2=vmask, clinfo3='dyn' )
212      !
213   END SUBROUTINE zdf_drg_exp
214
215
216   SUBROUTINE zdf_drg_init
217      !!----------------------------------------------------------------------
218      !!                  ***  ROUTINE zdf_brg_init  ***
219      !!
220      !! ** Purpose :   Initialization of the bottom friction
221      !!
222      !! ** Method  :   Read the namdrg namelist and check their consistency
223      !!                called at the first timestep (nit000)
224      !!----------------------------------------------------------------------
225      INTEGER   ::   ji, jj      ! dummy loop indexes
226      INTEGER   ::   ios, ioptio   ! local integers
227      !!
228      NAMELIST/namdrg/ ln_drg_OFF, ln_lin, ln_non_lin, ln_loglayer, ln_drgimp
229      !!----------------------------------------------------------------------
230      !
231      !                     !==  drag nature  ==!
232      !
233      READ  ( numnam_ref, namdrg, IOSTAT = ios, ERR = 901)
234901   IF( ios /= 0 )   CALL ctl_nam( ios , 'namdrg in reference namelist' )
235      READ  ( numnam_cfg, namdrg, IOSTAT = ios, ERR = 902 )
236902   IF( ios >  0 )   CALL ctl_nam( ios , 'namdrg in configuration namelist' )
237      IF(lwm) WRITE ( numond, namdrg )
238      !
239      IF(lwp) THEN
240         WRITE(numout,*)
241         WRITE(numout,*) 'zdf_drg_init : top and/or bottom drag setting'
242         WRITE(numout,*) '~~~~~~~~~~~~'
243         WRITE(numout,*) '   Namelist namdrg : top/bottom friction choices'
244         WRITE(numout,*) '      free-slip       : Cd = 0                  ln_drg_OFF  = ', ln_drg_OFF 
245         WRITE(numout,*) '      linear  drag    : Cd = Cd0                ln_lin      = ', ln_lin
246         WRITE(numout,*) '      non-linear  drag: Cd = Cd0_nl |U|         ln_non_lin  = ', ln_non_lin
247         WRITE(numout,*) '      logarithmic drag: Cd = vkarmn/log(z/z0)   ln_loglayer = ', ln_loglayer
248         WRITE(numout,*) '      implicit friction                         ln_drgimp   = ', ln_drgimp
249      ENDIF
250      !
251      ioptio = 0                       ! set ndrg and control check
252      IF( ln_drg_OFF  ) THEN   ;   ndrg = np_OFF        ;   ioptio = ioptio + 1   ;   ENDIF
253      IF( ln_lin      ) THEN   ;   ndrg = np_lin        ;   ioptio = ioptio + 1   ;   ENDIF
254      IF( ln_non_lin  ) THEN   ;   ndrg = np_non_lin    ;   ioptio = ioptio + 1   ;   ENDIF
255      IF( ln_loglayer ) THEN   ;   ndrg = np_loglayer   ;   ioptio = ioptio + 1   ;   ENDIF
256      !
257      IF( ioptio /= 1 )   CALL ctl_stop( 'zdf_drg_init: Choose ONE type of drag coef in namdrg' )
258      !
259      !
260      !                     !==  BOTTOM drag setting  ==!   (applied at seafloor)
261      !
262      ALLOCATE( rCd0_bot(jpi,jpj), rCdU_bot(jpi,jpj) )
263      CALL drg_init( 'BOTTOM'   , mbkt       ,                                         &   ! <== in
264         &           r_Cdmin_bot, r_Cdmax_bot, r_z0_bot, r_ke0_bot, rCd0_bot, rCdU_bot )   ! ==> out
265
266      !
267      !                     !==  TOP drag setting  ==!   (applied at the top of ocean cavities)
268      !
269      IF( ln_isfcav ) THEN              ! Ocean cavities: top friction setting
270         ALLOCATE( rCd0_top(jpi,jpj), rCdU_top(jpi,jpj) )
271         CALL drg_init( 'TOP   '   , mikt       ,                                         &   ! <== in
272            &           r_Cdmin_top, r_Cdmax_top, r_z0_top, r_ke0_top, rCd0_top, rCdU_top )   ! ==> out
273      ENDIF
274      !
275   END SUBROUTINE zdf_drg_init
276
277
278   SUBROUTINE drg_init( cd_topbot, k_mk,  &
279      &                 pCdmin, pCdmax, pz0, pke0, pCd0, pCdU ) 
280      !!----------------------------------------------------------------------
281      !!                  ***  ROUTINE drg_init  ***
282      !!
283      !! ** Purpose :   Initialization of the top/bottom friction CdO and Cd
284      !!              from namelist parameters
285      !!----------------------------------------------------------------------
286      CHARACTER(len=6)        , INTENT(in   ) ::   cd_topbot       ! top/ bot indicator
287      INTEGER , DIMENSION(:,:), INTENT(in   ) ::   k_mk            ! 1st/last  wet level
288      REAL(wp)                , INTENT(  out) ::   pCdmin, pCdmax  ! min and max drag coef. [-]
289      REAL(wp)                , INTENT(  out) ::   pz0             ! roughness              [m]
290      REAL(wp)                , INTENT(  out) ::   pke0            ! background KE          [m2/s2]
291      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   pCd0            ! masked precomputed part of the non-linear drag coefficient
292      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   pCdU            ! minus linear drag*|U| at t-points  [m/s]
293      !!
294      CHARACTER(len=40) ::   cl_namdrg, cl_file, cl_varname, cl_namref, cl_namcfg  ! local names
295      INTEGER ::   ji, jj              ! dummy loop indexes
296      LOGICAL ::   ll_top, ll_bot      ! local logical
297      INTEGER ::   ios, inum, imk      ! local integers
298      REAL(wp)::   zmsk, zzz, zcd      ! local scalars
299      REAL(wp), DIMENSION(jpi,jpj) ::   zmsk_boost   ! 2D workspace
300      !!
301      NAMELIST/namdrg_top/ rn_Cd0, rn_Uc0, rn_Cdmax, rn_ke0, rn_z0, ln_boost, rn_boost
302      NAMELIST/namdrg_bot/ rn_Cd0, rn_Uc0, rn_Cdmax, rn_ke0, rn_z0, ln_boost, rn_boost
303      !!----------------------------------------------------------------------
304      !
305      !                          !==  set TOP / BOTTOM specificities  ==!
306      ll_top = .FALSE.
307      ll_bot = .FALSE.
308      !
309      SELECT CASE (cd_topbot)
310      CASE( 'TOP   ' )
311         ll_top = .TRUE.
312         cl_namdrg  = 'namdrg_top'
313         cl_namref  = 'namdrg_top in reference     namelist'
314         cl_namcfg  = 'namdrg_top in configuration namelist'
315         cl_file    = 'tfr_coef.nc'
316         cl_varname = 'tfr_coef'
317      CASE( 'BOTTOM' )
318         ll_bot = .TRUE.
319         cl_namdrg  = 'namdrg_bot'
320         cl_namref  = 'namdrg_bot  in reference     namelist'
321         cl_namcfg  = 'namdrg_bot  in configuration namelist'
322         cl_file    = 'bfr_coef.nc'
323         cl_varname = 'bfr_coef'
324      CASE DEFAULT
325         CALL ctl_stop( 'drg_init: bad value for cd_topbot ' )
326      END SELECT
327      !
328      !                          !==  read namlist  ==!
329      !
330      IF(ll_top)   READ  ( numnam_ref, namdrg_top, IOSTAT = ios, ERR = 901)
331      IF(ll_bot)   READ  ( numnam_ref, namdrg_bot, IOSTAT = ios, ERR = 901)
332901   IF( ios /= 0 )   CALL ctl_nam( ios , TRIM(cl_namref) )
333      IF(ll_top)   READ  ( numnam_cfg, namdrg_top, IOSTAT = ios, ERR = 902 )
334      IF(ll_bot)   READ  ( numnam_cfg, namdrg_bot, IOSTAT = ios, ERR = 902 )
335902   IF( ios >  0 )   CALL ctl_nam( ios , TRIM(cl_namcfg) )
336      IF(lwm .AND. ll_top)   WRITE ( numond, namdrg_top )
337      IF(lwm .AND. ll_bot)   WRITE ( numond, namdrg_bot )
338      !
339      IF(lwp) THEN
340         WRITE(numout,*)
341         WRITE(numout,*) '   Namelist ',TRIM(cl_namdrg),' : set ',TRIM(cd_topbot),' friction parameters'
342         WRITE(numout,*) '      drag coefficient                        rn_Cd0   = ', rn_Cd0
343         WRITE(numout,*) '      characteristic velocity (linear case)   rn_Uc0   = ', rn_Uc0, ' m/s'
344         WRITE(numout,*) '      non-linear drag maximum                 rn_Cdmax = ', rn_Cdmax
345         WRITE(numout,*) '      background kinetic energy  (n-l case)   rn_ke0   = ', rn_ke0
346         WRITE(numout,*) '      bottom roughness           (n-l case)   rn_z0    = ', rn_z0
347         WRITE(numout,*) '      set a regional boost of Cd0             ln_boost = ', ln_boost
348         WRITE(numout,*) '         associated boost factor              rn_boost = ', rn_boost
349      ENDIF
350      !
351      !                          !==  return some namelist parametres  ==!   (used in non_lin and loglayer cases)
352      pCdmin = rn_Cd0
353      pCdmax = rn_Cdmax
354      pz0    = rn_z0
355      pke0   = rn_ke0
356      !
357      !                          !==  mask * boost factor  ==!
358      !
359      IF( ln_boost ) THEN           !* regional boost:   boost factor = 1 + regional boost
360         IF(lwp) WRITE(numout,*)
361         IF(lwp) WRITE(numout,*) '   ==>>>   use a regional boost read in ', TRIM(cl_file), ' file'
362         IF(lwp) WRITE(numout,*) '           using enhancement factor of ', rn_boost
363         ! cl_varname is a coefficient in [0,1] giving where to apply the regional boost
364         CALL iom_open ( TRIM(cl_file), inum )
365         CALL iom_get  ( inum, jpdom_global, TRIM(cl_varname), zmsk_boost, 1 )
366         CALL iom_close( inum)
367         zmsk_boost(:,:) = 1._wp + rn_boost * zmsk_boost(:,:)
368         !
369      ELSE                          !* no boost:   boost factor = 1
370         zmsk_boost(:,:) = 1._wp
371      ENDIF
372      !                             !* mask outside ocean cavities area (top) or land area (bot)
373      IF(ll_top)   zmsk_boost(:,:) = zmsk_boost(:,:) * ssmask(:,:) * (1. - tmask(:,:,1) )  ! none zero in ocean cavities only
374      IF(ll_bot)   zmsk_boost(:,:) = zmsk_boost(:,:) * ssmask(:,:)                         ! x seafloor mask
375      !
376      !
377      SELECT CASE( ndrg )
378      !
379      CASE( np_OFF  )            !==  No top/bottom friction  ==!   (pCdU = 0)
380         IF(lwp) WRITE(numout,*)
381         IF(lwp) WRITE(numout,*) '   ==>>>   ',TRIM(cd_topbot),' free-slip, friction set to zero'
382         !
383         l_zdfdrg = .FALSE.         ! no time variation of the drag: set it one for all
384         !
385         pCdU(:,:) = 0._wp
386         pCd0(:,:) = 0._wp
387         !
388      CASE( np_lin )             !==  linear friction  ==!   (pCdU = Cd0 * Uc0)
389         IF(lwp) WRITE(numout,*)
390         IF(lwp) WRITE(numout,*) '   ==>>>   linear ',TRIM(cd_topbot),' friction (constant coef = Cd0*Uc0 = ', rn_Cd0*rn_Uc0, ')'
391         !
392         l_zdfdrg = .FALSE.         ! no time variation of the Cd*|U| : set it one for all
393         !                     
394         pCd0(:,:) = rn_Cd0 * zmsk_boost(:,:)  !* constant in time drag coefficient (= mask (and boost) Cd0)
395         pCdU(:,:) = - pCd0(:,:) * rn_Uc0      !  using a constant velocity
396         !
397      CASE( np_non_lin )         !== non-linear friction  ==!   (pCd0 = Cd0 )
398         IF(lwp) WRITE(numout,*)
399         IF(lwp) WRITE(numout,*) '   ==>>>   quadratic ',TRIM(cd_topbot),' friction (propotional to module of the velocity)'
400         IF(lwp) WRITE(numout,*) '   with    a drag coefficient Cd0 = ', rn_Cd0, ', and'
401         IF(lwp) WRITE(numout,*) '           a background velocity module of (rn_ke0)^1/2 = ', SQRT(rn_ke0), 'm/s)'
402         !
403         l_zdfdrg = .TRUE.          !* Cd*|U| updated at each time-step (it depends on ocean velocity)
404         !
405         pCd0(:,:) = rn_Cd0 * zmsk_boost(:,:)  !* constant in time proportionality coefficient (= mask (and boost) Cd0)
406         pCdU(:,:) = 0._wp                     
407         !
408      CASE( np_loglayer )       !== logarithmic layer formulation of friction  ==!   (CdU = (vkarman log(z/z0))^2 |U| )
409         IF(lwp) WRITE(numout,*)
410         IF(lwp) WRITE(numout,*) '   ==>>>   quadratic ',TRIM(cd_topbot),' drag (propotional to module of the velocity)'
411         IF(lwp) WRITE(numout,*) '   with   a logarithmic Cd0 formulation Cd0 = ( vkarman log(z/z0) )^2 ,'
412         IF(lwp) WRITE(numout,*) '          a background velocity module of (rn_ke0)^1/2 = ', SQRT(pke0), 'm/s), '
413         IF(lwp) WRITE(numout,*) '          a logarithmic formulation: a roughness of ', pz0, ' meters,   and '
414         IF(lwp) WRITE(numout,*) '          a proportionality factor bounded by min/max values of ', pCdmin, pCdmax
415         !
416         l_zdfdrg = .TRUE.          !* Cd*|U| updated at each time-step (it depends on ocean velocity)
417         !
418         IF( ln_linssh ) THEN       !* pCd0 = (v log(z/z0))^2   as velocity points have a fixed z position
419            IF(lwp) WRITE(numout,*)
420            IF(lwp) WRITE(numout,*) '   N.B.   linear free surface case, Cd0 computed one for all'
421            !
422            l_log_not_linssh = .FALSE.    !- don't update Cd at each time step
423            !
424            DO_2D( 1, 1, 1, 1 )
425               zzz =  0.5_wp * e3t_0(ji,jj,k_mk(ji,jj))
426               zcd = (  vkarmn / LOG( zzz / rn_z0 )  )**2
427               pCd0(ji,jj) = zmsk_boost(ji,jj) * MIN(  MAX( rn_Cd0 , zcd ) , rn_Cdmax  )  ! rn_Cd0 < Cd0 < rn_Cdmax
428            END_2D
429         ELSE                       !* Cd updated at each time-step ==> pCd0 = mask * boost
430            IF(lwp) WRITE(numout,*)
431            IF(lwp) WRITE(numout,*) '   N.B.   non-linear free surface case, Cd0 updated at each time-step '
432            !
433            l_log_not_linssh = .TRUE.     ! compute the drag coef. at each time-step
434            !
435            pCd0(:,:) = zmsk_boost(:,:)
436         ENDIF
437         pCdU(:,:) = 0._wp          ! initialisation to zero (will be updated at each time step)
438         !
439      CASE DEFAULT
440         CALL ctl_stop( 'drg_init: bad flag value for ndrg ' )
441      END SELECT
442      !
443   END SUBROUTINE drg_init
444
445   !!======================================================================
446END MODULE zdfdrg
Note: See TracBrowser for help on using the repository browser.