1 | MODULE 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 license (see ./LICENSE) |
---|
60 | !!---------------------------------------------------------------------- |
---|
61 | CONTAINS |
---|
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 | & ln_zad_Aimp, & ! apdative-implicit vertical advection |
---|
83 | & rn_avm0, rn_avt0, nn_avb, nn_havtb ! coefficients |
---|
84 | !!---------------------------------------------------------------------- |
---|
85 | ! |
---|
86 | IF(lwp) THEN |
---|
87 | WRITE(numout,*) |
---|
88 | WRITE(numout,*) 'zdf_phy_init: ocean vertical physics' |
---|
89 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
90 | ENDIF |
---|
91 | ! |
---|
92 | ! !== Namelist ==! |
---|
93 | REWIND( numnam_ref ) ! Namelist namzdf in reference namelist : Vertical mixing parameters |
---|
94 | READ ( numnam_ref, namzdf, IOSTAT = ios, ERR = 901) |
---|
95 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf in reference namelist', lwp ) |
---|
96 | ! |
---|
97 | REWIND( numnam_cfg ) ! Namelist namzdf in reference namelist : Vertical mixing parameters |
---|
98 | READ ( numnam_cfg, namzdf, IOSTAT = ios, ERR = 902 ) |
---|
99 | 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'namzdf in configuration namelist', lwp ) |
---|
100 | IF(lwm) WRITE ( numond, namzdf ) |
---|
101 | ! |
---|
102 | IF(lwp) THEN ! Parameter print |
---|
103 | WRITE(numout,*) ' Namelist namzdf : set vertical mixing mixing parameters' |
---|
104 | WRITE(numout,*) ' adaptive-implicit vertical advection' |
---|
105 | WRITE(numout,*) ' Courant number targeted application ln_zad_Aimp = ', ln_zad_Aimp |
---|
106 | WRITE(numout,*) ' vertical closure scheme' |
---|
107 | WRITE(numout,*) ' constant vertical mixing coefficient ln_zdfcst = ', ln_zdfcst |
---|
108 | WRITE(numout,*) ' Richardson number dependent closure ln_zdfric = ', ln_zdfric |
---|
109 | WRITE(numout,*) ' Turbulent Kinetic Energy closure (TKE) ln_zdftke = ', ln_zdftke |
---|
110 | WRITE(numout,*) ' Generic Length Scale closure (GLS) ln_zdfgls = ', ln_zdfgls |
---|
111 | WRITE(numout,*) ' OSMOSIS-OBL closure (OSM) ln_zdfosm = ', ln_zdfosm |
---|
112 | WRITE(numout,*) ' convection: ' |
---|
113 | WRITE(numout,*) ' enhanced vertical diffusion ln_zdfevd = ', ln_zdfevd |
---|
114 | WRITE(numout,*) ' applied on momentum (=1/0) nn_evdm = ', nn_evdm |
---|
115 | WRITE(numout,*) ' vertical coefficient for evd rn_evd = ', rn_evd |
---|
116 | WRITE(numout,*) ' non-penetrative convection (npc) ln_zdfnpc = ', ln_zdfnpc |
---|
117 | WRITE(numout,*) ' npc call frequency nn_npc = ', nn_npc |
---|
118 | WRITE(numout,*) ' npc print frequency nn_npcp = ', nn_npcp |
---|
119 | WRITE(numout,*) ' double diffusive mixing ln_zdfddm = ', ln_zdfddm |
---|
120 | WRITE(numout,*) ' maximum avs for dd mixing rn_avts = ', rn_avts |
---|
121 | WRITE(numout,*) ' heat/salt buoyancy flux ratio rn_hsbfr= ', rn_hsbfr |
---|
122 | WRITE(numout,*) ' gravity wave-induced mixing' |
---|
123 | WRITE(numout,*) ' surface wave (Qiao et al 2010) ln_zdfswm = ', ln_zdfswm ! surface wave induced mixing |
---|
124 | WRITE(numout,*) ' internal wave (de Lavergne et al 2017) ln_zdfiwm = ', ln_zdfiwm |
---|
125 | WRITE(numout,*) ' coefficients : ' |
---|
126 | WRITE(numout,*) ' vertical eddy viscosity rn_avm0 = ', rn_avm0 |
---|
127 | WRITE(numout,*) ' vertical eddy diffusivity rn_avt0 = ', rn_avt0 |
---|
128 | WRITE(numout,*) ' constant background or profile nn_avb = ', nn_avb |
---|
129 | WRITE(numout,*) ' horizontal variation for avtb nn_havtb = ', nn_havtb |
---|
130 | ENDIF |
---|
131 | |
---|
132 | IF( ln_zad_Aimp ) THEN |
---|
133 | IF( zdf_phy_alloc() /= 0 ) & |
---|
134 | & CALL ctl_stop( 'STOP', 'zdf_phy_init : unable to allocate adaptive-implicit z-advection arrays' ) |
---|
135 | wi(:,:,:) = 0._wp |
---|
136 | ENDIF |
---|
137 | ! !== Background eddy viscosity and diffusivity ==! |
---|
138 | IF( nn_avb == 0 ) THEN ! Define avmb, avtb from namelist parameter |
---|
139 | avmb(:) = rn_avm0 |
---|
140 | avtb(:) = rn_avt0 |
---|
141 | ELSE ! Background profile of avt (fit a theoretical/observational profile (Krauss 1990) |
---|
142 | avmb(:) = rn_avm0 |
---|
143 | avtb(:) = rn_avt0 + ( 3.e-4_wp - 2._wp * rn_avt0 ) * 1.e-4_wp * gdepw_1d(:) ! m2/s |
---|
144 | IF(ln_sco .AND. lwp) CALL ctl_warn( 'avtb profile not valid in sco' ) |
---|
145 | ENDIF |
---|
146 | ! ! 2D shape of the avtb |
---|
147 | avtb_2d(:,:) = 1._wp ! uniform |
---|
148 | ! |
---|
149 | IF( nn_havtb == 1 ) THEN ! decrease avtb by a factor of ten in the equatorial band |
---|
150 | ! ! -15S -5S : linear decrease from avt0 to avt0/10. |
---|
151 | ! ! -5S +5N : cst value avt0/10. |
---|
152 | ! ! 5N 15N : linear increase from avt0/10, to avt0 |
---|
153 | WHERE(-15. <= gphit .AND. gphit < -5 ) avtb_2d = (1. - 0.09 * (gphit + 15.)) |
---|
154 | WHERE( -5. <= gphit .AND. gphit < 5 ) avtb_2d = 0.1 |
---|
155 | WHERE( 5. <= gphit .AND. gphit < 15 ) avtb_2d = (0.1 + 0.09 * (gphit - 5.)) |
---|
156 | ENDIF |
---|
157 | ! |
---|
158 | DO jk = 1, jpk ! set turbulent closure Kz to the background value (avt_k, avm_k) |
---|
159 | avt_k(:,:,jk) = avtb_2d(:,:) * avtb(jk) * wmask (:,:,jk) |
---|
160 | avm_k(:,:,jk) = avmb(jk) * wmask (:,:,jk) |
---|
161 | END DO |
---|
162 | !!gm to be tested only the 1st & last levels |
---|
163 | ! avt (:,:, 1 ) = 0._wp ; avs(:,:, 1 ) = 0._wp ; avm (:,:, 1 ) = 0._wp |
---|
164 | ! avt (:,:,jpk) = 0._wp ; avs(:,:,jpk) = 0._wp ; avm (:,:,jpk) = 0._wp |
---|
165 | !!gm |
---|
166 | avt (:,:,:) = 0._wp ; avs(:,:,:) = 0._wp ; avm (:,:,:) = 0._wp |
---|
167 | |
---|
168 | ! !== Convection ==! |
---|
169 | ! |
---|
170 | IF( ln_zdfnpc .AND. ln_zdfevd ) CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfnpc and ln_zdfevd' ) |
---|
171 | IF( ln_zdfosm .AND. ln_zdfevd ) CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfosm and ln_zdfevd' ) |
---|
172 | IF( lk_top .AND. ln_zdfnpc ) CALL ctl_stop( 'zdf_phy_init: npc scheme is not working with key_top' ) |
---|
173 | IF( lk_top .AND. ln_zdfosm ) CALL ctl_stop( 'zdf_phy_init: osmosis scheme is not working with key_top' ) |
---|
174 | IF(lwp) THEN |
---|
175 | WRITE(numout,*) |
---|
176 | IF ( ln_zdfnpc ) THEN ; WRITE(numout,*) ' ==>>> convection: use non penetrative convective scheme' |
---|
177 | ELSEIF( ln_zdfevd ) THEN ; WRITE(numout,*) ' ==>>> convection: use enhanced vertical diffusion scheme' |
---|
178 | ELSE ; WRITE(numout,*) ' ==>>> convection: no specific scheme used' |
---|
179 | ENDIF |
---|
180 | ENDIF |
---|
181 | |
---|
182 | IF(lwp) THEN !== Double Diffusion Mixing parameterization ==! (ddm) |
---|
183 | WRITE(numout,*) |
---|
184 | IF( ln_zdfddm ) THEN ; WRITE(numout,*) ' ==>>> use double diffusive mixing: avs /= avt' |
---|
185 | ELSE ; WRITE(numout,*) ' ==>>> No double diffusive mixing: avs = avt' |
---|
186 | ENDIF |
---|
187 | ENDIF |
---|
188 | |
---|
189 | ! !== type of vertical turbulent closure ==! (set nzdf_phy) |
---|
190 | ioptio = 0 |
---|
191 | IF( ln_zdfcst ) THEN ; ioptio = ioptio + 1 ; nzdf_phy = np_CST ; ENDIF |
---|
192 | IF( ln_zdfric ) THEN ; ioptio = ioptio + 1 ; nzdf_phy = np_RIC ; CALL zdf_ric_init ; ENDIF |
---|
193 | IF( ln_zdftke ) THEN ; ioptio = ioptio + 1 ; nzdf_phy = np_TKE ; CALL zdf_tke_init ; ENDIF |
---|
194 | IF( ln_zdfgls ) THEN ; ioptio = ioptio + 1 ; nzdf_phy = np_GLS ; CALL zdf_gls_init ; ENDIF |
---|
195 | IF( ln_zdfosm ) THEN ; ioptio = ioptio + 1 ; nzdf_phy = np_OSM ; CALL zdf_osm_init ; ENDIF |
---|
196 | ! |
---|
197 | IF( ioptio /= 1 ) CALL ctl_stop( 'zdf_phy_init: one and only one vertical diffusion option has to be defined ' ) |
---|
198 | IF( ln_isfcav ) THEN |
---|
199 | IF( ln_zdfric .OR. ln_zdfgls ) CALL ctl_stop( 'zdf_phy_init: zdfric and zdfgls never tested with ice shelves cavities ' ) |
---|
200 | ENDIF |
---|
201 | ! ! shear production term flag |
---|
202 | IF( ln_zdfcst ) THEN ; l_zdfsh2 = .FALSE. |
---|
203 | ELSE ; l_zdfsh2 = .TRUE. |
---|
204 | ENDIF |
---|
205 | |
---|
206 | ! !== gravity wave-driven mixing ==! |
---|
207 | IF( ln_zdfiwm ) CALL zdf_iwm_init ! internal wave-driven mixing |
---|
208 | IF( ln_zdfswm ) CALL zdf_swm_init ! surface wave-driven mixing |
---|
209 | |
---|
210 | ! !== top/bottom friction ==! |
---|
211 | CALL zdf_drg_init |
---|
212 | ! |
---|
213 | ! !== time-stepping ==! |
---|
214 | ! Check/update of time stepping done in dynzdf_init/trazdf_init |
---|
215 | !!gm move it here ? |
---|
216 | ! |
---|
217 | END SUBROUTINE zdf_phy_init |
---|
218 | |
---|
219 | |
---|
220 | SUBROUTINE zdf_phy( kt ) |
---|
221 | !!---------------------------------------------------------------------- |
---|
222 | !! *** ROUTINE zdf_phy *** |
---|
223 | !! |
---|
224 | !! ** Purpose : Update ocean physics at each time-step |
---|
225 | !! |
---|
226 | !! ** Method : |
---|
227 | !! |
---|
228 | !! ** Action : avm, avt vertical eddy viscosity and diffusivity at w-points |
---|
229 | !! nmld ??? mixed layer depth in level and meters <<<<====verifier ! |
---|
230 | !! bottom stress..... <<<<====verifier ! |
---|
231 | !!---------------------------------------------------------------------- |
---|
232 | INTEGER, INTENT(in) :: kt ! ocean time-step index |
---|
233 | ! |
---|
234 | INTEGER :: ji, jj, jk ! dummy loop indice |
---|
235 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: zsh2 ! shear production |
---|
236 | !! --------------------------------------------------------------------- |
---|
237 | ! |
---|
238 | IF( ln_timing ) CALL timing_start('zdf_phy') |
---|
239 | ! |
---|
240 | IF( l_zdfdrg ) THEN !== update top/bottom drag ==! (non-linear cases) |
---|
241 | ! |
---|
242 | ! !* bottom drag |
---|
243 | CALL zdf_drg( kt, mbkt , r_Cdmin_bot, r_Cdmax_bot, & ! <<== in |
---|
244 | & r_z0_bot, r_ke0_bot, rCd0_bot, & |
---|
245 | & rCdU_bot ) ! ==>> out : bottom drag [m/s] |
---|
246 | IF( ln_isfcav ) THEN !* top drag (ocean cavities) |
---|
247 | CALL zdf_drg( kt, mikt , r_Cdmin_top, r_Cdmax_top, & ! <<== in |
---|
248 | & r_z0_top, r_ke0_top, rCd0_top, & |
---|
249 | & rCdU_top ) ! ==>> out : bottom drag [m/s] |
---|
250 | ENDIF |
---|
251 | ENDIF |
---|
252 | ! |
---|
253 | ! !== Kz from chosen turbulent closure ==! (avm_k, avt_k) |
---|
254 | ! |
---|
255 | IF( l_zdfsh2 ) & !* shear production at w-points (energy conserving form) |
---|
256 | CALL zdf_sh2( ub, vb, un, vn, avm_k, & ! <<== in |
---|
257 | & zsh2 ) ! ==>> out : shear production |
---|
258 | ! |
---|
259 | SELECT CASE ( nzdf_phy ) !* Vertical eddy viscosity and diffusivity coefficients at w-points |
---|
260 | CASE( np_RIC ) ; CALL zdf_ric( kt, gdept_n, zsh2, avm_k, avt_k ) ! Richardson number dependent Kz |
---|
261 | CASE( np_TKE ) ; CALL zdf_tke( kt , zsh2, avm_k, avt_k ) ! TKE closure scheme for Kz |
---|
262 | CASE( np_GLS ) ; CALL zdf_gls( kt , zsh2, avm_k, avt_k ) ! GLS closure scheme for Kz |
---|
263 | CASE( np_OSM ) ; CALL zdf_osm( kt , avm_k, avt_k ) ! OSMOSIS closure scheme for Kz |
---|
264 | ! CASE( np_CST ) ! Constant Kz (reset avt, avm to the background value) |
---|
265 | ! ! avt_k and avm_k set one for all at initialisation phase |
---|
266 | !!gm avt(2:jpim1,2:jpjm1,1:jpkm1) = rn_avt0 * wmask(2:jpim1,2:jpjm1,1:jpkm1) |
---|
267 | !!gm avm(2:jpim1,2:jpjm1,1:jpkm1) = rn_avm0 * wmask(2:jpim1,2:jpjm1,1:jpkm1) |
---|
268 | END SELECT |
---|
269 | ! |
---|
270 | ! !== ocean Kz ==! (avt, avs, avm) |
---|
271 | ! |
---|
272 | ! !* start from turbulent closure values |
---|
273 | avt(:,:,2:jpkm1) = avt_k(:,:,2:jpkm1) |
---|
274 | avm(:,:,2:jpkm1) = avm_k(:,:,2:jpkm1) |
---|
275 | ! |
---|
276 | IF( ln_rnf_mouth ) THEN !* increase diffusivity at rivers mouths |
---|
277 | DO jk = 2, nkrnf |
---|
278 | avt(:,:,jk) = avt(:,:,jk) + 2._wp * rn_avt_rnf * rnfmsk(:,:) * wmask(:,:,jk) |
---|
279 | END DO |
---|
280 | ENDIF |
---|
281 | ! |
---|
282 | IF( ln_zdfevd ) CALL zdf_evd( kt, avm, avt ) !* convection: enhanced vertical eddy diffusivity |
---|
283 | ! |
---|
284 | ! !* double diffusive mixing |
---|
285 | IF( ln_zdfddm ) THEN ! update avt and compute avs |
---|
286 | CALL zdf_ddm( kt, avm, avt, avs ) |
---|
287 | ELSE ! same mixing on all tracers |
---|
288 | avs(2:jpim1,2:jpjm1,1:jpkm1) = avt(2:jpim1,2:jpjm1,1:jpkm1) |
---|
289 | ENDIF |
---|
290 | ! |
---|
291 | ! !* wave-induced mixing |
---|
292 | IF( ln_zdfswm ) CALL zdf_swm( kt, avm, avt, avs ) ! surface wave (Qiao et al. 2004) |
---|
293 | IF( ln_zdfiwm ) CALL zdf_iwm( kt, avm, avt, avs ) ! internal wave (de Lavergne et al 2017) |
---|
294 | |
---|
295 | #if defined key_agrif |
---|
296 | ! interpolation parent grid => child grid for avm_k ( ex : at west border: update column 1 and 2) |
---|
297 | IF( l_zdfsh2 ) CALL Agrif_avm |
---|
298 | #endif |
---|
299 | |
---|
300 | ! !* Lateral boundary conditions (sign unchanged) |
---|
301 | IF( l_zdfsh2 ) THEN |
---|
302 | CALL lbc_lnk_multi( 'zdfphy', avm_k, 'W', 1. , avt_k, 'W', 1., & |
---|
303 | & avm , 'W', 1. , avt , 'W', 1. , avs , 'W', 1. ) |
---|
304 | ELSE |
---|
305 | CALL lbc_lnk_multi( 'zdfphy', avm , 'W', 1. , avt , 'W', 1. , avs , 'W', 1. ) |
---|
306 | ENDIF |
---|
307 | ! |
---|
308 | IF( l_zdfdrg ) THEN ! drag have been updated (non-linear cases) |
---|
309 | IF( ln_isfcav ) THEN ; CALL lbc_lnk_multi( 'zdfphy', rCdU_top, 'T', 1. , rCdU_bot, 'T', 1. ) ! top & bot drag |
---|
310 | ELSE ; CALL lbc_lnk ( 'zdfphy', rCdU_bot, 'T', 1. ) ! bottom drag only |
---|
311 | ENDIF |
---|
312 | ENDIF |
---|
313 | ! |
---|
314 | CALL zdf_mxl( kt ) !* mixed layer depth, and level |
---|
315 | ! |
---|
316 | IF( lrst_oce ) THEN !* write TKE, GLS or RIC fields in the restart file |
---|
317 | IF( ln_zdftke ) CALL tke_rst( kt, 'WRITE' ) |
---|
318 | IF( ln_zdfgls ) CALL gls_rst( kt, 'WRITE' ) |
---|
319 | IF( ln_zdfric ) CALL ric_rst( kt, 'WRITE' ) |
---|
320 | ! NB. OSMOSIS restart (osm_rst) will be called in step.F90 after wn has been updated |
---|
321 | ENDIF |
---|
322 | ! |
---|
323 | IF( ln_timing ) CALL timing_stop('zdf_phy') |
---|
324 | ! |
---|
325 | END SUBROUTINE zdf_phy |
---|
326 | INTEGER FUNCTION zdf_phy_alloc() |
---|
327 | !!---------------------------------------------------------------------- |
---|
328 | !! *** FUNCTION zdf_phy_alloc *** |
---|
329 | !!---------------------------------------------------------------------- |
---|
330 | ! Allocate wi array (declared in oce.F90) for use with the adaptive-implicit vertical velocity option |
---|
331 | ALLOCATE( wi(jpi,jpj,jpk), Cu_adv(jpi,jpj,jpk), STAT= zdf_phy_alloc ) |
---|
332 | IF( zdf_phy_alloc /= 0 ) CALL ctl_warn('zdf_phy_alloc: failed to allocate ln_zad_Aimp=T required arrays') |
---|
333 | CALL mpp_sum ( 'zdfphy', zdf_phy_alloc ) |
---|
334 | END FUNCTION zdf_phy_alloc |
---|
335 | |
---|
336 | !!====================================================================== |
---|
337 | END MODULE zdfphy |
---|