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