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/2020/dev_r13648_ASINTER-04_laurent_bulk_ice/src/OCE/ZDF – NEMO

source: NEMO/branches/2020/dev_r13648_ASINTER-04_laurent_bulk_ice/src/OCE/ZDF/zdfphy.F90 @ 14021

Last change on this file since 14021 was 14021, checked in by laurent, 3 years ago

Caught up with trunk rev 14020...

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