[3] | 1 | MODULE zdfmxl |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE zdfmxl *** |
---|
| 4 | !! Ocean physics: mixed layer depth |
---|
| 5 | !!====================================================================== |
---|
[1559] | 6 | !! History : 1.0 ! 2003-08 (G. Madec) original code |
---|
[1585] | 7 | !! 3.2 ! 2009-07 (S. Masson, G. Madec) IOM + merge of DO-loop |
---|
[4990] | 8 | !! 3.7 ! 2012-03 (G. Madec) make public the density criteria for trdmxl |
---|
| 9 | !! - ! 2014-02 (F. Roquet) mixed layer depth calculated using N2 instead of rhop |
---|
[3] | 10 | !!---------------------------------------------------------------------- |
---|
[1585] | 11 | !! zdf_mxl : Compute the turbocline and mixed layer depths. |
---|
[3] | 12 | !!---------------------------------------------------------------------- |
---|
| 13 | USE oce ! ocean dynamics and tracers variables |
---|
| 14 | USE dom_oce ! ocean space and time domain variables |
---|
[1585] | 15 | USE zdf_oce ! ocean vertical physics |
---|
[3] | 16 | USE in_out_manager ! I/O manager |
---|
[258] | 17 | USE prtctl ! Print control |
---|
[4990] | 18 | USE phycst ! physical constants |
---|
[2715] | 19 | USE iom ! I/O library |
---|
| 20 | USE lib_mpp ! MPP library |
---|
[3294] | 21 | USE wrk_nemo ! work arrays |
---|
| 22 | USE timing ! Timing |
---|
[2758] | 23 | USE trc_oce, ONLY : lk_offline ! offline flag |
---|
[7591] | 24 | USE lbclnk |
---|
| 25 | USE eosbn2 ! Equation of state |
---|
[3] | 26 | |
---|
| 27 | IMPLICIT NONE |
---|
| 28 | PRIVATE |
---|
| 29 | |
---|
[2715] | 30 | PUBLIC zdf_mxl ! called by step.F90 |
---|
[6661] | 31 | PUBLIC zdf_mxl_tref ! called by asminc.F90 |
---|
[6626] | 32 | PUBLIC zdf_mxl_alloc ! Used in zdf_tke_init |
---|
[3] | 33 | |
---|
[2715] | 34 | INTEGER , PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: nmln !: number of level in the mixed layer (used by TOP) |
---|
| 35 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmld !: mixing layer depth (turbocline) [m] |
---|
| 36 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmlp !: mixed layer depth (rho=rho0+zdcrit) [m] |
---|
| 37 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmlpt !: mixed layer depth at t-points [m] |
---|
[6630] | 38 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmld_tref !: mixed layer depth at t-points - temperature criterion [m] |
---|
[7591] | 39 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmld_kara !: mixed layer depth of Kara et al [m] |
---|
[3] | 40 | |
---|
[4990] | 41 | REAL(wp), PUBLIC :: rho_c = 0.01_wp !: density criterion for mixed layer depth |
---|
| 42 | REAL(wp) :: avt_c = 5.e-4_wp ! Kz criterion for the turbocline depth |
---|
| 43 | |
---|
[7591] | 44 | ! Namelist variables for namzdf_karaml |
---|
| 45 | |
---|
| 46 | LOGICAL :: ln_kara ! Logical switch for calculating kara mixed |
---|
| 47 | ! layer |
---|
| 48 | LOGICAL :: ln_kara_write25h ! Logical switch for 25 hour outputs |
---|
| 49 | INTEGER :: jpmld_type ! mixed layer type |
---|
| 50 | REAL(wp) :: ppz_ref ! depth of initial T_ref |
---|
| 51 | REAL(wp) :: ppdT_crit ! Critical temp diff |
---|
| 52 | REAL(wp) :: ppiso_frac ! Fraction of ppdT_crit used |
---|
| 53 | |
---|
| 54 | !Used for 25h mean |
---|
| 55 | LOGICAL, PRIVATE :: kara_25h_init = .TRUE. !Logical used to initalise 25h |
---|
| 56 | !outputs. Necissary, because we need to |
---|
| 57 | !initalise the kara_25h on the zeroth |
---|
| 58 | !timestep (i.e in the nemogcm_init call) |
---|
| 59 | REAL, PRIVATE, ALLOCATABLE, DIMENSION(:,:) :: hmld_kara_25h |
---|
| 60 | |
---|
[3] | 61 | !! * Substitutions |
---|
| 62 | # include "domzgr_substitute.h90" |
---|
| 63 | !!---------------------------------------------------------------------- |
---|
[2715] | 64 | !! NEMO/OPA 4.0 , NEMO Consortium (2011) |
---|
[6613] | 65 | !! $Id$ |
---|
[2715] | 66 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
[3] | 67 | !!---------------------------------------------------------------------- |
---|
| 68 | CONTAINS |
---|
| 69 | |
---|
[2715] | 70 | INTEGER FUNCTION zdf_mxl_alloc() |
---|
| 71 | !!---------------------------------------------------------------------- |
---|
| 72 | !! *** FUNCTION zdf_mxl_alloc *** |
---|
| 73 | !!---------------------------------------------------------------------- |
---|
[2787] | 74 | zdf_mxl_alloc = 0 ! set to zero if no array to be allocated |
---|
[2758] | 75 | IF( .NOT. ALLOCATED( nmln ) ) THEN |
---|
[6630] | 76 | ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), & |
---|
| 77 | & hmld_tref(jpi,jpj), STAT= zdf_mxl_alloc ) |
---|
[2758] | 78 | ! |
---|
| 79 | IF( lk_mpp ) CALL mpp_sum ( zdf_mxl_alloc ) |
---|
| 80 | IF( zdf_mxl_alloc /= 0 ) CALL ctl_warn('zdf_mxl_alloc: failed to allocate arrays.') |
---|
| 81 | ! |
---|
| 82 | ENDIF |
---|
[2715] | 83 | END FUNCTION zdf_mxl_alloc |
---|
| 84 | |
---|
| 85 | |
---|
[3] | 86 | SUBROUTINE zdf_mxl( kt ) |
---|
| 87 | !!---------------------------------------------------------------------- |
---|
| 88 | !! *** ROUTINE zdfmxl *** |
---|
| 89 | !! |
---|
[1585] | 90 | !! ** Purpose : Compute the turbocline depth and the mixed layer depth |
---|
| 91 | !! with density criteria. |
---|
[3] | 92 | !! |
---|
[1577] | 93 | !! ** Method : The mixed layer depth is the shallowest W depth with |
---|
| 94 | !! the density of the corresponding T point (just bellow) bellow a |
---|
[4245] | 95 | !! given value defined locally as rho(10m) + rho_c |
---|
[1585] | 96 | !! The turbocline depth is the depth at which the vertical |
---|
| 97 | !! eddy diffusivity coefficient (resulting from the vertical physics |
---|
| 98 | !! alone, not the isopycnal part, see trazdf.F) fall below a given |
---|
[4990] | 99 | !! value defined locally (avt_c here taken equal to 5 cm/s2 by default) |
---|
[3] | 100 | !! |
---|
[1585] | 101 | !! ** Action : nmln, hmld, hmlp, hmlpt |
---|
[1559] | 102 | !!---------------------------------------------------------------------- |
---|
[2715] | 103 | INTEGER, INTENT(in) :: kt ! ocean time-step index |
---|
[4990] | 104 | ! |
---|
[6625] | 105 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 106 | INTEGER :: iikn, iiki, ikt, imkt ! local integer |
---|
| 107 | REAL(wp) :: zN2_c ! local scalar |
---|
[4990] | 108 | INTEGER, POINTER, DIMENSION(:,:) :: imld ! 2D workspace |
---|
[3] | 109 | !!---------------------------------------------------------------------- |
---|
[3294] | 110 | ! |
---|
| 111 | IF( nn_timing == 1 ) CALL timing_start('zdf_mxl') |
---|
| 112 | ! |
---|
| 113 | CALL wrk_alloc( jpi,jpj, imld ) |
---|
[3] | 114 | |
---|
| 115 | IF( kt == nit000 ) THEN |
---|
| 116 | IF(lwp) WRITE(numout,*) |
---|
| 117 | IF(lwp) WRITE(numout,*) 'zdf_mxl : mixed layer depth' |
---|
| 118 | IF(lwp) WRITE(numout,*) '~~~~~~~ ' |
---|
[2715] | 119 | ! ! allocate zdfmxl arrays |
---|
| 120 | IF( zdf_mxl_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'zdf_mxl : unable to allocate arrays' ) |
---|
[3] | 121 | ENDIF |
---|
| 122 | |
---|
[1559] | 123 | ! w-level of the mixing and mixed layers |
---|
[4990] | 124 | nmln(:,:) = nlb10 ! Initialization to the number of w ocean point |
---|
| 125 | hmlp(:,:) = 0._wp ! here hmlp used as a dummy variable, integrating vertically N^2 |
---|
| 126 | zN2_c = grav * rho_c * r1_rau0 ! convert density criteria into N^2 criteria |
---|
| 127 | DO jk = nlb10, jpkm1 |
---|
| 128 | DO jj = 1, jpj ! Mixed layer level: w-level |
---|
| 129 | DO ji = 1, jpi |
---|
| 130 | ikt = mbkt(ji,jj) |
---|
| 131 | hmlp(ji,jj) = hmlp(ji,jj) + MAX( rn2b(ji,jj,jk) , 0._wp ) * fse3w(ji,jj,jk) |
---|
| 132 | IF( hmlp(ji,jj) < zN2_c ) nmln(ji,jj) = MIN( jk , ikt ) + 1 ! Mixed layer level |
---|
| 133 | END DO |
---|
| 134 | END DO |
---|
| 135 | END DO |
---|
| 136 | ! |
---|
| 137 | ! w-level of the turbocline |
---|
| 138 | imld(:,:) = mbkt(:,:) + 1 ! Initialization to the number of w ocean point |
---|
| 139 | DO jk = jpkm1, nlb10, -1 ! from the bottom to nlb10 |
---|
[3] | 140 | DO jj = 1, jpj |
---|
| 141 | DO ji = 1, jpi |
---|
[6625] | 142 | imkt = mikt(ji,jj) |
---|
| 143 | IF( avt (ji,jj,jk) < avt_c ) imld(ji,jj) = MAX( imkt, jk ) ! Turbocline |
---|
[3] | 144 | END DO |
---|
| 145 | END DO |
---|
| 146 | END DO |
---|
[1559] | 147 | ! depth of the mixing and mixed layers |
---|
[7591] | 148 | |
---|
| 149 | CALL zdf_mxl_kara( kt ) ! kara MLD |
---|
| 150 | |
---|
[3] | 151 | DO jj = 1, jpj |
---|
| 152 | DO ji = 1, jpi |
---|
[1585] | 153 | iiki = imld(ji,jj) |
---|
[1577] | 154 | iikn = nmln(ji,jj) |
---|
[6625] | 155 | imkt = mikt(ji,jj) |
---|
[6626] | 156 | hmld (ji,jj) = ( fsdepw(ji,jj,iiki ) - fsdepw(ji,jj,imkt ) ) * ssmask(ji,jj) ! Turbocline depth |
---|
| 157 | hmlp (ji,jj) = ( fsdepw(ji,jj,iikn ) - fsdepw(ji,jj,imkt ) ) * ssmask(ji,jj) ! Mixed layer depth |
---|
| 158 | hmlpt(ji,jj) = ( fsdept(ji,jj,iikn-1) - fsdepw(ji,jj,imkt ) ) * ssmask(ji,jj) ! depth of the last T-point inside the mixed layer |
---|
[3] | 159 | END DO |
---|
| 160 | END DO |
---|
[6625] | 161 | IF( .NOT.lk_offline ) THEN ! no need to output in offline mode |
---|
| 162 | CALL iom_put( "mldr10_1", hmlp ) ! mixed layer depth |
---|
| 163 | CALL iom_put( "mldkz5" , hmld ) ! turbocline depth |
---|
[2758] | 164 | ENDIF |
---|
[6661] | 165 | |
---|
| 166 | IF(ln_ctl) CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ', ovlap=1 ) |
---|
| 167 | ! |
---|
| 168 | CALL wrk_dealloc( jpi,jpj, imld ) |
---|
| 169 | ! |
---|
| 170 | IF( nn_timing == 1 ) CALL timing_stop('zdf_mxl') |
---|
| 171 | ! |
---|
| 172 | END SUBROUTINE zdf_mxl |
---|
[6630] | 173 | |
---|
[6661] | 174 | |
---|
| 175 | SUBROUTINE zdf_mxl_tref() |
---|
| 176 | !!---------------------------------------------------------------------- |
---|
| 177 | !! *** ROUTINE zdf_mxl_tref *** |
---|
| 178 | !! |
---|
| 179 | !! ** Purpose : Compute the mixed layer depth with temperature criteria. |
---|
| 180 | !! |
---|
| 181 | !! ** Method : The temperature-defined mixed layer depth is required |
---|
| 182 | !! when assimilating SST in a 2D analysis. |
---|
| 183 | !! |
---|
| 184 | !! ** Action : hmld_tref |
---|
| 185 | !!---------------------------------------------------------------------- |
---|
| 186 | ! |
---|
| 187 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 188 | REAL(wp) :: t_ref ! Reference temperature |
---|
| 189 | REAL(wp) :: temp_c = 0.2 ! temperature criterion for mixed layer depth |
---|
| 190 | !!---------------------------------------------------------------------- |
---|
| 191 | ! |
---|
| 192 | ! Initialise array |
---|
| 193 | IF( zdf_mxl_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'zdf_mxl_tref : unable to allocate arrays' ) |
---|
| 194 | |
---|
[6630] | 195 | !For the AMM model assimiation uses a temperature based mixed layer depth |
---|
| 196 | !This is defined here |
---|
| 197 | DO jj = 1, jpj |
---|
| 198 | DO ji = 1, jpi |
---|
| 199 | hmld_tref(ji,jj)=fsdept(ji,jj,1 ) |
---|
| 200 | IF(ssmask(ji,jj) > 0.)THEN |
---|
| 201 | t_ref=tsn(ji,jj,1,jp_tem) |
---|
| 202 | DO jk=2,jpk |
---|
| 203 | IF(ssmask(ji,jj)==0.)THEN |
---|
| 204 | hmld_tref(ji,jj)=fsdept(ji,jj,jk ) |
---|
| 205 | EXIT |
---|
| 206 | ELSEIF( ABS(tsn(ji,jj,jk,jp_tem)-t_ref) < temp_c)THEN |
---|
| 207 | hmld_tref(ji,jj)=fsdept(ji,jj,jk ) |
---|
| 208 | ELSE |
---|
| 209 | EXIT |
---|
| 210 | ENDIF |
---|
| 211 | ENDDO |
---|
| 212 | ENDIF |
---|
| 213 | ENDDO |
---|
| 214 | ENDDO |
---|
[3] | 215 | |
---|
[6661] | 216 | END SUBROUTINE zdf_mxl_tref |
---|
| 217 | |
---|
[7591] | 218 | SUBROUTINE zdf_mxl_kara( kt ) |
---|
| 219 | !!---------------------------------------------------------------------------------- |
---|
| 220 | !! *** ROUTINE zdf_mxl_kara *** |
---|
| 221 | ! |
---|
| 222 | ! Calculate mixed layer depth according to the definition of |
---|
| 223 | ! Kara et al, 2000, JGR, 105, 16803. |
---|
| 224 | ! |
---|
| 225 | ! If mld_type=1 the mixed layer depth is calculated as the depth at which the |
---|
| 226 | ! density has increased by an amount equivalent to a temperature difference of |
---|
| 227 | ! 0.8C at the surface. |
---|
| 228 | ! |
---|
| 229 | ! For other values of mld_type the mixed layer is calculated as the depth at |
---|
| 230 | ! which the temperature differs by 0.8C from the surface temperature. |
---|
| 231 | ! |
---|
| 232 | ! Original version: David Acreman |
---|
| 233 | ! |
---|
| 234 | !!----------------------------------------------------------------------------------- |
---|
| 235 | |
---|
| 236 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
| 237 | |
---|
| 238 | NAMELIST/namzdf_karaml/ ln_kara,jpmld_type, ppz_ref, ppdT_crit, & |
---|
| 239 | & ppiso_frac, ln_kara_write25h |
---|
| 240 | |
---|
| 241 | ! Local variables |
---|
| 242 | REAL, DIMENSION(jpi,jpk) :: ppzdep ! depth for use in calculating d(rho) |
---|
| 243 | REAL(wp), DIMENSION(jpi,jpj,jpts) :: ztsn_2d !Local version of tsn |
---|
| 244 | LOGICAL :: ll_found(jpi,jpj) ! Is T_b to be found by interpolation ? |
---|
| 245 | LOGICAL :: ll_belowml(jpi,jpj,jpk) ! Flag points below mixed layer when ll_found=F |
---|
| 246 | INTEGER :: ji, jj, jk ! loop counter |
---|
| 247 | INTEGER :: ik_ref(jpi,jpj) ! index of reference level |
---|
| 248 | INTEGER :: ik_iso(jpi,jpj) ! index of last uniform temp level |
---|
| 249 | REAL :: zT(jpi,jpj,jpk) ! Temperature or denisty |
---|
| 250 | REAL :: zT_ref(jpi,jpj) ! reference temperature |
---|
| 251 | REAL :: zT_b ! base temperature |
---|
| 252 | REAL :: zdTdz(jpi,jpj,jpk-2) ! gradient of zT |
---|
| 253 | REAL :: zmoddT(jpi,jpj,jpk-2) ! Absolute temperature difference |
---|
| 254 | REAL :: zdz ! depth difference |
---|
| 255 | REAL :: zdT ! temperature difference |
---|
| 256 | REAL :: zdelta_T(jpi,jpj) ! difference critereon |
---|
| 257 | REAL :: zRHO1(jpi,jpj), zRHO2(jpi,jpj) ! Densities |
---|
| 258 | INTEGER, SAVE :: idt, i_steps |
---|
| 259 | INTEGER, SAVE :: i_cnt_25h ! Counter for 25 hour means |
---|
| 260 | INTEGER :: ios ! Local integer output status for namelist read |
---|
| 261 | |
---|
| 262 | |
---|
| 263 | |
---|
| 264 | !!------------------------------------------------------------------------------------- |
---|
| 265 | |
---|
| 266 | IF( kt == nit000 ) THEN |
---|
| 267 | ! Default values |
---|
| 268 | ln_kara = .FALSE. |
---|
| 269 | ln_kara_write25h = .FALSE. |
---|
| 270 | jpmld_type = 1 |
---|
| 271 | ppz_ref = 10.0 |
---|
| 272 | ppdT_crit = 0.2 |
---|
| 273 | ppiso_frac = 0.1 |
---|
| 274 | ! Read namelist |
---|
| 275 | REWIND ( numnam_ref ) |
---|
| 276 | READ ( numnam_ref, namzdf_karaml, IOSTAT=ios, ERR= 901 ) |
---|
| 277 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_karaml in reference namelist', lwp ) |
---|
| 278 | |
---|
| 279 | REWIND( numnam_cfg ) ! Namelist nam_diadiaop in configuration namelist 3D hourly diagnostics |
---|
| 280 | READ ( numnam_cfg, namzdf_karaml, IOSTAT = ios, ERR = 902 ) |
---|
| 281 | 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_karaml in configuration namelist', lwp ) |
---|
| 282 | IF(lwm) WRITE ( numond, namzdf_karaml ) |
---|
| 283 | |
---|
| 284 | |
---|
| 285 | WRITE(numout,*) '===== Kara mixed layer =====' |
---|
| 286 | WRITE(numout,*) 'ln_kara = ', ln_kara |
---|
| 287 | WRITE(numout,*) 'jpmld_type = ', jpmld_type |
---|
| 288 | WRITE(numout,*) 'ppz_ref = ', ppz_ref |
---|
| 289 | WRITE(numout,*) 'ppdT_crit = ', ppdT_crit |
---|
| 290 | WRITE(numout,*) 'ppiso_frac = ', ppiso_frac |
---|
| 291 | WRITE(numout,*) 'ln_kara_write25h = ', ln_kara_write25h |
---|
| 292 | WRITE(numout,*) '============================' |
---|
| 293 | |
---|
| 294 | IF ( .NOT.ln_kara ) THEN |
---|
| 295 | WRITE(numout,*) "ln_kara not set; Kara mixed layer not calculated" |
---|
| 296 | ELSE |
---|
| 297 | IF (.NOT. ALLOCATED(hmld_kara) ) ALLOCATE( hmld_kara(jpi,jpj) ) |
---|
| 298 | IF ( ln_kara_write25h .AND. kara_25h_init ) THEN |
---|
| 299 | i_cnt_25h = 0 |
---|
| 300 | IF (.NOT. ALLOCATED(hmld_kara_25h) ) & |
---|
| 301 | & ALLOCATE( hmld_kara_25h(jpi,jpj) ) |
---|
| 302 | hmld_kara_25h = 0._wp |
---|
| 303 | IF( nacc == 1 ) THEN |
---|
| 304 | idt = INT(rdtmin) |
---|
| 305 | ELSE |
---|
| 306 | idt = INT(rdt) |
---|
| 307 | ENDIF |
---|
| 308 | IF( MOD( 3600,idt) == 0 ) THEN |
---|
| 309 | i_steps = 3600 / idt |
---|
| 310 | ELSE |
---|
| 311 | CALL ctl_stop('STOP', & |
---|
| 312 | & 'zdf_mxl_kara: timestep must give MOD(3600,rdt)'// & |
---|
| 313 | & ' = 0 otherwise no hourly values are possible') |
---|
| 314 | ENDIF |
---|
| 315 | ENDIF |
---|
| 316 | ENDIF |
---|
| 317 | ENDIF |
---|
| 318 | |
---|
| 319 | IF ( ln_kara ) THEN |
---|
| 320 | |
---|
| 321 | !set critical ln_kara |
---|
| 322 | ztsn_2d = tsn(:,:,1,:) |
---|
| 323 | ztsn_2d(:,:,jp_tem) = ztsn_2d(:,:,jp_tem) + ppdT_crit |
---|
| 324 | |
---|
| 325 | ! Set the mixed layer depth criterion at each grid point |
---|
| 326 | ppzdep = 0._wp |
---|
| 327 | IF( jpmld_type == 1 ) THEN |
---|
| 328 | CALL eos ( tsn(:,:,1,:), & |
---|
| 329 | & ppzdep(:,:), zRHO1(:,:) ) |
---|
| 330 | CALL eos ( ztsn_2d(:,:,:), & |
---|
| 331 | & ppzdep(:,:), zRHO2(:,:) ) |
---|
| 332 | zdelta_T(:,:) = abs( zRHO1(:,:) - zRHO2(:,:) ) * rau0 |
---|
| 333 | ! RHO from eos (2d version) doesn't calculate north or east halo: |
---|
| 334 | CALL lbc_lnk( zdelta_T, 'T', 1. ) |
---|
| 335 | zT(:,:,:) = rhop(:,:,:) |
---|
| 336 | ELSE |
---|
| 337 | zdelta_T(:,:) = ppdT_crit |
---|
| 338 | zT(:,:,:) = tsn(:,:,:,jp_tem) |
---|
| 339 | ENDIF |
---|
| 340 | |
---|
| 341 | ! Calculate the gradient of zT and absolute difference for use later |
---|
| 342 | DO jk = 1 ,jpk-2 |
---|
| 343 | zdTdz(:,:,jk) = ( zT(:,:,jk+1) - zT(:,:,jk) ) / fse3w(:,:,jk+1) |
---|
| 344 | zmoddT(:,:,jk) = abs( zT(:,:,jk+1) - zT(:,:,jk) ) |
---|
| 345 | END DO |
---|
| 346 | |
---|
| 347 | ! Find density/temperature at the reference level (Kara et al use 10m). |
---|
| 348 | ! ik_ref is the index of the box centre immediately above or at the reference level |
---|
| 349 | ! Find ppz_ref in the array of model level depths and find the ref |
---|
| 350 | ! density/temperature by linear interpolation. |
---|
| 351 | ik_ref = -1 |
---|
| 352 | DO jk = jpkm1, 2, -1 |
---|
| 353 | WHERE( fsdept(:,:,jk) > ppz_ref ) |
---|
| 354 | ik_ref(:,:) = jk - 1 |
---|
| 355 | zT_ref(:,:) = zT(:,:,jk-1) + & |
---|
| 356 | & zdTdz(:,:,jk-1) * ( ppz_ref - fsdept(:,:,jk-1) ) |
---|
| 357 | ENDWHERE |
---|
| 358 | END DO |
---|
| 359 | IF ( ANY( ik_ref < 0 ) .OR. ANY( ik_ref > jpkm1 ) ) THEN |
---|
| 360 | CALL ctl_stop( "STOP", & |
---|
| 361 | & "zdf_mxl_kara: unable to find reference level for kara ML" ) |
---|
| 362 | ELSE |
---|
| 363 | ! If the first grid box centre is below the reference level then use the |
---|
| 364 | ! top model level to get zT_ref |
---|
| 365 | WHERE( fsdept(:,:,1) > ppz_ref ) |
---|
| 366 | zT_ref = zT(:,:,1) |
---|
| 367 | ik_ref = 1 |
---|
| 368 | ENDWHERE |
---|
| 369 | |
---|
| 370 | ! Search for a uniform density/temperature region where adjacent levels |
---|
| 371 | ! differ by less than ppiso_frac * deltaT. |
---|
| 372 | ! ik_iso is the index of the last level in the uniform layer |
---|
| 373 | ! ll_found indicates whether the mixed layer depth can be found by interpolation |
---|
| 374 | ik_iso(:,:) = ik_ref(:,:) |
---|
| 375 | ll_found(:,:) = .false. |
---|
| 376 | DO jj = 1, nlcj |
---|
| 377 | DO ji = 1, nlci |
---|
| 378 | !CDIR NOVECTOR |
---|
| 379 | DO jk = ik_ref(ji,jj), mbathy(ji,jj)-1 |
---|
| 380 | IF( zmoddT(ji,jj,jk) > ( ppiso_frac * zdelta_T(ji,jj) ) ) THEN |
---|
| 381 | ik_iso(ji,jj) = jk |
---|
| 382 | ll_found(ji,jj) = ( zmoddT(ji,jj,jk) > zdelta_T(ji,jj) ) |
---|
| 383 | EXIT |
---|
| 384 | ENDIF |
---|
| 385 | END DO |
---|
| 386 | END DO |
---|
| 387 | END DO |
---|
| 388 | |
---|
| 389 | ! Use linear interpolation to find depth of mixed layer base where possible |
---|
| 390 | hmld_kara(:,:) = ppz_ref |
---|
| 391 | DO jj = 1, jpj |
---|
| 392 | DO ji = 1, jpi |
---|
| 393 | IF( ll_found(ji,jj) .and. tmask(ji,jj,1) == 1.0 ) THEN |
---|
| 394 | zdz = abs( zdelta_T(ji,jj) / zdTdz(ji,jj,ik_iso(ji,jj)) ) |
---|
| 395 | hmld_kara(ji,jj) = fsdept(ji,jj,ik_iso(ji,jj)) + zdz |
---|
| 396 | ENDIF |
---|
| 397 | END DO |
---|
| 398 | END DO |
---|
| 399 | |
---|
| 400 | ! If ll_found = .false. then calculate MLD using difference of zdelta_T |
---|
| 401 | ! from the reference density/temperature |
---|
| 402 | |
---|
| 403 | ! Prevent this section from working on land points |
---|
| 404 | WHERE( tmask(:,:,1) /= 1.0 ) |
---|
| 405 | ll_found = .true. |
---|
| 406 | ENDWHERE |
---|
| 407 | |
---|
| 408 | DO jk = 1, jpk |
---|
| 409 | ll_belowml(:,:,jk) = abs( zT(:,:,jk) - zT_ref(:,:) ) >= & |
---|
| 410 | & zdelta_T(:,:) |
---|
| 411 | END DO |
---|
| 412 | |
---|
| 413 | ! Set default value where interpolation cannot be used (ll_found=false) |
---|
| 414 | DO jj = 1, jpj |
---|
| 415 | DO ji = 1, jpi |
---|
| 416 | IF( .NOT. ll_found(ji,jj) ) & |
---|
| 417 | & hmld_kara(ji,jj) = fsdept(ji,jj,mbathy(ji,jj)) |
---|
| 418 | END DO |
---|
| 419 | END DO |
---|
| 420 | |
---|
| 421 | DO jj = 1, jpj |
---|
| 422 | DO ji = 1, jpi |
---|
| 423 | !CDIR NOVECTOR |
---|
| 424 | DO jk = ik_ref(ji,jj)+1, mbathy(ji,jj) |
---|
| 425 | IF( ll_found(ji,jj) ) EXIT |
---|
| 426 | IF( ll_belowml(ji,jj,jk) ) THEN |
---|
| 427 | zT_b = zT_ref(ji,jj) + zdelta_T(ji,jj) * & |
---|
| 428 | & SIGN(1.0, zdTdz(ji,jj,jk-1) ) |
---|
| 429 | zdT = zT_b - zT(ji,jj,jk-1) |
---|
| 430 | zdz = zdT / zdTdz(ji,jj,jk-1) |
---|
| 431 | hmld_kara(ji,jj) = fsdept(ji,jj,jk-1) + zdz |
---|
| 432 | EXIT |
---|
| 433 | ENDIF |
---|
| 434 | END DO |
---|
| 435 | END DO |
---|
| 436 | END DO |
---|
| 437 | |
---|
| 438 | hmld_kara(:,:) = hmld_kara(:,:) * tmask(:,:,1) |
---|
| 439 | |
---|
| 440 | IF( ln_kara_write25h ) THEN |
---|
| 441 | !Double IF required as i_steps not defined if ln_kara_write25h = |
---|
| 442 | ! FALSE |
---|
| 443 | IF ( ( MOD( kt, i_steps ) == 0 ) .OR. kara_25h_init ) THEN |
---|
| 444 | hmld_kara_25h = hmld_kara_25h + hmld_kara |
---|
| 445 | i_cnt_25h = i_cnt_25h + 1 |
---|
| 446 | IF ( kara_25h_init ) kara_25h_init = .FALSE. |
---|
| 447 | ENDIF |
---|
| 448 | ENDIF |
---|
| 449 | |
---|
| 450 | #if defined key_iomput |
---|
| 451 | IF( kt /= nit000 ) THEN |
---|
| 452 | CALL iom_put( "mldkara" , hmld_kara ) |
---|
| 453 | IF( ( MOD( i_cnt_25h, 25) == 0 ) .AND. ln_kara_write25h ) & |
---|
| 454 | CALL iom_put( "kara25h" , ( hmld_kara_25h / 25._wp ) ) |
---|
| 455 | ENDIF |
---|
| 456 | #endif |
---|
| 457 | |
---|
| 458 | ENDIF |
---|
| 459 | ENDIF |
---|
| 460 | |
---|
| 461 | END SUBROUTINE zdf_mxl_kara |
---|
| 462 | |
---|
[3] | 463 | !!====================================================================== |
---|
| 464 | END MODULE zdfmxl |
---|