Changeset 1577 for trunk/NEMO/OPA_SRC/ZDF/zdfmxl.F90
- Timestamp:
- 2009-08-04T16:56:59+02:00 (15 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/NEMO/OPA_SRC/ZDF/zdfmxl.F90
r1559 r1577 5 5 !!====================================================================== 6 6 !! History : 1.0 ! 2003-08 (G. Madec) original code 7 !! 3.2 ! 2009-07 (S. Masson, G. Madec) IOM + merge of DO-loop7 !! 3.2 ! 2009-07 (S. Masson, G. Madec) IOM + 10m ref. 8 8 !!---------------------------------------------------------------------- 9 !! zdf_mxl : Compute the turbocline and mixed layer depths.9 !! zdf_mxl : Compute mixed layer depth and level 10 10 !!---------------------------------------------------------------------- 11 11 USE oce ! ocean dynamics and tracers variables 12 12 USE dom_oce ! ocean space and time domain variables 13 USE zdf_oce ! ocean vertical physics14 13 USE in_out_manager ! I/O manager 15 14 USE prtctl ! Print control … … 22 21 23 22 INTEGER , PUBLIC, DIMENSION(jpi,jpj) :: nmln !: number of level in the mixed layer 24 REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: hmld !: mixing layer depth (turbocline) [m]25 23 REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: hmlp !: mixed layer depth (rho=rho0+zdcrit) [m] 26 24 REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: hmlpt !: mixed layer depth at t-points [m] 27 28 REAL(wp) :: avt_c = 5.e-4_wp ! Kz criterion for the turbocline depth29 REAL(wp) :: rho_c = 0.01_wp ! density criterion for mixed layer depth30 25 31 26 !! * Substitutions … … 43 38 !! *** ROUTINE zdfmxl *** 44 39 !! 45 !! ** Purpose : Compute the turbocline depth and the mixed layer depth 46 !! with density criteria. 40 !! ** Purpose : the mixed layer depth with density criteria. 47 41 !! 48 !! ** Method : The turbocline depth is the depth at which the vertical 49 !! eddy diffusivity coefficient (resulting from the vertical physics 50 !! alone, not the isopycnal part, see trazdf.F) fall below a given 51 !! value defined locally (avt_c here taken equal to 5 cm/s2) 42 !! ** Method : The mixed layer depth is the shallowest W depth with 43 !! the density of the corresponding T point (just bellow) bellow a 44 !! given value defined locally as rho(10m) + zrho_c 52 45 !! 53 !! ** Action : nmln, hml d, hmlp, hmlpt46 !! ** Action : nmln, hmlp, hmlpt 54 47 !!---------------------------------------------------------------------- 55 48 INTEGER, INTENT( in ) :: kt ! ocean time-step index 56 49 !! 57 INTEGER :: ji, jj, jk! dummy loop indices58 INTEGER :: iki, ikn ! temporary integer59 INTEGER, DIMENSION(jpi,jpj) :: imld ! temporary workspace50 INTEGER :: ji, jj, jk ! dummy loop indices 51 INTEGER :: iikn ! temporary integer within a do loop 52 REAL(wp) :: zrho_c = 0.01_wp ! density criterion for mixed layer depth 60 53 !!---------------------------------------------------------------------- 61 54 … … 67 60 68 61 ! w-level of the mixing and mixed layers 69 nmln(:,:) = mbathy(:,:) ! Initialization to the number of w ocean point mbathy 70 imld(:,:) = mbathy(:,:) 71 DO jk = jpkm1, 2, -1 ! from the bottom to the surface 62 nmln(:,:) = mbathy(:,:) ! Initialization to the number of w ocean point mbathy 63 DO jk = jpkm1, nlb10, -1 ! from the bottom to nlb10 72 64 DO jj = 1, jpj 73 65 DO ji = 1, jpi 74 IF( avt (ji,jj,jk) < avt_c ) imld(ji,jj) = jk ! Turbocline 75 IF( rhop(ji,jj,jk) > rhop(ji,jj,1) + rho_c ) nmln(ji,jj) = jk ! Mixed layer 66 IF( rhop(ji,jj,jk) > rhop(ji,jj,nla10) + zrho_c ) nmln(ji,jj) = jk ! Mixed layer 76 67 END DO 77 68 END DO … … 80 71 DO jj = 1, jpj 81 72 DO ji = 1, jpi 82 iki = imld(ji,jj) 83 ikn = nmln(ji,jj) 84 hmld (ji,jj) = fsdepw(ji,jj,iki) * tmask(ji,jj,1) ! Turbocline depth 85 hmlp (ji,jj) = fsdepw(ji,jj,ikn) * tmask(ji,jj,1) ! Mixed layer depth 86 hmlpt(ji,jj) = fsdept(ji,jj,ikn-1) ! depth of the last T-point inside the mixed layer 73 iikn = nmln(ji,jj) 74 hmlp (ji,jj) = fsdepw(ji,jj,iikn ) * tmask(ji,jj,1) ! Mixed layer depth 75 hmlpt(ji,jj) = fsdept(ji,jj,iikn-1) ! depth of the last T-point inside the mixed layer 87 76 END DO 88 77 END DO 89 CALL iom_put( "mldturb", hmld ) ! turbocline depth 90 CALL iom_put( "mld010" , hmlp ) ! mixed layer depth 78 CALL iom_put( "mldr10_1" , hmlp ) ! mixed layer depth 91 79 92 IF(ln_ctl) CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hml d, clinfo2=' hmld: ', ovlap=1 )80 IF(ln_ctl) CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ', ovlap=1 ) 93 81 ! 94 82 END SUBROUTINE zdf_mxl
Note: See TracChangeset
for help on using the changeset viewer.