Changeset 7976
- Timestamp:
- 2017-04-26T17:50:32+02:00 (7 years ago)
- Location:
- branches/UKMO/AMM15_v3_6_STABLE/NEMOGCM
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/UKMO/AMM15_v3_6_STABLE/NEMOGCM/CONFIG/SHARED/namelist_ref
r5578 r7976 904 904 !! namzdf_ddm double diffusive mixing parameterization ("key_zdfddm") 905 905 !! namzdf_tmx tidal mixing parameterization ("key_zdftmx") 906 !! namzdf_mldzint vertically-interpolated mixed-layer depth parameters 906 907 !!====================================================================== 907 908 ! … … 1009 1010 ln_tmx_itf = .true. ! ITF specific parameterisation 1010 1011 rn_tfe_itf = 1. ! ITF tidal dissipation efficiency 1012 / 1013 !------------------------------------------------------------------------------------------ 1014 &namzdf_mldzint ! Parameters for vertically-interpolated mixed-layer depth diagnostic 1015 !------------------------------------------------------------------------------------------ 1016 nn_mld_diag = 0 ! Number of MLD diagnostics to use from below 1017 1018 ! ! MLD criterion ! Reference ! Finite difference ! Gradient layer ! 1019 ! ! type ! depth ! criterion ! criterion ! 1020 sn_mld1 = 1 , 10.0 , 0.2 , 0.1 1021 sn_mld2 = 0 , 0.0 , 0.0 , 0.0 1022 sn_mld3 = 0 , 0.0 , 0.0 , 0.0 1023 sn_mld4 = 0 , 0.0 , 0.0 , 0.0 1024 sn_mld5 = 0 , 0.0 , 0.0 , 0.0 1011 1025 / 1012 1026 -
branches/UKMO/AMM15_v3_6_STABLE/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90
r6928 r7976 47 47 USE ioipsl 48 48 USE dynspg_oce, ONLY: un_adv, vn_adv ! barotropic velocities 49 USE eosbn2 ! equation of state (eos_bn2 routine) 50 49 51 50 52 #if defined key_lim2 … … 132 134 REAL(wp), POINTER, DIMENSION(:,:) :: z2d ! 2D workspace 133 135 REAL(wp), POINTER, DIMENSION(:,:,:) :: z3d ! 3D workspace 136 REAL(wp), POINTER, DIMENSION(:,:,:) :: zrhd , zrhop ! 3D workspace 134 137 !!---------------------------------------------------------------------- 135 138 ! … … 138 141 CALL wrk_alloc( jpi , jpj , z2d ) 139 142 CALL wrk_alloc( jpi , jpj, jpk , z3d ) 143 CALL wrk_alloc( jpi , jpj, jpk , zrhd , zrhop ) 140 144 ! 141 145 ! Output the initial state and forcings … … 376 380 CALL iom_put( "v_salttr", 0.5 * z2d ) ! heat transport in j-direction 377 381 ENDIF 382 383 IF( iom_use("rhop") ) THEN 384 CALL eos( tsn, zrhd, zrhop, fsdept_n(:,:,:) ) ! now in situ and potential density 385 zrhop(:,:,jpk) = 0._wp 386 CALL iom_put( 'rhop', zrhop ) 387 ENDIF 388 378 389 ! 379 390 CALL wrk_dealloc( jpi , jpj , z2d ) 380 391 CALL wrk_dealloc( jpi , jpj, jpk , z3d ) 392 CALL wrk_dealloc( jpi , jpj, jpk , zrhd , zrhop ) 381 393 ! 382 394 IF( nn_timing == 1 ) CALL timing_stop('dia_wri') -
branches/UKMO/AMM15_v3_6_STABLE/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfmxl.F90
r6928 r7976 18 18 USE phycst ! physical constants 19 19 USE iom ! I/O library 20 USE eosbn2 ! for zdf_mxl_zint 20 21 USE lib_mpp ! MPP library 21 22 USE wrk_nemo ! work arrays … … 27 28 28 29 PUBLIC zdf_mxl ! called by step.F90 29 PUBLIC zdf_mxl_alloc ! Used in zdf_tke_init30 30 31 31 INTEGER , PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: nmln !: number of level in the mixed layer (used by TOP) … … 33 33 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmlp !: mixed layer depth (rho=rho0+zdcrit) [m] 34 34 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmlpt !: mixed layer depth at t-points [m] 35 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: hmld_zint !: vertically-interpolated mixed layer depth [m] 36 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: htc_mld ! Heat content of hmld_zint 37 LOGICAL, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ll_found ! Is T_b to be found by interpolation ? 38 LOGICAL, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: ll_belowml ! Flag points below mixed layer when ll_found=F 35 39 36 40 REAL(wp), PUBLIC :: rho_c = 0.01_wp !: density criterion for mixed layer depth 37 41 REAL(wp) :: avt_c = 5.e-4_wp ! Kz criterion for the turbocline depth 42 43 TYPE, PUBLIC :: MXL_ZINT !: Structure for MLD defs 44 INTEGER :: mld_type ! mixed layer type 45 REAL(wp) :: zref ! depth of initial T_ref 46 REAL(wp) :: dT_crit ! Critical temp diff 47 REAL(wp) :: iso_frac ! Fraction of rn_dT_crit used 48 END TYPE MXL_ZINT 38 49 39 50 !! * Substitutions … … 52 63 zdf_mxl_alloc = 0 ! set to zero if no array to be allocated 53 64 IF( .NOT. ALLOCATED( nmln ) ) THEN 54 ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), STAT= zdf_mxl_alloc ) 65 ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), hmld_zint(jpi,jpj), & 66 & htc_mld(jpi,jpj), & 67 & ll_found(jpi,jpj), ll_belowml(jpi,jpj,jpk), STAT= zdf_mxl_alloc ) 55 68 ! 56 69 IF( lk_mpp ) CALL mpp_sum ( zdf_mxl_alloc ) … … 128 141 iikn = nmln(ji,jj) 129 142 imkt = mikt(ji,jj) 130 hmld (ji,jj) = ( fsdepw(ji,jj,iiki ) - fsdepw(ji,jj,imkt ) )* ssmask(ji,jj) ! Turbocline depth131 hmlp (ji,jj) = ( fsdepw(ji,jj,iikn ) - fsdepw(ji,jj, imkt) ) * ssmask(ji,jj) ! Mixed layer depth132 hmlpt(ji,jj) = ( fsdept(ji,jj,iikn-1) - fsdepw(ji,jj,imkt ) )* ssmask(ji,jj) ! depth of the last T-point inside the mixed layer143 hmld (ji,jj) = ( fsdepw(ji,jj,iiki ) - fsdepw(ji,jj,imkt ) ) * ssmask(ji,jj) ! Turbocline depth 144 hmlp (ji,jj) = ( fsdepw(ji,jj,iikn ) - fsdepw(ji,jj,MAX( imkt,nla10 ) ) ) * ssmask(ji,jj) ! Mixed layer depth 145 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 133 146 END DO 134 147 END DO … … 138 151 ENDIF 139 152 153 ! Vertically-interpolated mixed-layer depth diagnostic 154 CALL zdf_mxl_zint( kt ) 155 140 156 IF(ln_ctl) CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ', ovlap=1 ) 141 157 ! … … 145 161 ! 146 162 END SUBROUTINE zdf_mxl 163 164 SUBROUTINE zdf_mxl_zint_mld( sf ) 165 !!---------------------------------------------------------------------------------- 166 !! *** ROUTINE zdf_mxl_zint_mld *** 167 ! 168 ! Calculate vertically-interpolated mixed layer depth diagnostic. 169 ! 170 ! This routine can calculate the mixed layer depth diagnostic suggested by 171 ! Kara et al, 2000, JGR, 105, 16803, but is more general and can calculate 172 ! vertically-interpolated mixed-layer depth diagnostics with other parameter 173 ! settings set in the namzdf_mldzint namelist. 174 ! 175 ! If mld_type=1 the mixed layer depth is calculated as the depth at which the 176 ! density has increased by an amount equivalent to a temperature difference of 177 ! 0.8C at the surface. 178 ! 179 ! For other values of mld_type the mixed layer is calculated as the depth at 180 ! which the temperature differs by 0.8C from the surface temperature. 181 ! 182 ! David Acreman, Daley Calvert 183 ! 184 !!----------------------------------------------------------------------------------- 185 186 TYPE(MXL_ZINT), INTENT(in) :: sf 187 188 ! Diagnostic criteria 189 INTEGER :: nn_mld_type ! mixed layer type 190 REAL(wp) :: rn_zref ! depth of initial T_ref 191 REAL(wp) :: rn_dT_crit ! Critical temp diff 192 REAL(wp) :: rn_iso_frac ! Fraction of rn_dT_crit used 193 194 ! Local variables 195 REAL(wp), PARAMETER :: zepsilon = 1.e-30 ! local small value 196 INTEGER, POINTER, DIMENSION(:,:) :: ikmt ! number of active tracer levels 197 INTEGER, POINTER, DIMENSION(:,:) :: ik_ref ! index of reference level 198 INTEGER, POINTER, DIMENSION(:,:) :: ik_iso ! index of last uniform temp level 199 REAL, POINTER, DIMENSION(:,:,:) :: zT ! Temperature or density 200 REAL, POINTER, DIMENSION(:,:) :: ppzdep ! depth for use in calculating d(rho) 201 REAL, POINTER, DIMENSION(:,:) :: zT_ref ! reference temperature 202 REAL :: zT_b ! base temperature 203 REAL, POINTER, DIMENSION(:,:,:) :: zdTdz ! gradient of zT 204 REAL, POINTER, DIMENSION(:,:,:) :: zmoddT ! Absolute temperature difference 205 REAL :: zdz ! depth difference 206 REAL :: zdT ! temperature difference 207 REAL, POINTER, DIMENSION(:,:) :: zdelta_T ! difference critereon 208 REAL, POINTER, DIMENSION(:,:) :: zRHO1, zRHO2 ! Densities 209 INTEGER :: ji, jj, jk ! loop counter 210 211 !!------------------------------------------------------------------------------------- 212 ! 213 CALL wrk_alloc( jpi, jpj, ikmt, ik_ref, ik_iso) 214 CALL wrk_alloc( jpi, jpj, ppzdep, zT_ref, zdelta_T, zRHO1, zRHO2 ) 215 CALL wrk_alloc( jpi, jpj, jpk, zT, zdTdz, zmoddT ) 216 217 ! Unpack structure 218 nn_mld_type = sf%mld_type 219 rn_zref = sf%zref 220 rn_dT_crit = sf%dT_crit 221 rn_iso_frac = sf%iso_frac 222 223 ! Set the mixed layer depth criterion at each grid point 224 IF( nn_mld_type == 0 ) THEN 225 zdelta_T(:,:) = rn_dT_crit 226 zT(:,:,:) = rhop(:,:,:) 227 ELSE IF( nn_mld_type == 1 ) THEN 228 ppzdep(:,:)=0.0 229 call eos ( tsn(:,:,1,:), ppzdep(:,:), zRHO1(:,:) ) 230 ! Use zT temporarily as a copy of tsn with rn_dT_crit added to SST 231 ! [assumes number of tracers less than number of vertical levels] 232 zT(:,:,1:jpts)=tsn(:,:,1,1:jpts) 233 zT(:,:,jp_tem)=zT(:,:,1)+rn_dT_crit 234 CALL eos( zT(:,:,1:jpts), ppzdep(:,:), zRHO2(:,:) ) 235 zdelta_T(:,:) = abs( zRHO1(:,:) - zRHO2(:,:) ) * rau0 236 ! RHO from eos (2d version) doesn't calculate north or east halo: 237 CALL lbc_lnk( zdelta_T, 'T', 1. ) 238 zT(:,:,:) = rhop(:,:,:) 239 ELSE 240 zdelta_T(:,:) = rn_dT_crit 241 zT(:,:,:) = tsn(:,:,:,jp_tem) 242 END IF 243 244 ! Calculate the gradient of zT and absolute difference for use later 245 DO jk = 1 ,jpk-2 246 zdTdz(:,:,jk) = ( zT(:,:,jk+1) - zT(:,:,jk) ) / fse3w(:,:,jk+1) 247 zmoddT(:,:,jk) = abs( zT(:,:,jk+1) - zT(:,:,jk) ) 248 END DO 249 250 ! Find density/temperature at the reference level (Kara et al use 10m). 251 ! ik_ref is the index of the box centre immediately above or at the reference level 252 ! Find rn_zref in the array of model level depths and find the ref 253 ! density/temperature by linear interpolation. 254 DO jk = jpkm1, 2, -1 255 WHERE ( fsdept(:,:,jk) > rn_zref ) 256 ik_ref(:,:) = jk - 1 257 zT_ref(:,:) = zT(:,:,jk-1) + zdTdz(:,:,jk-1) * ( rn_zref - fsdept(:,:,jk-1) ) 258 END WHERE 259 END DO 260 261 ! If the first grid box centre is below the reference level then use the 262 ! top model level to get zT_ref 263 WHERE ( fsdept(:,:,1) > rn_zref ) 264 zT_ref = zT(:,:,1) 265 ik_ref = 1 266 END WHERE 267 268 ! The number of active tracer levels is 1 less than the number of active w levels 269 ikmt(:,:) = mbathy(:,:) - 1 270 271 ! Initialize / reset 272 ll_found(:,:) = .false. 273 274 IF ( rn_iso_frac - zepsilon > 0. ) THEN 275 ! Search for a uniform density/temperature region where adjacent levels 276 ! differ by less than rn_iso_frac * deltaT. 277 ! ik_iso is the index of the last level in the uniform layer 278 ! ll_found indicates whether the mixed layer depth can be found by interpolation 279 ik_iso(:,:) = ik_ref(:,:) 280 DO jj = 1, nlcj 281 DO ji = 1, nlci 282 !CDIR NOVECTOR 283 DO jk = ik_ref(ji,jj), ikmt(ji,jj)-1 284 IF ( zmoddT(ji,jj,jk) > ( rn_iso_frac * zdelta_T(ji,jj) ) ) THEN 285 ik_iso(ji,jj) = jk 286 ll_found(ji,jj) = ( zmoddT(ji,jj,jk) > zdelta_T(ji,jj) ) 287 EXIT 288 END IF 289 END DO 290 END DO 291 END DO 292 293 ! Use linear interpolation to find depth of mixed layer base where possible 294 hmld_zint(:,:) = rn_zref 295 DO jj = 1, jpj 296 DO ji = 1, jpi 297 IF (ll_found(ji,jj) .and. tmask(ji,jj,1) == 1.0) THEN 298 zdz = abs( zdelta_T(ji,jj) / zdTdz(ji,jj,ik_iso(ji,jj)) ) 299 hmld_zint(ji,jj) = fsdept(ji,jj,ik_iso(ji,jj)) + zdz 300 END IF 301 END DO 302 END DO 303 END IF 304 305 ! If ll_found = .false. then calculate MLD using difference of zdelta_T 306 ! from the reference density/temperature 307 308 ! Prevent this section from working on land points 309 WHERE ( tmask(:,:,1) /= 1.0 ) 310 ll_found = .true. 311 END WHERE 312 313 DO jk=1, jpk 314 ll_belowml(:,:,jk) = abs( zT(:,:,jk) - zT_ref(:,:) ) >= zdelta_T(:,:) 315 END DO 316 317 ! Set default value where interpolation cannot be used (ll_found=false) 318 DO jj = 1, jpj 319 DO ji = 1, jpi 320 IF ( .not. ll_found(ji,jj) ) hmld_zint(ji,jj) = fsdept(ji,jj,ikmt(ji,jj)) 321 END DO 322 END DO 323 324 DO jj = 1, jpj 325 DO ji = 1, jpi 326 !CDIR NOVECTOR 327 DO jk = ik_ref(ji,jj)+1, ikmt(ji,jj) 328 IF ( ll_found(ji,jj) ) EXIT 329 IF ( ll_belowml(ji,jj,jk) ) THEN 330 zT_b = zT_ref(ji,jj) + zdelta_T(ji,jj) * SIGN(1.0, zdTdz(ji,jj,jk-1) ) 331 zdT = zT_b - zT(ji,jj,jk-1) 332 zdz = zdT / zdTdz(ji,jj,jk-1) 333 hmld_zint(ji,jj) = fsdept(ji,jj,jk-1) + zdz 334 EXIT 335 END IF 336 END DO 337 END DO 338 END DO 339 340 hmld_zint(:,:) = hmld_zint(:,:)*tmask(:,:,1) 341 ! 342 CALL wrk_dealloc( jpi, jpj, ikmt, ik_ref, ik_iso) 343 CALL wrk_dealloc( jpi, jpj, ppzdep, zT_ref, zdelta_T, zRHO1, zRHO2 ) 344 CALL wrk_dealloc( jpi,jpj, jpk, zT, zdTdz, zmoddT ) 345 ! 346 END SUBROUTINE zdf_mxl_zint_mld 347 348 SUBROUTINE zdf_mxl_zint_htc( kt ) 349 !!---------------------------------------------------------------------- 350 !! *** ROUTINE zdf_mxl_zint_htc *** 351 !! 352 !! ** Purpose : 353 !! 354 !! ** Method : 355 !!---------------------------------------------------------------------- 356 357 INTEGER, INTENT(in) :: kt ! ocean time-step index 358 359 INTEGER :: ji, jj, jk 360 INTEGER :: ikmax 361 REAL(wp) :: zc, zcoef 362 ! 363 INTEGER, ALLOCATABLE, DIMENSION(:,:) :: ilevel 364 REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zthick_0, zthick 365 366 !!---------------------------------------------------------------------- 367 368 IF( .NOT. ALLOCATED(ilevel) ) THEN 369 ALLOCATE( ilevel(jpi,jpj), zthick_0(jpi,jpj), & 370 & zthick(jpi,jpj), STAT=ji ) 371 IF( lk_mpp ) CALL mpp_sum(ji) 372 IF( ji /= 0 ) CALL ctl_stop( 'STOP', 'zdf_mxl_zint_htc : unable to allocate arrays' ) 373 ENDIF 374 375 ! Find last whole model T level above the MLD 376 ilevel(:,:) = 0 377 zthick_0(:,:) = 0._wp 378 379 DO jk = 1, jpkm1 380 DO jj = 1, jpj 381 DO ji = 1, jpi 382 zthick_0(ji,jj) = zthick_0(ji,jj) + fse3t(ji,jj,jk) 383 IF( zthick_0(ji,jj) < hmld_zint(ji,jj) ) ilevel(ji,jj) = jk 384 END DO 385 END DO 386 WRITE(numout,*) 'zthick_0(jk =',jk,') =',zthick_0(2,2) 387 WRITE(numout,*) 'fsdepw(jk+1 =',jk+1,') =',fsdepw(2,2,jk+1) 388 END DO 389 390 ! Surface boundary condition 391 IF( lk_vvl ) THEN ; zthick(:,:) = 0._wp ; htc_mld(:,:) = 0._wp 392 ELSE ; zthick(:,:) = sshn(:,:) ; htc_mld(:,:) = tsn(:,:,1,jp_tem) * sshn(:,:) * tmask(:,:,1) 393 ENDIF 394 395 ! Deepest whole T level above the MLD 396 ikmax = MIN( MAXVAL( ilevel(:,:) ), jpkm1 ) 397 398 ! Integration down to last whole model T level 399 DO jk = 1, ikmax 400 DO jj = 1, jpj 401 DO ji = 1, jpi 402 zc = fse3t(ji,jj,jk) * REAL( MIN( MAX( 0, ilevel(ji,jj) - jk + 1 ) , 1 ) ) ! 0 below ilevel 403 zthick(ji,jj) = zthick(ji,jj) + zc 404 htc_mld(ji,jj) = htc_mld(ji,jj) + zc * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) 405 END DO 406 END DO 407 END DO 408 409 ! Subsequent partial T level 410 zthick(:,:) = hmld_zint(:,:) - zthick(:,:) ! remaining thickness to reach MLD 411 412 DO jj = 1, jpj 413 DO ji = 1, jpi 414 htc_mld(ji,jj) = htc_mld(ji,jj) + tsn(ji,jj,ilevel(ji,jj)+1,jp_tem) & 415 & * MIN( fse3t(ji,jj,ilevel(ji,jj)+1), zthick(ji,jj) ) * tmask(ji,jj,ilevel(ji,jj)+1) 416 END DO 417 END DO 418 419 WRITE(numout,*) 'htc_mld(after) =',htc_mld(2,2) 420 421 ! Convert to heat content 422 zcoef = rau0 * rcp 423 htc_mld(:,:) = zcoef * htc_mld(:,:) 424 425 END SUBROUTINE zdf_mxl_zint_htc 426 427 SUBROUTINE zdf_mxl_zint( kt ) 428 !!---------------------------------------------------------------------- 429 !! *** ROUTINE zdf_mxl_zint *** 430 !! 431 !! ** Purpose : 432 !! 433 !! ** Method : 434 !!---------------------------------------------------------------------- 435 436 INTEGER, INTENT(in) :: kt ! ocean time-step index 437 438 INTEGER :: ios 439 INTEGER :: jn 440 441 INTEGER :: nn_mld_diag = 0 ! number of diagnostics 442 443 CHARACTER(len=1) :: cmld 444 445 TYPE(MXL_ZINT) :: sn_mld1, sn_mld2, sn_mld3, sn_mld4, sn_mld5 446 TYPE(MXL_ZINT), SAVE, DIMENSION(5) :: mld_diags 447 448 NAMELIST/namzdf_mldzint/ nn_mld_diag, sn_mld1, sn_mld2, sn_mld3, sn_mld4, sn_mld5 449 450 !!---------------------------------------------------------------------- 451 452 IF( kt == nit000 ) THEN 453 REWIND( numnam_ref ) ! Namelist namzdf_mldzint in reference namelist 454 READ ( numnam_ref, namzdf_mldzint, IOSTAT = ios, ERR = 901) 455 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_mldzint in reference namelist', lwp ) 456 457 REWIND( numnam_cfg ) ! Namelist namzdf_mldzint in configuration namelist 458 READ ( numnam_cfg, namzdf_mldzint, IOSTAT = ios, ERR = 902 ) 459 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_mldzint in configuration namelist', lwp ) 460 IF(lwm) WRITE ( numond, namzdf_mldzint ) 461 462 IF( nn_mld_diag > 5 ) CALL ctl_stop( 'STOP', 'zdf_mxl_ini: Specify no more than 5 MLD definitions' ) 463 464 mld_diags(1) = sn_mld1 465 mld_diags(2) = sn_mld2 466 mld_diags(3) = sn_mld3 467 mld_diags(4) = sn_mld4 468 mld_diags(5) = sn_mld5 469 470 IF( nn_mld_diag > 0 ) THEN 471 WRITE(numout,*) '=============== Vertically-interpolated mixed layer ================' 472 WRITE(numout,*) '(Diagnostic number, nn_mld_type, rn_zref, rn_dT_crit, rn_iso_frac)' 473 DO jn = 1, nn_mld_diag 474 WRITE(numout,*) 'MLD criterion',jn,':' 475 WRITE(numout,*) ' nn_mld_type =', mld_diags(jn)%mld_type 476 WRITE(numout,*) ' rn_zref =' , mld_diags(jn)%zref 477 WRITE(numout,*) ' rn_dT_crit =' , mld_diags(jn)%dT_crit 478 WRITE(numout,*) ' rn_iso_frac =', mld_diags(jn)%iso_frac 479 END DO 480 WRITE(numout,*) '====================================================================' 481 ENDIF 482 ENDIF 483 484 IF( nn_mld_diag > 0 ) THEN 485 DO jn = 1, nn_mld_diag 486 WRITE(cmld,'(I1)') jn 487 IF( iom_use( "mldzint_"//cmld ) .OR. iom_use( "mldhtc_"//cmld ) ) THEN 488 CALL zdf_mxl_zint_mld( mld_diags(jn) ) 489 490 IF( iom_use( "mldzint_"//cmld ) ) THEN 491 CALL iom_put( "mldzint_"//cmld, hmld_zint(:,:) ) 492 ENDIF 493 494 IF( iom_use( "mldhtc_"//cmld ) ) THEN 495 CALL zdf_mxl_zint_htc( kt ) 496 CALL iom_put( "mldhtc_"//cmld , htc_mld(:,:) ) 497 ENDIF 498 ENDIF 499 END DO 500 ENDIF 501 502 END SUBROUTINE zdf_mxl_zint 147 503 148 504 !!======================================================================
Note: See TracChangeset
for help on using the changeset viewer.