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/2021/dev_r14116_HPC-10_mcastril_Mixed_Precision_implementation/src/OCE/ZDF – NEMO

source: NEMO/branches/2021/dev_r14116_HPC-10_mcastril_Mixed_Precision_implementation/src/OCE/ZDF/zdfphy.F90 @ 15660

Last change on this file since 15660 was 15540, checked in by sparonuz, 3 years ago

Mixed precision version, tested up to 30 years on ORCA2.

  • Property svn:keywords set to Id
File size: 22.8 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   ! TEMP: [tiling] This change not necessary after finalisation of zdf_osm (not yet tiled)
15   USE domtile
16   USE zdf_oce        ! vertical physics: shared variables
17   USE zdfdrg         ! vertical physics: top/bottom drag coef.
18   USE zdfsh2         ! vertical physics: shear production term of TKE
19   USE zdfric         ! vertical physics: RIChardson dependent vertical mixing
20   USE zdftke         ! vertical physics: TKE vertical mixing
21   USE zdfgls         ! vertical physics: GLS vertical mixing
22   USE zdfosm         ! vertical physics: OSMOSIS vertical mixing
23   USE zdfddm         ! vertical physics: double diffusion mixing
24   USE zdfevd         ! vertical physics: convection via enhanced vertical diffusion
25   USE zdfmfc         ! vertical physics: Mass Flux Convection
26   USE zdfiwm         ! vertical physics: internal wave-induced mixing
27   USE zdfswm         ! vertical physics: surface  wave-induced mixing
28   USE zdfmxl         ! vertical physics: mixed layer
29   USE tranpc         ! convection: non penetrative adjustment
30   USE trc_oce        ! variables shared between passive tracer & ocean
31   USE sbc_oce        ! surface module (only for nn_isf in the option compatibility test)
32   USE sbcrnf         ! surface boundary condition: runoff variables
33   USE sbc_ice        ! sea ice drag
34#if defined key_agrif
35   USE agrif_oce_interp   ! interpavm
36#endif
37   !
38   USE in_out_manager ! I/O manager
39   USE iom            ! IOM library
40   USE lbclnk         ! lateral boundary conditions
41   USE lib_mpp        ! distribued memory computing
42   USE timing         ! Timing
43
44   IMPLICIT NONE
45   PRIVATE
46
47   PUBLIC   zdf_phy_init  ! called by nemogcm.F90
48   PUBLIC   zdf_phy       ! called by step.F90
49
50   INTEGER ::   nzdf_phy   ! type of vertical closure used
51   !                       ! associated indicators
52   INTEGER, PARAMETER ::   np_CST = 1   ! Constant Kz
53   INTEGER, PARAMETER ::   np_RIC = 2   ! Richardson number dependent Kz
54   INTEGER, PARAMETER ::   np_TKE = 3   ! Turbulente Kinetic Eenergy closure scheme for Kz
55   INTEGER, PARAMETER ::   np_GLS = 4   ! Generic Length Scale closure scheme for Kz
56   INTEGER, PARAMETER ::   np_OSM = 5   ! OSMOSIS-OBL closure scheme for Kz
57
58   LOGICAL, PUBLIC ::   l_zdfsh2   ! shear production term flag (=F for CST, =T otherwise (i.e. TKE, GLS, RIC))
59
60   REAL(wp), SAVE, ALLOCATABLE, DIMENSION(:,:,:) ::   avm_k_n !: "Now" avm_k used for calculation of zsh2 with tiling
61
62   !! * Substitutions
63#  include "single_precision_substitute.h90"
64#  include "do_loop_substitute.h90"
65   !!----------------------------------------------------------------------
66   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
67   !! $Id$
68   !! Software governed by the CeCILL license (see ./LICENSE)
69   !!----------------------------------------------------------------------
70CONTAINS
71
72   SUBROUTINE zdf_phy_init( Kmm )
73      !!----------------------------------------------------------------------
74      !!                  ***  ROUTINE zdf_phy_init  ***
75      !!
76      !! ** Purpose :   initializations of the vertical ocean physics
77      !!
78      !! ** Method  :   Read namelist namzdf, control logicals
79      !!                set horizontal shape and vertical profile of background mixing coef.
80      !!----------------------------------------------------------------------
81      INTEGER, INTENT(in)    :: Kmm ! time level index (middle)
82      !
83      INTEGER ::   jk            ! dummy loop indices
84      INTEGER ::   ioptio, ios   ! local integers
85      !!
86      NAMELIST/namzdf/ ln_zdfcst, ln_zdfric, ln_zdftke, ln_zdfgls,   &     ! type of closure scheme
87         &             ln_zdfosm,                                    &     ! type of closure scheme
88         &             ln_zdfmfc,                                    &     ! convection : mass flux
89         &             ln_zdfevd, nn_evdm, rn_evd ,                  &     ! convection : evd
90         &             ln_zdfnpc, nn_npc , nn_npcp,                  &     ! convection : npc
91         &             ln_zdfddm, rn_avts, rn_hsbfr,                 &     ! double diffusion
92         &             ln_zdfswm,                                    &     ! surface  wave-induced mixing
93         &             ln_zdfiwm,                                    &     ! internal  -      -      -
94         &             ln_zad_Aimp,                                  &     ! apdative-implicit vertical advection
95         &             rn_avm0, rn_avt0, nn_avb, nn_havtb                  ! coefficients
96      !!----------------------------------------------------------------------
97      !
98      IF(lwp) THEN
99         WRITE(numout,*)
100         WRITE(numout,*) 'zdf_phy_init: ocean vertical physics'
101         WRITE(numout,*) '~~~~~~~~~~~~'
102      ENDIF
103      !
104      !                           !==  Namelist  ==!
105      READ  ( numnam_ref, namzdf, IOSTAT = ios, ERR = 901)
106901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namzdf in reference namelist' )
107      !
108      READ  ( numnam_cfg, namzdf, IOSTAT = ios, ERR = 902 )
109902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namzdf in configuration namelist' )
110      IF(lwm)   WRITE ( numond, namzdf )
111      !
112      IF(lwp) THEN                      ! Parameter print
113         WRITE(numout,*) '   Namelist namzdf : set vertical mixing mixing parameters'
114         WRITE(numout,*) '      adaptive-implicit vertical advection'
115         WRITE(numout,*) '         Courant number targeted application   ln_zad_Aimp = ', ln_zad_Aimp
116         WRITE(numout,*) '      vertical closure scheme'
117         WRITE(numout,*) '         constant vertical mixing coefficient    ln_zdfcst = ', ln_zdfcst
118         WRITE(numout,*) '         Richardson number dependent closure     ln_zdfric = ', ln_zdfric
119         WRITE(numout,*) '         Turbulent Kinetic Energy closure (TKE)  ln_zdftke = ', ln_zdftke
120         WRITE(numout,*) '         Generic Length Scale closure (GLS)      ln_zdfgls = ', ln_zdfgls
121         WRITE(numout,*) '         OSMOSIS-OBL closure (OSM)               ln_zdfosm = ', ln_zdfosm
122         WRITE(numout,*) '      convection: '
123         WRITE(numout,*) '         convection mass flux (mfc)              ln_zdfmfc = ', ln_zdfmfc
124         WRITE(numout,*) '         enhanced vertical diffusion             ln_zdfevd = ', ln_zdfevd
125         WRITE(numout,*) '            applied on momentum (=1/0)             nn_evdm = ', nn_evdm
126         WRITE(numout,*) '            vertical coefficient for evd           rn_evd  = ', rn_evd
127         WRITE(numout,*) '         non-penetrative convection (npc)        ln_zdfnpc = ', ln_zdfnpc
128         WRITE(numout,*) '            npc call  frequency                    nn_npc  = ', nn_npc
129         WRITE(numout,*) '            npc print frequency                    nn_npcp = ', nn_npcp
130         WRITE(numout,*) '      double diffusive mixing                    ln_zdfddm = ', ln_zdfddm
131         WRITE(numout,*) '         maximum avs for dd mixing                 rn_avts = ', rn_avts
132         WRITE(numout,*) '         heat/salt buoyancy flux ratio             rn_hsbfr= ', rn_hsbfr
133         WRITE(numout,*) '      gravity wave-induced mixing'
134         WRITE(numout,*) '         surface  wave (Qiao et al 2010)         ln_zdfswm = ', ln_zdfswm                                          ! surface wave induced mixing
135         WRITE(numout,*) '         internal wave (de Lavergne et al 2017)  ln_zdfiwm = ', ln_zdfiwm
136         WRITE(numout,*) '      coefficients : '
137         WRITE(numout,*) '         vertical eddy viscosity                 rn_avm0   = ', rn_avm0
138         WRITE(numout,*) '         vertical eddy diffusivity               rn_avt0   = ', rn_avt0
139         WRITE(numout,*) '         constant background or profile          nn_avb    = ', nn_avb
140         WRITE(numout,*) '         horizontal variation for avtb           nn_havtb  = ', nn_havtb
141      ENDIF
142
143      IF( ln_zad_Aimp ) THEN
144         IF( zdf_phy_alloc() /= 0 )   &
145            &       CALL ctl_stop( 'STOP', 'zdf_phy_init : unable to allocate adaptive-implicit z-advection arrays' )
146         Cu_adv(:,:,:) = 0._wp
147         wi    (:,:,:) = 0._wp
148      ENDIF
149      !                          !==  Background eddy viscosity and diffusivity  ==!
150      IF( nn_avb == 0 ) THEN             ! Define avmb, avtb from namelist parameter
151         avmb(:) = rn_avm0
152         avtb(:) = rn_avt0
153      ELSE                               ! Background profile of avt (fit a theoretical/observational profile (Krauss 1990)
154         avmb(:) = rn_avm0
155         avtb(:) = rn_avt0 + ( 3.e-4_wp - 2._wp * rn_avt0 ) * 1.e-4_wp * gdepw_1d(:)   ! m2/s
156         IF(ln_sco .AND. lwp)   CALL ctl_warn( 'avtb profile not valid in sco' )
157      ENDIF
158      !                                  ! 2D shape of the avtb
159      avtb_2d(:,:) = 1._wp                   ! uniform
160      !
161      IF( nn_havtb == 1 ) THEN               ! decrease avtb by a factor of ten in the equatorial band
162           !                                 !   -15S -5S : linear decrease from avt0 to avt0/10.
163           !                                 !   -5S  +5N : cst value avt0/10.
164           !                                 !    5N  15N : linear increase from avt0/10, to avt0
165           WHERE(-15. <= gphit .AND. gphit < -5 )   avtb_2d = (1.  - 0.09 * (gphit + 15.))
166           WHERE( -5. <= gphit .AND. gphit <  5 )   avtb_2d =  0.1
167           WHERE(  5. <= gphit .AND. gphit < 15 )   avtb_2d = (0.1 + 0.09 * (gphit -  5.))
168      ENDIF
169      !
170      DO jk = 1, jpk                      ! set turbulent closure Kz to the background value (avt_k, avm_k)
171         avt_k(:,:,jk) = avtb_2d(:,:) * avtb(jk) * wmask (:,:,jk)
172         avm_k(:,:,jk) =                avmb(jk) * wmask (:,:,jk)
173      END DO
174!!gm  to be tested only the 1st & last levels
175!      avt  (:,:, 1 ) = 0._wp   ;   avs(:,:, 1 ) = 0._wp   ;   avm  (:,:, 1 ) = 0._wp
176!      avt  (:,:,jpk) = 0._wp   ;   avs(:,:,jpk) = 0._wp   ;   avm  (:,:,jpk) = 0._wp
177!!gm
178      avt  (:,:,:) = 0._wp   ;   avs(:,:,:) = 0._wp   ;   avm  (:,:,:) = 0._wp
179
180      !                          !==  Convection  ==!
181      !
182      IF( ln_zdfnpc .AND. ln_zdfevd )   CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfnpc and ln_zdfevd' )
183      IF( ln_zdfosm .AND. ln_zdfevd )   CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfosm and ln_zdfevd' )
184      IF( ln_zdfmfc .AND. ln_zdfevd )   CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfmfc and ln_zdfevd' )
185      IF( ln_zdfmfc .AND. ln_zdfnpc )   CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfmfc and ln_zdfnpc' )
186      IF( ln_zdfmfc .AND. ln_zdfosm )   CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfmfc and ln_zdfosm' )
187      IF( lk_top    .AND. ln_zdfnpc )   CALL ctl_stop( 'zdf_phy_init: npc scheme is not working with key_top' )
188      IF( lk_top    .AND. ln_zdfosm )   CALL ctl_warn( 'zdf_phy_init: osmosis gives no non-local fluxes for TOP tracers yet' )
189      IF( lk_top    .AND. ln_zdfmfc )   CALL ctl_stop( 'zdf_phy_init: Mass Flux scheme is not working with key_top' )
190      IF(lwp) THEN
191         WRITE(numout,*)
192         IF    ( ln_zdfnpc ) THEN  ;   WRITE(numout,*) '   ==>>>   convection: use non penetrative convective scheme'
193         ELSEIF( ln_zdfevd ) THEN  ;   WRITE(numout,*) '   ==>>>   convection: use enhanced vertical diffusion scheme'
194         ELSEIF( ln_zdfmfc ) THEN  ;   WRITE(numout,*) '   ==>>>   convection: use Mass Flux scheme'
195         ELSE                      ;   WRITE(numout,*) '   ==>>>   convection: no specific scheme used'
196         ENDIF
197      ENDIF
198
199      IF(lwp) THEN               !==  Double Diffusion Mixing parameterization  ==!   (ddm)
200         WRITE(numout,*)
201         IF( ln_zdfddm ) THEN   ;   WRITE(numout,*) '   ==>>>   use double diffusive mixing: avs /= avt'
202         ELSE                   ;   WRITE(numout,*) '   ==>>>   No  double diffusive mixing: avs = avt'
203         ENDIF
204      ENDIF
205
206      !                          !==  type of vertical turbulent closure  ==!    (set nzdf_phy)
207      ioptio = 0
208      IF( ln_zdfcst ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_CST   ;   ENDIF
209      IF( ln_zdfric ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_RIC   ;   CALL zdf_ric_init          ;   ENDIF
210      IF( ln_zdftke ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_TKE   ;   CALL zdf_tke_init( Kmm )   ;   ENDIF
211      IF( ln_zdfgls ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_GLS   ;   CALL zdf_gls_init          ;   ENDIF
212      IF( ln_zdfosm ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_OSM   ;   CALL zdf_osm_init( Kmm )   ;   ENDIF
213      !
214      IF( ioptio /= 1 )    CALL ctl_stop( 'zdf_phy_init: one and only one vertical diffusion option has to be defined ' )
215      IF( ln_isfcav ) THEN
216      IF( ln_zdfric )      CALL ctl_stop( 'zdf_phy_init: zdfric never tested with ice shelves cavities ' )
217      ENDIF
218      !                                ! shear production term flag
219      IF( ln_zdfcst .OR. ln_zdfosm ) THEN   ;   l_zdfsh2 = .FALSE.
220      ELSE                                  ;   l_zdfsh2 = .TRUE.
221      ENDIF
222      IF( ln_tile .AND. l_zdfsh2 ) ALLOCATE( avm_k_n(jpi,jpj,jpk) )
223      !                          !== Mass Flux Convectiive algorithm  ==!
224      IF( ln_zdfmfc )   CALL zdf_mfc_init       ! Convection computed with eddy diffusivity mass flux
225      !
226      !                          !== gravity wave-driven mixing  ==!
227      IF( ln_zdfiwm )   CALL zdf_iwm_init       ! internal wave-driven mixing
228      IF( ln_zdfswm )   CALL zdf_swm_init       ! surface  wave-driven mixing
229
230      !                          !== top/bottom friction  ==!
231      CALL zdf_drg_init
232      !
233      !                          !== time-stepping  ==!
234      ! Check/update of time stepping done in dynzdf_init/trazdf_init
235      !!gm move it here ?
236      !
237   END SUBROUTINE zdf_phy_init
238
239
240   SUBROUTINE zdf_phy( kt, Kbb, Kmm, Krhs )
241      !!----------------------------------------------------------------------
242      !!                     ***  ROUTINE zdf_phy  ***
243      !!
244      !! ** Purpose :  Update ocean physics at each time-step
245      !!
246      !! ** Method  :
247      !!
248      !! ** Action  :   avm, avt vertical eddy viscosity and diffusivity at w-points
249      !!                nmld ??? mixed layer depth in level and meters   <<<<====verifier !
250      !!                bottom stress.....                               <<<<====verifier !
251      !!----------------------------------------------------------------------
252      INTEGER, INTENT(in) ::   kt         ! ocean time-step index
253      INTEGER, INTENT(in) ::   Kbb, Kmm, Krhs   ! ocean time level indices
254      !
255      INTEGER ::   ji, jj, jk   ! dummy loop indice
256      REAL(dp), DIMENSION(A2D(nn_hls),jpk) ::   zsh2   ! shear production
257      !! ---------------------------------------------------------------------
258      !
259      IF( ln_timing )   CALL timing_start('zdf_phy')
260      !
261      IF( l_zdfdrg ) THEN     !==  update top/bottom drag  ==!   (non-linear cases)
262         !
263         !                       !* bottom drag
264         CALL zdf_drg( kt, Kmm, mbkt , r_Cdmin_bot, r_Cdmax_bot,   &   ! <<== in
265            &              r_z0_bot,   r_ke0_bot,    rCd0_bot,   &
266            &                                        rCdU_bot  )     ! ==>> out : bottom drag [m/s]
267         IF( ln_isfcav ) THEN    !* top drag   (ocean cavities)
268            CALL zdf_drg( kt, Kmm, mikt , r_Cdmin_top, r_Cdmax_top,   &   ! <<== in
269               &              r_z0_top,   r_ke0_top,    rCd0_top,   &
270               &                                        rCdU_top  )     ! ==>> out : bottom drag [m/s]
271         ENDIF
272      ENDIF
273      !
274#if defined key_si3
275      IF ( ln_drgice_imp) THEN
276         IF ( ln_isfcav ) THEN
277            DO_2D_OVR( 1, 1, 1, 1 )
278               rCdU_top(ji,jj) = rCdU_top(ji,jj) + ssmask(ji,jj) * tmask(ji,jj,1) * rCdU_ice(ji,jj)
279            END_2D
280         ELSE
281            DO_2D_OVR( 1, 1, 1, 1 )
282               rCdU_top(ji,jj) = rCdU_ice(ji,jj)
283            END_2D
284         ENDIF
285      ENDIF
286#endif
287      !
288      CALL zdf_mxl( kt, Kmm )                        !* mixed layer depth, and level
289      !
290      !                       !==  Kz from chosen turbulent closure  ==!   (avm_k, avt_k)
291      !
292      ! NOTE: [tiling] the closure schemes (zdf_tke etc) will update avm_k. With tiling, the calculation of zsh2 on adjacent tiles then uses both updated (next timestep) and non-updated (current timestep) values of avm_k. To preserve results, we save a read-only copy of the "now" avm_k to use in the calculation of zsh2.
293      IF( l_zdfsh2 ) THEN        !* shear production at w-points (energy conserving form)
294         IF( ln_tile ) THEN
295            IF( ntile == 1 ) avm_k_n(:,:,:) = avm_k(:,:,:)     ! Preserve "now" avm_k for calculation of zsh2
296            CALL zdf_sh2( Kbb, Kmm, avm_k_n, &     ! <<== in
297               &                     zsh2    )     ! ==>> out : shear production
298         ELSE
299            CALL zdf_sh2( Kbb, Kmm, CASTSP(avm_k),   &     ! <<== in
300               &                     zsh2    )     ! ==>> out : shear production
301         ENDIF
302      ENDIF
303      !
304      SELECT CASE ( nzdf_phy )                  !* Vertical eddy viscosity and diffusivity coefficients at w-points
305      CASE( np_RIC )   ;   CALL zdf_ric( kt,      Kmm, zsh2, avm_k, avt_k )    ! Richardson number dependent Kz
306      CASE( np_TKE )   ;   CALL zdf_tke( kt, Kbb, Kmm, zsh2, avm_k, avt_k )    ! TKE closure scheme for Kz
307      CASE( np_GLS )   ;   CALL zdf_gls( kt, Kbb, Kmm, zsh2, avm_k, avt_k )    ! GLS closure scheme for Kz
308      CASE( np_OSM )   ;   CALL zdf_osm( kt, Kbb, Kmm, Krhs, avm_k, avt_k )    ! OSMOSIS closure scheme for Kz
309   !     CASE( np_CST )                                  ! Constant Kz (reset avt, avm to the background value)
310   !         ! avt_k and avm_k set one for all at initialisation phase
311!!gm         avt(2:jpim1,2:jpjm1,1:jpkm1) = rn_avt0 * wmask(2:jpim1,2:jpjm1,1:jpkm1)
312!!gm         avm(2:jpim1,2:jpjm1,1:jpkm1) = rn_avm0 * wmask(2:jpim1,2:jpjm1,1:jpkm1)
313      END SELECT
314      !
315      !                          !==  ocean Kz  ==!   (avt, avs, avm)
316      !
317      !                                         !* start from turbulent closure values
318      DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
319         avt(ji,jj,jk) = avt_k(ji,jj,jk)
320         avm(ji,jj,jk) = avm_k(ji,jj,jk)
321      END_3D
322      !
323      IF( ln_rnf_mouth ) THEN                   !* increase diffusivity at rivers mouths
324         DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, nkrnf )
325            avt(ji,jj,jk) = avt(ji,jj,jk) + 2._wp * rn_avt_rnf * rnfmsk(ji,jj) * wmask(ji,jj,jk)
326         END_3D
327      ENDIF
328      !
329      IF( ln_zdfevd )   CALL zdf_evd( kt, Kmm, Krhs, avm, avt )  !* convection: enhanced vertical eddy diffusivity
330      !
331      !                                         !* double diffusive mixing
332      IF( ln_zdfddm ) THEN                            ! update avt and compute avs
333                        CALL zdf_ddm( kt, Kmm,  avm, avt, avs )
334      ELSE                                            ! same mixing on all tracers
335         DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 )
336            avs(ji,jj,jk) = avt(ji,jj,jk)
337         END_3D
338      ENDIF
339      !
340      !                                         !* wave-induced mixing
341      IF( ln_zdfswm )   CALL zdf_swm( kt, Kmm, avm, avt, avs )   ! surface  wave (Qiao et al. 2004)
342      IF( ln_zdfiwm )   CALL zdf_iwm( kt, Kmm, avm, avt, avs )   ! internal wave (de Lavergne et al 2017)
343
344#if defined key_agrif
345      ! interpolation parent grid => child grid for avm_k ( ex : at west border: update column 1 and 2)
346      IF( .NOT. l_istiled .OR. ntile == nijtile ) THEN                       ! Do only on the last tile
347         IF( l_zdfsh2 )   CALL Agrif_avm
348      ENDIF
349#endif
350
351      !                                         !* Lateral boundary conditions (sign unchanged)
352      IF(nn_hls==1) THEN
353         IF( l_zdfsh2 ) THEN
354            CALL lbc_lnk( 'zdfphy', avm_k, 'W', 1.0_wp , avt_k, 'W', 1.0_wp,   &
355                  &                 avm  , 'W', 1.0_wp , avt  , 'W', 1.0_wp , avs , 'W', 1.0_wp )
356         ELSE
357            CALL lbc_lnk( 'zdfphy', avm  , 'W', 1.0_wp , avt  , 'W', 1.0_wp , avs , 'W', 1.0_wp )
358         ENDIF
359         !
360         IF( l_zdfdrg ) THEN     ! drag  have been updated (non-linear cases)
361            IF( ln_isfcav ) THEN   ;  CALL lbc_lnk( 'zdfphy', rCdU_top, 'T', 1.0_wp , rCdU_bot, 'T', 1.0_wp )   ! top & bot drag
362            ELSE                   ;  CALL lbc_lnk( 'zdfphy', rCdU_bot, 'T', 1.0_wp )                           ! bottom drag only
363            ENDIF
364         ENDIF
365      ENDIF
366      !
367      CALL zdf_mxl_turb( kt, Kmm )                   !* turbocline depth
368      !
369      IF( .NOT. l_istiled .OR. ntile == nijtile ) THEN                       ! Do only on the last tile
370         IF( lrst_oce ) THEN                       !* write TKE, GLS or RIC fields in the restart file
371            IF( ln_zdftke )   CALL tke_rst( kt, 'WRITE' )
372            IF( ln_zdfgls )   CALL gls_rst( kt, 'WRITE' )
373            IF( ln_zdfric )   CALL ric_rst( kt, 'WRITE' )
374            ! NB. OSMOSIS restart (osm_rst) will be called in step.F90 after ww has been updated
375         ENDIF
376      ENDIF
377      !
378      ! diagnostics of energy dissipation
379      IF( iom_use('avt_k') .OR. iom_use('avm_k') .OR. iom_use('eshear_k') .OR. iom_use('estrat_k') ) THEN
380         IF( l_zdfsh2 ) THEN
381            zsh2(:,:,1  ) = 0._wp
382            zsh2(:,:,jpk) = 0._wp
383            CALL lbc_lnk( 'zdfphy', zsh2, 'W', 1._wp )
384            CALL iom_put( 'avt_k'   ,   avt_k       * wmask )
385            CALL iom_put( 'avm_k'   ,   avm_k       * wmask )
386            CALL iom_put( 'eshear_k',   zsh2        * wmask )
387            CALL iom_put( 'estrat_k', - avt_k * rn2 * wmask )
388         ENDIF
389      ENDIF
390      !
391      IF( ln_timing )   CALL timing_stop('zdf_phy')
392      !
393   END SUBROUTINE zdf_phy
394
395
396   INTEGER FUNCTION zdf_phy_alloc()
397      !!----------------------------------------------------------------------
398      !!                 ***  FUNCTION zdf_phy_alloc  ***
399      !!----------------------------------------------------------------------
400     ! Allocate wi array (declared in oce.F90) for use with the adaptive-implicit vertical velocity option
401     ALLOCATE(     wi(jpi,jpj,jpk), Cu_adv(jpi,jpj,jpk),  STAT= zdf_phy_alloc )
402     IF( zdf_phy_alloc /= 0 )   CALL ctl_warn('zdf_phy_alloc: failed to allocate ln_zad_Aimp=T required arrays')
403     CALL mpp_sum ( 'zdfphy', zdf_phy_alloc )
404   END FUNCTION zdf_phy_alloc
405
406   !!======================================================================
407END MODULE zdfphy
Note: See TracBrowser for help on using the repository browser.