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.
zdfphy.F90 in NEMO/branches/UKMO/dev_r9922_GC3_cpl_pkg/src/OCE/ZDF – NEMO

source: NEMO/branches/UKMO/dev_r9922_GC3_cpl_pkg/src/OCE/ZDF/zdfphy.F90 @ 9947

Last change on this file since 9947 was 9947, checked in by timgraham, 6 years ago

Clear svn keywords

File size: 18.1 KB
Line 
1MODULE zdfphy
2   !!======================================================================
3   !!                      ***  MODULE  zdfphy  ***
4   !! Vertical ocean physics :   manager of all vertical physics packages
5   !!======================================================================
6   !! History :  4.0  !  2017-04  (G. Madec)  original code
7   !!----------------------------------------------------------------------
8
9   !!----------------------------------------------------------------------
10   !!   zdf_phy_init  : initialization of all vertical physics packages
11   !!   zdf_phy       : upadate at each time-step the vertical mixing coeff.
12   !!----------------------------------------------------------------------
13   USE oce            ! ocean dynamics and tracers variables
14   USE zdf_oce        ! vertical physics: shared variables         
15   USE zdfdrg         ! vertical physics: top/bottom drag coef.
16   USE zdfsh2         ! vertical physics: shear production term of TKE
17   USE zdfric         ! vertical physics: RIChardson dependent vertical mixing   
18   USE zdftke         ! vertical physics: TKE vertical mixing
19   USE zdfgls         ! vertical physics: GLS vertical mixing
20   USE zdfosm         ! vertical physics: OSMOSIS vertical mixing
21   USE zdfddm         ! vertical physics: double diffusion mixing     
22   USE zdfevd         ! vertical physics: convection via enhanced vertical diffusion 
23   USE zdfiwm         ! vertical physics: internal wave-induced mixing 
24   USE zdfswm         ! vertical physics: surface  wave-induced mixing
25   USE zdfmxl         ! vertical physics: mixed layer
26   USE tranpc         ! convection: non penetrative adjustment
27   USE trc_oce        ! variables shared between passive tracer & ocean           
28   USE sbc_oce        ! surface module (only for nn_isf in the option compatibility test)
29   USE sbcrnf         ! surface boundary condition: runoff variables
30#if defined key_agrif
31   USE agrif_oce_interp   ! interpavm
32#endif
33   !
34   USE in_out_manager ! I/O manager
35   USE iom            ! IOM library
36   USE lbclnk         ! lateral boundary conditions
37   USE lib_mpp        ! distribued memory computing
38   USE timing         ! Timing
39
40   IMPLICIT NONE
41   PRIVATE
42
43   PUBLIC   zdf_phy_init  ! called by nemogcm.F90
44   PUBLIC   zdf_phy       ! called by step.F90
45
46   INTEGER ::   nzdf_phy   ! type of vertical closure used
47   !                       ! associated indicators
48   INTEGER, PARAMETER ::   np_CST = 1   ! Constant Kz
49   INTEGER, PARAMETER ::   np_RIC = 2   ! Richardson number dependent Kz
50   INTEGER, PARAMETER ::   np_TKE = 3   ! Turbulente Kinetic Eenergy closure scheme for Kz
51   INTEGER, PARAMETER ::   np_GLS = 4   ! Generic Length Scale closure scheme for Kz
52   INTEGER, PARAMETER ::   np_OSM = 5   ! OSMOSIS-OBL closure scheme for Kz
53
54   LOGICAL ::   l_zdfsh2   ! shear production term flag (=F for CST, =T otherwise (i.e. TKE, GLS, RIC))
55
56   !!----------------------------------------------------------------------
57   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
58   !! $Id$
59   !! Software governed by the CeCILL licence     (./LICENSE)
60   !!----------------------------------------------------------------------
61CONTAINS
62
63   SUBROUTINE zdf_phy_init
64      !!----------------------------------------------------------------------
65      !!                  ***  ROUTINE zdf_phy_init  ***
66      !!
67      !! ** Purpose :   initializations of the vertical ocean physics
68      !!
69      !! ** Method  :   Read namelist namzdf, control logicals
70      !!                set horizontal shape and vertical profile of background mixing coef.
71      !!----------------------------------------------------------------------
72      INTEGER ::   jk            ! dummy loop indices
73      INTEGER ::   ioptio, ios   ! local integers
74      !!
75      NAMELIST/namzdf/ ln_zdfcst, ln_zdfric, ln_zdftke, ln_zdfgls,   &     ! type of closure scheme
76         &             ln_zdfosm,                                    &     ! type of closure scheme
77         &             ln_zdfevd, nn_evdm, rn_evd ,                  &     ! convection : evd
78         &             ln_zdfnpc, nn_npc , nn_npcp,                  &     ! convection : npc
79         &             ln_zdfddm, rn_avts, rn_hsbfr,                 &     ! double diffusion
80         &             ln_zdfswm,                                    &     ! surface  wave-induced mixing
81         &             ln_zdfiwm,                                    &     ! internal  -      -      -
82         &             rn_avm0, rn_avt0, nn_avb, nn_havtb                  ! coefficients
83      !!----------------------------------------------------------------------
84      !
85      IF(lwp) THEN
86         WRITE(numout,*)
87         WRITE(numout,*) 'zdf_phy_init: ocean vertical physics'
88         WRITE(numout,*) '~~~~~~~~~~~~'
89      ENDIF
90      !
91      !                           !==  Namelist  ==!
92      REWIND( numnam_ref )              ! Namelist namzdf in reference namelist : Vertical mixing parameters
93      READ  ( numnam_ref, namzdf, IOSTAT = ios, ERR = 901)
94901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namzdf in reference namelist', lwp )
95      !
96      REWIND( numnam_cfg )              ! Namelist namzdf in reference namelist : Vertical mixing parameters
97      READ  ( numnam_cfg, namzdf, IOSTAT = ios, ERR = 902 )
98902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namzdf in configuration namelist', lwp )
99      IF(lwm)   WRITE ( numond, namzdf )
100      !
101      IF(lwp) THEN                      ! Parameter print
102         WRITE(numout,*) '   Namelist namzdf : set vertical mixing mixing parameters'
103         WRITE(numout,*) '      vertical closure scheme'
104         WRITE(numout,*) '         constant vertical mixing coefficient    ln_zdfcst = ', ln_zdfcst
105         WRITE(numout,*) '         Richardson number dependent closure     ln_zdfric = ', ln_zdfric
106         WRITE(numout,*) '         Turbulent Kinetic Energy closure (TKE)  ln_zdftke = ', ln_zdftke
107         WRITE(numout,*) '         Generic Length Scale closure (GLS)      ln_zdfgls = ', ln_zdfgls
108         WRITE(numout,*) '         OSMOSIS-OBL closure (OSM)               ln_zdfosm = ', ln_zdfosm
109         WRITE(numout,*) '      convection: '
110         WRITE(numout,*) '         enhanced vertical diffusion             ln_zdfevd = ', ln_zdfevd
111         WRITE(numout,*) '            applied on momentum (=1/0)             nn_evdm = ', nn_evdm
112         WRITE(numout,*) '            vertical coefficient for evd           rn_evd  = ', rn_evd
113         WRITE(numout,*) '         non-penetrative convection (npc)        ln_zdfnpc = ', ln_zdfnpc
114         WRITE(numout,*) '            npc call  frequency                    nn_npc  = ', nn_npc
115         WRITE(numout,*) '            npc print frequency                    nn_npcp = ', nn_npcp
116         WRITE(numout,*) '      double diffusive mixing                    ln_zdfddm = ', ln_zdfddm
117         WRITE(numout,*) '         maximum avs for dd mixing                 rn_avts = ', rn_avts
118         WRITE(numout,*) '         heat/salt buoyancy flux ratio             rn_hsbfr= ', rn_hsbfr
119         WRITE(numout,*) '      gravity wave-induced mixing'
120         WRITE(numout,*) '         surface  wave (Qiao et al 2010)         ln_zdfswm = ', ln_zdfswm                                          ! surface wave induced mixing
121         WRITE(numout,*) '         internal wave (de Lavergne et al 2017)  ln_zdfiwm = ', ln_zdfiwm
122         WRITE(numout,*) '      coefficients : '
123         WRITE(numout,*) '         vertical eddy viscosity                 rn_avm0   = ', rn_avm0
124         WRITE(numout,*) '         vertical eddy diffusivity               rn_avt0   = ', rn_avt0
125         WRITE(numout,*) '         constant background or profile          nn_avb    = ', nn_avb
126         WRITE(numout,*) '         horizontal variation for avtb           nn_havtb  = ', nn_havtb
127      ENDIF
128
129      !                          !==  Background eddy viscosity and diffusivity  ==!
130      IF( nn_avb == 0 ) THEN             ! Define avmb, avtb from namelist parameter
131         avmb(:) = rn_avm0
132         avtb(:) = rn_avt0                     
133      ELSE                               ! Background profile of avt (fit a theoretical/observational profile (Krauss 1990)
134         avmb(:) = rn_avm0
135         avtb(:) = rn_avt0 + ( 3.e-4_wp - 2._wp * rn_avt0 ) * 1.e-4_wp * gdepw_1d(:)   ! m2/s
136         IF(ln_sco .AND. lwp)   CALL ctl_warn( 'avtb profile not valid in sco' )
137      ENDIF
138      !                                  ! 2D shape of the avtb
139      avtb_2d(:,:) = 1._wp                   ! uniform
140      !
141      IF( nn_havtb == 1 ) THEN               ! decrease avtb by a factor of ten in the equatorial band
142           !                                 !   -15S -5S : linear decrease from avt0 to avt0/10.
143           !                                 !   -5S  +5N : cst value avt0/10.
144           !                                 !    5N  15N : linear increase from avt0/10, to avt0
145           WHERE(-15. <= gphit .AND. gphit < -5 )   avtb_2d = (1.  - 0.09 * (gphit + 15.))
146           WHERE( -5. <= gphit .AND. gphit <  5 )   avtb_2d =  0.1
147           WHERE(  5. <= gphit .AND. gphit < 15 )   avtb_2d = (0.1 + 0.09 * (gphit -  5.))
148      ENDIF
149      !
150      DO jk = 1, jpk                      ! set turbulent closure Kz to the background value (avt_k, avm_k)
151         avt_k(:,:,jk) = avtb_2d(:,:) * avtb(jk) * wmask (:,:,jk)
152         avm_k(:,:,jk) =                avmb(jk) * wmask (:,:,jk)
153      END DO
154!!gm  to be tested only the 1st & last levels
155!      avt  (:,:, 1 ) = 0._wp   ;   avs(:,:, 1 ) = 0._wp   ;   avm  (:,:, 1 ) = 0._wp
156!      avt  (:,:,jpk) = 0._wp   ;   avs(:,:,jpk) = 0._wp   ;   avm  (:,:,jpk) = 0._wp
157!!gm
158      avt  (:,:,:) = 0._wp   ;   avs(:,:,:) = 0._wp   ;   avm  (:,:,:) = 0._wp
159
160      !                          !==  Convection  ==!
161      !
162      IF( ln_zdfnpc .AND. ln_zdfevd )   CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfnpc and ln_zdfevd' )
163      IF( ln_zdfosm .AND. ln_zdfevd )   CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfosm and ln_zdfevd' )
164      IF( lk_top    .AND. ln_zdfnpc )   CALL ctl_stop( 'zdf_phy_init: npc scheme is not working with key_top' )
165      IF( lk_top    .AND. ln_zdfosm )   CALL ctl_stop( 'zdf_phy_init: osmosis scheme is not working with key_top' )
166      IF(lwp) THEN
167         WRITE(numout,*)
168         IF    ( ln_zdfnpc ) THEN  ;   WRITE(numout,*) '   ==>>>   convection: use non penetrative convective scheme'
169         ELSEIF( ln_zdfevd ) THEN  ;   WRITE(numout,*) '   ==>>>   convection: use enhanced vertical diffusion scheme'
170         ELSE                      ;   WRITE(numout,*) '   ==>>>   convection: no specific scheme used'
171         ENDIF
172      ENDIF
173
174      IF(lwp) THEN               !==  Double Diffusion Mixing parameterization  ==!   (ddm)
175         WRITE(numout,*)
176         IF( ln_zdfddm ) THEN   ;   WRITE(numout,*) '   ==>>>   use double diffusive mixing: avs /= avt'
177         ELSE                   ;   WRITE(numout,*) '   ==>>>   No  double diffusive mixing: avs = avt'
178         ENDIF
179      ENDIF
180
181      !                          !==  type of vertical turbulent closure  ==!    (set nzdf_phy)
182      ioptio = 0 
183      IF( ln_zdfcst ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_CST   ;   ENDIF
184      IF( ln_zdfric ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_RIC   ;   CALL zdf_ric_init   ;   ENDIF
185      IF( ln_zdftke ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_TKE   ;   CALL zdf_tke_init   ;   ENDIF
186      IF( ln_zdfgls ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_GLS   ;   CALL zdf_gls_init   ;   ENDIF
187      IF( ln_zdfosm ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_OSM   ;   CALL zdf_osm_init   ;   ENDIF
188      !
189      IF( ioptio /= 1 )    CALL ctl_stop( 'zdf_phy_init: one and only one vertical diffusion option has to be defined ' )
190      IF( ln_isfcav ) THEN
191      IF( ln_zdfric .OR. ln_zdfgls )    CALL ctl_stop( 'zdf_phy_init: zdfric and zdfgls never tested with ice shelves cavities ' )
192      ENDIF
193      !                                ! shear production term flag
194      IF( ln_zdfcst ) THEN   ;   l_zdfsh2 = .FALSE.
195      ELSE                   ;   l_zdfsh2 = .TRUE.
196      ENDIF
197
198      !                          !== gravity wave-driven mixing  ==!
199      IF( ln_zdfiwm )   CALL zdf_iwm_init       ! internal wave-driven mixing
200      IF( ln_zdfswm )   CALL zdf_swm_init       ! surface  wave-driven mixing
201
202      !                          !== top/bottom friction  ==!
203      CALL zdf_drg_init
204      !
205      !                          !== time-stepping  ==!
206      ! Check/update of time stepping done in dynzdf_init/trazdf_init
207      !!gm move it here ?
208      !
209   END SUBROUTINE zdf_phy_init
210
211
212   SUBROUTINE zdf_phy( kt )
213      !!----------------------------------------------------------------------
214      !!                     ***  ROUTINE zdf_phy  ***
215      !!
216      !! ** Purpose :  Update ocean physics at each time-step
217      !!
218      !! ** Method  :
219      !!
220      !! ** Action  :   avm, avt vertical eddy viscosity and diffusivity at w-points
221      !!                nmld ??? mixed layer depth in level and meters   <<<<====verifier !
222      !!                bottom stress.....                               <<<<====verifier !
223      !!----------------------------------------------------------------------
224      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
225      !
226      INTEGER ::   ji, jj, jk   ! dummy loop indice
227      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zsh2   ! shear production
228      !! ---------------------------------------------------------------------
229      !
230      IF( ln_timing )   CALL timing_start('zdf_phy')
231      !
232      IF( l_zdfdrg ) THEN     !==  update top/bottom drag  ==!   (non-linear cases)
233         !
234         !                       !* bottom drag
235         CALL zdf_drg( kt, mbkt    , r_Cdmin_bot, r_Cdmax_bot,   &   ! <<== in
236            &              r_z0_bot,   r_ke0_bot,    rCd0_bot,   &
237            &                                        rCdU_bot  )     ! ==>> out : bottom drag [m/s]
238         IF( ln_isfcav ) THEN    !* top drag   (ocean cavities)
239            CALL zdf_drg( kt, mikt    , r_Cdmin_top, r_Cdmax_top,   &   ! <<== in
240               &              r_z0_top,   r_ke0_top,    rCd0_top,   &
241               &                                        rCdU_top  )     ! ==>> out : bottom drag [m/s]
242         ENDIF
243      ENDIF
244      !
245      !                       !==  Kz from chosen turbulent closure  ==!   (avm_k, avt_k)
246      !
247      IF( l_zdfsh2 )   &         !* shear production at w-points (energy conserving form)
248         CALL zdf_sh2( ub, vb, un, vn, avm_k,   &     ! <<== in
249            &                           zsh2    )     ! ==>> out : shear production
250      !
251      SELECT CASE ( nzdf_phy )                  !* Vertical eddy viscosity and diffusivity coefficients at w-points
252      CASE( np_RIC )   ;   CALL zdf_ric( kt, gdept_n, zsh2, avm_k, avt_k )    ! Richardson number dependent Kz
253      CASE( np_TKE )   ;   CALL zdf_tke( kt         , zsh2, avm_k, avt_k )    ! TKE closure scheme for Kz
254      CASE( np_GLS )   ;   CALL zdf_gls( kt         , zsh2, avm_k, avt_k )    ! GLS closure scheme for Kz
255      CASE( np_OSM )   ;   CALL zdf_osm( kt               , avm_k, avt_k )    ! OSMOSIS closure scheme for Kz
256!     CASE( np_CST )                                  ! Constant Kz (reset avt, avm to the background value)
257!         ! avt_k and avm_k set one for all at initialisation phase
258!!gm         avt(2:jpim1,2:jpjm1,1:jpkm1) = rn_avt0 * wmask(2:jpim1,2:jpjm1,1:jpkm1)
259!!gm         avm(2:jpim1,2:jpjm1,1:jpkm1) = rn_avm0 * wmask(2:jpim1,2:jpjm1,1:jpkm1)
260      END SELECT
261     
262      !                          !==  ocean Kz  ==!   (avt, avs, avm)
263      !
264      !                                         !* start from turbulent closure values
265      avt(:,:,2:jpkm1) = avt_k(:,:,2:jpkm1)
266      avm(:,:,2:jpkm1) = avm_k(:,:,2:jpkm1)
267      !
268      IF( ln_rnf_mouth ) THEN                   !* increase diffusivity at rivers mouths
269         DO jk = 2, nkrnf
270            avt(:,:,jk) = avt(:,:,jk) + 2._wp * rn_avt_rnf * rnfmsk(:,:) * wmask(:,:,jk)
271         END DO
272      ENDIF
273      !
274      IF( ln_zdfevd )   CALL zdf_evd( kt, avm, avt )  !* convection: enhanced vertical eddy diffusivity
275      !
276      !                                         !* double diffusive mixing
277      IF( ln_zdfddm ) THEN                            ! update avt and compute avs
278                        CALL zdf_ddm( kt, avm, avt, avs )
279      ELSE                                            ! same mixing on all tracers
280         avs(2:jpim1,2:jpjm1,1:jpkm1) = avt(2:jpim1,2:jpjm1,1:jpkm1)
281      ENDIF
282      !
283      !                                         !* wave-induced mixing
284      IF( ln_zdfswm )   CALL zdf_swm( kt, avm, avt, avs )   ! surface  wave (Qiao et al. 2004)
285      IF( ln_zdfiwm )   CALL zdf_iwm( kt, avm, avt, avs )   ! internal wave (de Lavergne et al 2017)
286
287#if defined key_agrif 
288      ! interpolation parent grid => child grid for avm_k ( ex : at west border: update column 1 and 2)
289      IF( l_zdfsh2 )   CALL Agrif_avm
290#endif
291
292      !                                         !* Lateral boundary conditions (sign unchanged)
293      IF( l_zdfsh2 ) THEN
294         CALL lbc_lnk_multi( avm_k, 'W', 1. , avt_k, 'W', 1.,   &
295            &                avm  , 'W', 1. , avt  , 'W', 1. , avs , 'W', 1. )
296      ELSE
297         CALL lbc_lnk_multi( avm  , 'W', 1. , avt  , 'W', 1. , avs , 'W', 1. )
298      ENDIF
299      !
300      IF( l_zdfdrg ) THEN     ! drag  have been updated (non-linear cases)
301         IF( ln_isfcav ) THEN   ;  CALL lbc_lnk_multi( rCdU_top, 'T', 1. , rCdU_bot, 'T', 1. )   ! top & bot drag
302         ELSE                   ;  CALL lbc_lnk      ( rCdU_bot, 'T', 1. )                       ! bottom drag only
303         ENDIF
304      ENDIF
305      !
306      CALL zdf_mxl( kt )                        !* mixed layer depth, and level
307      !
308      IF( lrst_oce ) THEN                       !* write TKE, GLS or RIC fields in the restart file
309         IF( ln_zdftke )   CALL tke_rst( kt, 'WRITE' )
310         IF( ln_zdfgls )   CALL gls_rst( kt, 'WRITE' )
311         IF( ln_zdfric )   CALL ric_rst( kt, 'WRITE' ) 
312         ! NB. OSMOSIS restart (osm_rst) will be called in step.F90 after wn has been updated
313      ENDIF
314      !
315      IF( ln_timing )   CALL timing_stop('zdf_phy')
316      !
317   END SUBROUTINE zdf_phy
318
319   !!======================================================================
320END MODULE zdfphy
Note: See TracBrowser for help on using the repository browser.