[8215] | 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 |
---|
[14072] | 11 | !! zdf_phy : upadate at each time-step the vertical mixing coeff. |
---|
[8215] | 12 | !!---------------------------------------------------------------------- |
---|
| 13 | USE oce ! ocean dynamics and tracers variables |
---|
[14834] | 14 | ! TEMP: [tiling] This change not necessary after finalisation of zdf_osm (not yet tiled) |
---|
| 15 | USE domtile |
---|
[14072] | 16 | USE zdf_oce ! vertical physics: shared variables |
---|
[8215] | 17 | USE zdfdrg ! vertical physics: top/bottom drag coef. |
---|
| 18 | USE zdfsh2 ! vertical physics: shear production term of TKE |
---|
[14072] | 19 | USE zdfric ! vertical physics: RIChardson dependent vertical mixing |
---|
[8215] | 20 | USE zdftke ! vertical physics: TKE vertical mixing |
---|
| 21 | USE zdfgls ! vertical physics: GLS vertical mixing |
---|
[8930] | 22 | USE zdfosm ! vertical physics: OSMOSIS vertical mixing |
---|
[14072] | 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 |
---|
[8215] | 27 | USE zdfswm ! vertical physics: surface wave-induced mixing |
---|
| 28 | USE zdfmxl ! vertical physics: mixed layer |
---|
| 29 | USE tranpc ! convection: non penetrative adjustment |
---|
[14072] | 30 | USE trc_oce ! variables shared between passive tracer & ocean |
---|
[8215] | 31 | USE sbc_oce ! surface module (only for nn_isf in the option compatibility test) |
---|
| 32 | USE sbcrnf ! surface boundary condition: runoff variables |
---|
[13472] | 33 | USE sbc_ice ! sea ice drag |
---|
[8863] | 34 | #if defined key_agrif |
---|
[9570] | 35 | USE agrif_oce_interp ! interpavm |
---|
[8863] | 36 | #endif |
---|
[8215] | 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 |
---|
[8568] | 42 | USE timing ! Timing |
---|
[8215] | 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 | |
---|
[14072] | 50 | INTEGER :: nzdf_phy ! type of vertical closure used |
---|
[8215] | 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 |
---|
[8930] | 56 | INTEGER, PARAMETER :: np_OSM = 5 ! OSMOSIS-OBL closure scheme for Kz |
---|
[8215] | 57 | |
---|
[14834] | 58 | LOGICAL, PUBLIC :: l_zdfsh2 ! shear production term flag (=F for CST, =T otherwise (i.e. TKE, GLS, RIC)) |
---|
[8215] | 59 | |
---|
[14834] | 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" |
---|
[8215] | 64 | !!---------------------------------------------------------------------- |
---|
[9598] | 65 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
[10069] | 66 | !! $Id$ |
---|
[10068] | 67 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
[8215] | 68 | !!---------------------------------------------------------------------- |
---|
| 69 | CONTAINS |
---|
| 70 | |
---|
[12377] | 71 | SUBROUTINE zdf_phy_init( Kmm ) |
---|
[8215] | 72 | !!---------------------------------------------------------------------- |
---|
| 73 | !! *** ROUTINE zdf_phy_init *** |
---|
[14072] | 74 | !! |
---|
[8215] | 75 | !! ** Purpose : initializations of the vertical ocean physics |
---|
| 76 | !! |
---|
[14072] | 77 | !! ** Method : Read namelist namzdf, control logicals |
---|
[8215] | 78 | !! set horizontal shape and vertical profile of background mixing coef. |
---|
| 79 | !!---------------------------------------------------------------------- |
---|
[12377] | 80 | INTEGER, INTENT(in) :: Kmm ! time level index (middle) |
---|
| 81 | ! |
---|
[8215] | 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 |
---|
[8930] | 86 | & ln_zdfosm, & ! type of closure scheme |
---|
[14010] | 87 | & ln_zdfmfc, & ! convection : mass flux |
---|
[8215] | 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 - - - |
---|
[10364] | 93 | & ln_zad_Aimp, & ! apdative-implicit vertical advection |
---|
[8215] | 94 | & rn_avm0, rn_avt0, nn_avb, nn_havtb ! coefficients |
---|
| 95 | !!---------------------------------------------------------------------- |
---|
| 96 | ! |
---|
[9169] | 97 | IF(lwp) THEN |
---|
| 98 | WRITE(numout,*) |
---|
| 99 | WRITE(numout,*) 'zdf_phy_init: ocean vertical physics' |
---|
| 100 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
| 101 | ENDIF |
---|
| 102 | ! |
---|
[8215] | 103 | ! !== Namelist ==! |
---|
| 104 | READ ( numnam_ref, namzdf, IOSTAT = ios, ERR = 901) |
---|
[11536] | 105 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf in reference namelist' ) |
---|
[8215] | 106 | ! |
---|
| 107 | READ ( numnam_cfg, namzdf, IOSTAT = ios, ERR = 902 ) |
---|
[11536] | 108 | 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'namzdf in configuration namelist' ) |
---|
[8215] | 109 | IF(lwm) WRITE ( numond, namzdf ) |
---|
| 110 | ! |
---|
| 111 | IF(lwp) THEN ! Parameter print |
---|
| 112 | WRITE(numout,*) ' Namelist namzdf : set vertical mixing mixing parameters' |
---|
[10364] | 113 | WRITE(numout,*) ' adaptive-implicit vertical advection' |
---|
| 114 | WRITE(numout,*) ' Courant number targeted application ln_zad_Aimp = ', ln_zad_Aimp |
---|
[8215] | 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 |
---|
[8930] | 120 | WRITE(numout,*) ' OSMOSIS-OBL closure (OSM) ln_zdfosm = ', ln_zdfosm |
---|
[8215] | 121 | WRITE(numout,*) ' convection: ' |
---|
[14010] | 122 | WRITE(numout,*) ' convection mass flux (mfc) ln_zdfmfc = ', ln_zdfmfc |
---|
[8215] | 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 | |
---|
[10364] | 142 | IF( ln_zad_Aimp ) THEN |
---|
| 143 | IF( zdf_phy_alloc() /= 0 ) & |
---|
[10907] | 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 |
---|
[10364] | 147 | ENDIF |
---|
[15249] | 148 | ! ! Initialise zdf_mxl arrays (only hmld as not set everywhere when nn_hls > 1) |
---|
| 149 | IF( zdf_mxl_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'zdf_mxl : unable to allocate arrays' ) |
---|
| 150 | hmld(:,:) = 0._wp |
---|
[8215] | 151 | ! !== Background eddy viscosity and diffusivity ==! |
---|
| 152 | IF( nn_avb == 0 ) THEN ! Define avmb, avtb from namelist parameter |
---|
| 153 | avmb(:) = rn_avm0 |
---|
[14072] | 154 | avtb(:) = rn_avt0 |
---|
[8215] | 155 | ELSE ! Background profile of avt (fit a theoretical/observational profile (Krauss 1990) |
---|
| 156 | avmb(:) = rn_avm0 |
---|
| 157 | avtb(:) = rn_avt0 + ( 3.e-4_wp - 2._wp * rn_avt0 ) * 1.e-4_wp * gdepw_1d(:) ! m2/s |
---|
| 158 | IF(ln_sco .AND. lwp) CALL ctl_warn( 'avtb profile not valid in sco' ) |
---|
| 159 | ENDIF |
---|
| 160 | ! ! 2D shape of the avtb |
---|
[14072] | 161 | avtb_2d(:,:) = 1._wp ! uniform |
---|
[8215] | 162 | ! |
---|
| 163 | IF( nn_havtb == 1 ) THEN ! decrease avtb by a factor of ten in the equatorial band |
---|
| 164 | ! ! -15S -5S : linear decrease from avt0 to avt0/10. |
---|
| 165 | ! ! -5S +5N : cst value avt0/10. |
---|
| 166 | ! ! 5N 15N : linear increase from avt0/10, to avt0 |
---|
| 167 | WHERE(-15. <= gphit .AND. gphit < -5 ) avtb_2d = (1. - 0.09 * (gphit + 15.)) |
---|
| 168 | WHERE( -5. <= gphit .AND. gphit < 5 ) avtb_2d = 0.1 |
---|
| 169 | WHERE( 5. <= gphit .AND. gphit < 15 ) avtb_2d = (0.1 + 0.09 * (gphit - 5.)) |
---|
| 170 | ENDIF |
---|
| 171 | ! |
---|
| 172 | DO jk = 1, jpk ! set turbulent closure Kz to the background value (avt_k, avm_k) |
---|
| 173 | avt_k(:,:,jk) = avtb_2d(:,:) * avtb(jk) * wmask (:,:,jk) |
---|
| 174 | avm_k(:,:,jk) = avmb(jk) * wmask (:,:,jk) |
---|
| 175 | END DO |
---|
| 176 | !!gm to be tested only the 1st & last levels |
---|
| 177 | ! avt (:,:, 1 ) = 0._wp ; avs(:,:, 1 ) = 0._wp ; avm (:,:, 1 ) = 0._wp |
---|
| 178 | ! avt (:,:,jpk) = 0._wp ; avs(:,:,jpk) = 0._wp ; avm (:,:,jpk) = 0._wp |
---|
| 179 | !!gm |
---|
| 180 | avt (:,:,:) = 0._wp ; avs(:,:,:) = 0._wp ; avm (:,:,:) = 0._wp |
---|
| 181 | |
---|
| 182 | ! !== Convection ==! |
---|
| 183 | ! |
---|
| 184 | IF( ln_zdfnpc .AND. ln_zdfevd ) CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfnpc and ln_zdfevd' ) |
---|
[8930] | 185 | IF( ln_zdfosm .AND. ln_zdfevd ) CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfosm and ln_zdfevd' ) |
---|
[14010] | 186 | IF( ln_zdfmfc .AND. ln_zdfevd ) CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfmfc and ln_zdfevd' ) |
---|
| 187 | IF( ln_zdfmfc .AND. ln_zdfnpc ) CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfmfc and ln_zdfnpc' ) |
---|
| 188 | IF( ln_zdfmfc .AND. ln_zdfosm ) CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfmfc and ln_zdfosm' ) |
---|
[8215] | 189 | IF( lk_top .AND. ln_zdfnpc ) CALL ctl_stop( 'zdf_phy_init: npc scheme is not working with key_top' ) |
---|
[14045] | 190 | IF( lk_top .AND. ln_zdfosm ) CALL ctl_warn( 'zdf_phy_init: osmosis gives no non-local fluxes for TOP tracers yet' ) |
---|
[14010] | 191 | IF( lk_top .AND. ln_zdfmfc ) CALL ctl_stop( 'zdf_phy_init: Mass Flux scheme is not working with key_top' ) |
---|
[8215] | 192 | IF(lwp) THEN |
---|
| 193 | WRITE(numout,*) |
---|
[9169] | 194 | IF ( ln_zdfnpc ) THEN ; WRITE(numout,*) ' ==>>> convection: use non penetrative convective scheme' |
---|
| 195 | ELSEIF( ln_zdfevd ) THEN ; WRITE(numout,*) ' ==>>> convection: use enhanced vertical diffusion scheme' |
---|
[14010] | 196 | ELSEIF( ln_zdfmfc ) THEN ; WRITE(numout,*) ' ==>>> convection: use Mass Flux scheme' |
---|
[9169] | 197 | ELSE ; WRITE(numout,*) ' ==>>> convection: no specific scheme used' |
---|
[8215] | 198 | ENDIF |
---|
| 199 | ENDIF |
---|
| 200 | |
---|
| 201 | IF(lwp) THEN !== Double Diffusion Mixing parameterization ==! (ddm) |
---|
| 202 | WRITE(numout,*) |
---|
[9169] | 203 | IF( ln_zdfddm ) THEN ; WRITE(numout,*) ' ==>>> use double diffusive mixing: avs /= avt' |
---|
| 204 | ELSE ; WRITE(numout,*) ' ==>>> No double diffusive mixing: avs = avt' |
---|
[8215] | 205 | ENDIF |
---|
| 206 | ENDIF |
---|
| 207 | |
---|
| 208 | ! !== type of vertical turbulent closure ==! (set nzdf_phy) |
---|
[14072] | 209 | ioptio = 0 |
---|
[8215] | 210 | IF( ln_zdfcst ) THEN ; ioptio = ioptio + 1 ; nzdf_phy = np_CST ; ENDIF |
---|
[12377] | 211 | IF( ln_zdfric ) THEN ; ioptio = ioptio + 1 ; nzdf_phy = np_RIC ; CALL zdf_ric_init ; ENDIF |
---|
| 212 | IF( ln_zdftke ) THEN ; ioptio = ioptio + 1 ; nzdf_phy = np_TKE ; CALL zdf_tke_init( Kmm ) ; ENDIF |
---|
| 213 | IF( ln_zdfgls ) THEN ; ioptio = ioptio + 1 ; nzdf_phy = np_GLS ; CALL zdf_gls_init ; ENDIF |
---|
| 214 | IF( ln_zdfosm ) THEN ; ioptio = ioptio + 1 ; nzdf_phy = np_OSM ; CALL zdf_osm_init( Kmm ) ; ENDIF |
---|
[8215] | 215 | ! |
---|
| 216 | IF( ioptio /= 1 ) CALL ctl_stop( 'zdf_phy_init: one and only one vertical diffusion option has to be defined ' ) |
---|
| 217 | IF( ln_isfcav ) THEN |
---|
[14967] | 218 | IF( ln_zdfric ) CALL ctl_stop( 'zdf_phy_init: zdfric never tested with ice shelves cavities ' ) |
---|
[8215] | 219 | ENDIF |
---|
| 220 | ! ! shear production term flag |
---|
[14834] | 221 | IF( ln_zdfcst .OR. ln_zdfosm ) THEN ; l_zdfsh2 = .FALSE. |
---|
| 222 | ELSE ; l_zdfsh2 = .TRUE. |
---|
[8215] | 223 | ENDIF |
---|
[14834] | 224 | IF( ln_tile .AND. l_zdfsh2 ) ALLOCATE( avm_k_n(jpi,jpj,jpk) ) |
---|
[14010] | 225 | ! !== Mass Flux Convectiive algorithm ==! |
---|
| 226 | IF( ln_zdfmfc ) CALL zdf_mfc_init ! Convection computed with eddy diffusivity mass flux |
---|
| 227 | ! |
---|
[8215] | 228 | ! !== gravity wave-driven mixing ==! |
---|
| 229 | IF( ln_zdfiwm ) CALL zdf_iwm_init ! internal wave-driven mixing |
---|
| 230 | IF( ln_zdfswm ) CALL zdf_swm_init ! surface wave-driven mixing |
---|
| 231 | |
---|
| 232 | ! !== top/bottom friction ==! |
---|
| 233 | CALL zdf_drg_init |
---|
| 234 | ! |
---|
| 235 | ! !== time-stepping ==! |
---|
| 236 | ! Check/update of time stepping done in dynzdf_init/trazdf_init |
---|
| 237 | !!gm move it here ? |
---|
| 238 | ! |
---|
| 239 | END SUBROUTINE zdf_phy_init |
---|
| 240 | |
---|
| 241 | |
---|
[12377] | 242 | SUBROUTINE zdf_phy( kt, Kbb, Kmm, Krhs ) |
---|
[8215] | 243 | !!---------------------------------------------------------------------- |
---|
| 244 | !! *** ROUTINE zdf_phy *** |
---|
| 245 | !! |
---|
| 246 | !! ** Purpose : Update ocean physics at each time-step |
---|
| 247 | !! |
---|
[14072] | 248 | !! ** Method : |
---|
[8215] | 249 | !! |
---|
| 250 | !! ** Action : avm, avt vertical eddy viscosity and diffusivity at w-points |
---|
| 251 | !! nmld ??? mixed layer depth in level and meters <<<<====verifier ! |
---|
| 252 | !! bottom stress..... <<<<====verifier ! |
---|
| 253 | !!---------------------------------------------------------------------- |
---|
[12377] | 254 | INTEGER, INTENT(in) :: kt ! ocean time-step index |
---|
| 255 | INTEGER, INTENT(in) :: Kbb, Kmm, Krhs ! ocean time level indices |
---|
[8215] | 256 | ! |
---|
| 257 | INTEGER :: ji, jj, jk ! dummy loop indice |
---|
[14834] | 258 | REAL(wp), DIMENSION(A2D(nn_hls),jpk) :: zsh2 ! shear production |
---|
[8215] | 259 | !! --------------------------------------------------------------------- |
---|
| 260 | ! |
---|
[8568] | 261 | IF( ln_timing ) CALL timing_start('zdf_phy') |
---|
| 262 | ! |
---|
[8215] | 263 | IF( l_zdfdrg ) THEN !== update top/bottom drag ==! (non-linear cases) |
---|
| 264 | ! |
---|
| 265 | ! !* bottom drag |
---|
[14072] | 266 | CALL zdf_drg( kt, Kmm, mbkt , r_Cdmin_bot, r_Cdmax_bot, & ! <<== in |
---|
[8215] | 267 | & r_z0_bot, r_ke0_bot, rCd0_bot, & |
---|
| 268 | & rCdU_bot ) ! ==>> out : bottom drag [m/s] |
---|
| 269 | IF( ln_isfcav ) THEN !* top drag (ocean cavities) |
---|
[14072] | 270 | CALL zdf_drg( kt, Kmm, mikt , r_Cdmin_top, r_Cdmax_top, & ! <<== in |
---|
[8215] | 271 | & r_z0_top, r_ke0_top, rCd0_top, & |
---|
| 272 | & rCdU_top ) ! ==>> out : bottom drag [m/s] |
---|
| 273 | ENDIF |
---|
| 274 | ENDIF |
---|
| 275 | ! |
---|
[13472] | 276 | #if defined key_si3 |
---|
| 277 | IF ( ln_drgice_imp) THEN |
---|
| 278 | IF ( ln_isfcav ) THEN |
---|
[14834] | 279 | DO_2D_OVR( 1, 1, 1, 1 ) |
---|
| 280 | rCdU_top(ji,jj) = rCdU_top(ji,jj) + ssmask(ji,jj) * tmask(ji,jj,1) * rCdU_ice(ji,jj) |
---|
| 281 | END_2D |
---|
[13472] | 282 | ELSE |
---|
[14834] | 283 | DO_2D_OVR( 1, 1, 1, 1 ) |
---|
| 284 | rCdU_top(ji,jj) = rCdU_ice(ji,jj) |
---|
| 285 | END_2D |
---|
[13472] | 286 | ENDIF |
---|
| 287 | ENDIF |
---|
| 288 | #endif |
---|
[14072] | 289 | ! |
---|
[14834] | 290 | CALL zdf_mxl( kt, Kmm ) !* mixed layer depth, and level |
---|
[14921] | 291 | ! |
---|
| 292 | ! !== Kz from chosen turbulent closure ==! (avm_k, avt_k) |
---|
| 293 | ! |
---|
| 294 | ! 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. |
---|
| 295 | IF( l_zdfsh2 ) THEN !* shear production at w-points (energy conserving form) |
---|
| 296 | IF( ln_tile ) THEN |
---|
| 297 | IF( ntile == 1 ) avm_k_n(:,:,:) = avm_k(:,:,:) ! Preserve "now" avm_k for calculation of zsh2 |
---|
| 298 | CALL zdf_sh2( Kbb, Kmm, avm_k_n, & ! <<== in |
---|
| 299 | & zsh2 ) ! ==>> out : shear production |
---|
| 300 | ELSE |
---|
| 301 | CALL zdf_sh2( Kbb, Kmm, avm_k, & ! <<== in |
---|
| 302 | & zsh2 ) ! ==>> out : shear production |
---|
[14834] | 303 | ENDIF |
---|
[14921] | 304 | ENDIF |
---|
| 305 | ! |
---|
| 306 | SELECT CASE ( nzdf_phy ) !* Vertical eddy viscosity and diffusivity coefficients at w-points |
---|
| 307 | CASE( np_RIC ) ; CALL zdf_ric( kt, Kmm, zsh2, avm_k, avt_k ) ! Richardson number dependent Kz |
---|
| 308 | CASE( np_TKE ) ; CALL zdf_tke( kt, Kbb, Kmm, zsh2, avm_k, avt_k ) ! TKE closure scheme for Kz |
---|
| 309 | CASE( np_GLS ) ; CALL zdf_gls( kt, Kbb, Kmm, zsh2, avm_k, avt_k ) ! GLS closure scheme for Kz |
---|
| 310 | CASE( np_OSM ) ; CALL zdf_osm( kt, Kbb, Kmm, Krhs, avm_k, avt_k ) ! OSMOSIS closure scheme for Kz |
---|
[14834] | 311 | ! CASE( np_CST ) ! Constant Kz (reset avt, avm to the background value) |
---|
| 312 | ! ! avt_k and avm_k set one for all at initialisation phase |
---|
[8215] | 313 | !!gm avt(2:jpim1,2:jpjm1,1:jpkm1) = rn_avt0 * wmask(2:jpim1,2:jpjm1,1:jpkm1) |
---|
| 314 | !!gm avm(2:jpim1,2:jpjm1,1:jpkm1) = rn_avm0 * wmask(2:jpim1,2:jpjm1,1:jpkm1) |
---|
[14921] | 315 | END SELECT |
---|
[14072] | 316 | ! |
---|
[8215] | 317 | ! !== ocean Kz ==! (avt, avs, avm) |
---|
[15553] | 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( .NOT. l_istiled .OR. ntile == nijtile ) THEN ! Do only on the last tile |
---|
| 321 | IF( l_zdfsh2 ) CALL Agrif_avm |
---|
| 322 | ENDIF |
---|
| 323 | #endif |
---|
[8215] | 324 | ! |
---|
| 325 | ! !* start from turbulent closure values |
---|
[14834] | 326 | DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 ) |
---|
| 327 | avt(ji,jj,jk) = avt_k(ji,jj,jk) |
---|
| 328 | avm(ji,jj,jk) = avm_k(ji,jj,jk) |
---|
| 329 | END_3D |
---|
[8215] | 330 | ! |
---|
| 331 | IF( ln_rnf_mouth ) THEN !* increase diffusivity at rivers mouths |
---|
[14834] | 332 | DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, nkrnf ) |
---|
| 333 | avt(ji,jj,jk) = avt(ji,jj,jk) + 2._wp * rn_avt_rnf * rnfmsk(ji,jj) * wmask(ji,jj,jk) |
---|
| 334 | END_3D |
---|
[8215] | 335 | ENDIF |
---|
| 336 | ! |
---|
[12377] | 337 | IF( ln_zdfevd ) CALL zdf_evd( kt, Kmm, Krhs, avm, avt ) !* convection: enhanced vertical eddy diffusivity |
---|
[8215] | 338 | ! |
---|
| 339 | ! !* double diffusive mixing |
---|
| 340 | IF( ln_zdfddm ) THEN ! update avt and compute avs |
---|
[12377] | 341 | CALL zdf_ddm( kt, Kmm, avm, avt, avs ) |
---|
[8215] | 342 | ELSE ! same mixing on all tracers |
---|
[14834] | 343 | DO_3D_OVR( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 ) |
---|
| 344 | avs(ji,jj,jk) = avt(ji,jj,jk) |
---|
| 345 | END_3D |
---|
[8215] | 346 | ENDIF |
---|
| 347 | ! |
---|
[14072] | 348 | ! !* wave-induced mixing |
---|
| 349 | IF( ln_zdfswm ) CALL zdf_swm( kt, Kmm, avm, avt, avs ) ! surface wave (Qiao et al. 2004) |
---|
[12377] | 350 | IF( ln_zdfiwm ) CALL zdf_iwm( kt, Kmm, avm, avt, avs ) ! internal wave (de Lavergne et al 2017) |
---|
[8215] | 351 | |
---|
| 352 | ! !* Lateral boundary conditions (sign unchanged) |
---|
[15182] | 353 | IF(nn_hls==1) THEN |
---|
[14834] | 354 | IF( l_zdfsh2 ) THEN |
---|
| 355 | CALL lbc_lnk( 'zdfphy', avm_k, 'W', 1.0_wp , avt_k, 'W', 1.0_wp, & |
---|
| 356 | & avm , 'W', 1.0_wp , avt , 'W', 1.0_wp , avs , 'W', 1.0_wp ) |
---|
| 357 | ELSE |
---|
| 358 | CALL lbc_lnk( 'zdfphy', avm , 'W', 1.0_wp , avt , 'W', 1.0_wp , avs , 'W', 1.0_wp ) |
---|
[9107] | 359 | ENDIF |
---|
[14834] | 360 | ! |
---|
| 361 | IF( l_zdfdrg ) THEN ! drag have been updated (non-linear cases) |
---|
| 362 | IF( ln_isfcav ) THEN ; CALL lbc_lnk( 'zdfphy', rCdU_top, 'T', 1.0_wp , rCdU_bot, 'T', 1.0_wp ) ! top & bot drag |
---|
| 363 | ELSE ; CALL lbc_lnk( 'zdfphy', rCdU_bot, 'T', 1.0_wp ) ! bottom drag only |
---|
| 364 | ENDIF |
---|
| 365 | ENDIF |
---|
[8215] | 366 | ENDIF |
---|
| 367 | ! |
---|
[14834] | 368 | CALL zdf_mxl_turb( kt, Kmm ) !* turbocline depth |
---|
[8866] | 369 | ! |
---|
[14834] | 370 | IF( .NOT. l_istiled .OR. ntile == nijtile ) THEN ! Do only on the last tile |
---|
| 371 | IF( lrst_oce ) THEN !* write TKE, GLS or RIC fields in the restart file |
---|
| 372 | IF( ln_zdftke ) CALL tke_rst( kt, 'WRITE' ) |
---|
| 373 | IF( ln_zdfgls ) CALL gls_rst( kt, 'WRITE' ) |
---|
| 374 | IF( ln_zdfric ) CALL ric_rst( kt, 'WRITE' ) |
---|
| 375 | ! NB. OSMOSIS restart (osm_rst) will be called in step.F90 after ww has been updated |
---|
| 376 | ENDIF |
---|
[8215] | 377 | ENDIF |
---|
| 378 | ! |
---|
[14983] | 379 | ! diagnostics of energy dissipation |
---|
| 380 | IF( iom_use('avt_k') .OR. iom_use('avm_k') .OR. iom_use('eshear_k') .OR. iom_use('estrat_k') ) THEN |
---|
| 381 | IF( l_zdfsh2 ) THEN |
---|
| 382 | CALL iom_put( 'avt_k' , avt_k * wmask ) |
---|
| 383 | CALL iom_put( 'avm_k' , avm_k * wmask ) |
---|
| 384 | CALL iom_put( 'eshear_k', zsh2 * wmask ) |
---|
| 385 | CALL iom_put( 'estrat_k', - avt_k * rn2 * wmask ) |
---|
| 386 | ENDIF |
---|
| 387 | ENDIF |
---|
| 388 | ! |
---|
[8568] | 389 | IF( ln_timing ) CALL timing_stop('zdf_phy') |
---|
| 390 | ! |
---|
[8215] | 391 | END SUBROUTINE zdf_phy |
---|
[13558] | 392 | |
---|
| 393 | |
---|
[10364] | 394 | INTEGER FUNCTION zdf_phy_alloc() |
---|
| 395 | !!---------------------------------------------------------------------- |
---|
| 396 | !! *** FUNCTION zdf_phy_alloc *** |
---|
| 397 | !!---------------------------------------------------------------------- |
---|
| 398 | ! Allocate wi array (declared in oce.F90) for use with the adaptive-implicit vertical velocity option |
---|
| 399 | ALLOCATE( wi(jpi,jpj,jpk), Cu_adv(jpi,jpj,jpk), STAT= zdf_phy_alloc ) |
---|
| 400 | IF( zdf_phy_alloc /= 0 ) CALL ctl_warn('zdf_phy_alloc: failed to allocate ln_zad_Aimp=T required arrays') |
---|
[10425] | 401 | CALL mpp_sum ( 'zdfphy', zdf_phy_alloc ) |
---|
[10364] | 402 | END FUNCTION zdf_phy_alloc |
---|
[8215] | 403 | |
---|
| 404 | !!====================================================================== |
---|
| 405 | END MODULE zdfphy |
---|