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/2020/temporary_r4_trunk/src/OCE/ZDF – NEMO

source: NEMO/branches/2020/temporary_r4_trunk/src/OCE/ZDF/zdfdrg.F90 @ 13470

Last change on this file since 13470 was 13470, checked in by smasson, 4 years ago

r4_trunk: second change of DO loops for routines to be merged, see #2523

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