New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 1577 for trunk/NEMO/OPA_SRC/ZDF/zdfmxl.F90 – NEMO

Ignore:
Timestamp:
2009-08-04T16:56:59+02:00 (15 years ago)
Author:
smasson
Message:

update diahth and zdfmxl, see ticket:468

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/ZDF/zdfmxl.F90

    r1559 r1577  
    55   !!====================================================================== 
    66   !! History :  1.0  ! 2003-08  (G. Madec)  original code 
    7    !!            3.2  ! 2009-07  (S. Masson, G. Madec)  IOM + merge of DO-loop 
     7   !!            3.2  ! 2009-07  (S. Masson, G. Madec)  IOM + 10m ref. 
    88   !!---------------------------------------------------------------------- 
    9    !!   zdf_mxl      : Compute the turbocline and mixed layer depths. 
     9   !!   zdf_mxl       : Compute mixed layer depth and level 
    1010   !!---------------------------------------------------------------------- 
    1111   USE oce             ! ocean dynamics and tracers variables 
    1212   USE dom_oce         ! ocean space and time domain variables 
    13    USE zdf_oce         ! ocean vertical physics 
    1413   USE in_out_manager  ! I/O manager 
    1514   USE prtctl          ! Print control 
     
    2221 
    2322   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] 
    2523   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   hmlp    !: mixed layer depth  (rho=rho0+zdcrit) [m] 
    2624   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 depth 
    29    REAL(wp) ::   rho_c = 0.01_wp    ! density criterion for mixed layer depth 
    3025 
    3126   !! * Substitutions 
     
    4338      !!                  ***  ROUTINE zdfmxl  *** 
    4439      !!                    
    45       !! ** Purpose :   Compute the turbocline depth and the mixed layer depth 
    46       !!              with density criteria. 
     40      !! ** Purpose :   the mixed layer depth with density criteria. 
    4741      !! 
    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 
    5245      !! 
    53       !! ** Action  :   nmln, hmld, hmlp, hmlpt 
     46      !! ** Action  :   nmln, hmlp, hmlpt 
    5447      !!---------------------------------------------------------------------- 
    5548      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index 
    5649      !! 
    57       INTEGER ::   ji, jj, jk     ! dummy loop indices 
    58       INTEGER ::   iki, ikn       ! temporary integer 
    59       INTEGER, DIMENSION(jpi,jpj) ::   imld   ! temporary workspace 
     50      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 
    6053      !!---------------------------------------------------------------------- 
    6154 
     
    6760 
    6861      ! 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 
    7264         DO jj = 1, jpj 
    7365            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 
    7667            END DO 
    7768         END DO 
     
    8071      DO jj = 1, jpj 
    8172         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 
    8776         END DO 
    8877      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 
    9179       
    92       IF(ln_ctl)   CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmld, 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 ) 
    9381      ! 
    9482   END SUBROUTINE zdf_mxl 
Note: See TracChangeset for help on using the changeset viewer.