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 1585 for trunk/NEMO – NEMO

Changeset 1585 for trunk/NEMO


Ignore:
Timestamp:
2009-08-06T09:29:43+02:00 (15 years ago)
Author:
smasson
Message:

move back turbocline computation in zdfmxl, see ticket:517

Location:
trunk/NEMO/OPA_SRC
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/DIA/diahth.F90

    r1577 r1585  
    2121   USE oce             ! ocean dynamics and tracers 
    2222   USE dom_oce         ! ocean space and time domain 
    23    USE zdf_oce         ! ocean vertical physics 
    2423   USE phycst          ! physical constants 
    2524   USE in_out_manager  ! I/O manager 
     
    3534   LOGICAL , PUBLIC, PARAMETER          ::   lk_diahth = .TRUE.   !: thermocline-20d depths flag 
    3635   ! note: following variables should move to local variables once iom_put is always used  
    37    REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   hmld                 !: mixing layer depth (turbocline)                [m] 
    3836   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   hth                  !: depth of the max vertical temperature gradient [m] 
    3937   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   hd20                 !: depth of 20 C isotherm                         [m] 
     
    112110            zztmp = bathy(ji,jj) 
    113111            hth     (ji,jj) = zztmp 
    114             hmld    (ji,jj) = zztmp 
    115112            zabs2   (ji,jj) = zztmp 
    116113            ztm2    (ji,jj) = zztmp 
     
    165162               ENDIF 
    166163                
    167                IF( avt (ji,jj,jk) < zavt5 )    hmld   (ji,jj) = zzdep                ! avt < zavt5 
    168  
    169164               IF( nla10 > 1 ) THEN  
    170165                  zztmp = rhop(ji,jj,jk) - rhop(ji,jj,1)                             ! delta rho(1) 
     
    178173       
    179174      CALL iom_put( "mlddzt", hth )            ! depth of the thermocline 
    180       CALL iom_put( "mldkz5", hmld )           ! turbocline depth 
    181175      IF( nla10 > 1 ) THEN  
    182176         CALL iom_put( "mldr0_3", zrho0_3 )   ! MLD delta rho(surf) = 0.03 
  • trunk/NEMO/OPA_SRC/DIA/diawri.F90

    r1581 r1585  
    300300         CALL histdef( nid_T, "soshfldo", "Shortwave Radiation"                , "W/m2"   ,   &  ! qsr 
    301301            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     302         CALL histdef( nid_T, "somixhgt", "Turbocline Depth"                   , "m"      ,   &  ! hmld 
     303            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    302304         CALL histdef( nid_T, "somxl010", "Mixed Layer Depth 0.01"             , "m"      ,   &  ! hmlp 
    303305            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    331333            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clmx, zsto, zout ) 
    332334#if defined key_diahth 
    333          CALL histdef( nid_T, "somixhgt", "Turbocline Depth"                   , "m"      ,   &  ! hmld 
    334             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    335335         CALL histdef( nid_T, "sothedep", "Thermocline Depth"                  , "m"      ,   & ! hth 
    336336            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    445445      CALL histwrite( nid_T, "sohefldo", it, qns + qsr     , ndim_hT, ndex_hT )   ! total heat flux 
    446446      CALL histwrite( nid_T, "soshfldo", it, qsr           , ndim_hT, ndex_hT )   ! solar heat flux 
     447      CALL histwrite( nid_T, "somixhgt", it, hmld          , ndim_hT, ndex_hT )   ! turbocline depth 
    447448      CALL histwrite( nid_T, "somxl010", it, hmlp          , ndim_hT, ndex_hT )   ! mixed layer depth 
    448449      CALL histwrite( nid_T, "soicecov", it, fr_i          , ndim_hT, ndex_hT )   ! ice fraction    
     
    466467 
    467468#if defined key_diahth 
    468       CALL histwrite( nid_T, "somixhgt", it, hmld          , ndim_hT, ndex_hT )   ! turbocline depth 
    469469      CALL histwrite( nid_T, "sothedep", it, hth           , ndim_hT, ndex_hT )   ! depth of the thermocline 
    470470      CALL histwrite( nid_T, "so20chgt", it, hd20          , ndim_hT, ndex_hT )   ! depth of the 20 isotherm 
  • trunk/NEMO/OPA_SRC/DIA/diawri_dimg.h90

    r1577 r1585  
    173173       fsel(:,:,8 ) = fsel(:,:,8 ) + qrp (:,:) 
    174174       fsel(:,:,9 ) = fsel(:,:,9 ) + erp (:,:) 
    175 #ifdef key_diahth    
    176175       fsel(:,:,10) = fsel(:,:,10) + hmld(:,:) 
    177 #endif 
    178176       fsel(:,:,11) = fsel(:,:,11) + hmlp(:,:) 
    179177       fsel(:,:,12) = fsel(:,:,12) + fr_i(:,:) 
     
    251249          fsel(:,:,8 ) = qrp (:,:) * tmask(:,:,1) 
    252250          fsel(:,:,9 ) = erp (:,:) * tmask(:,:,1) 
    253 #ifdef key_diahth            
    254251          fsel(:,:,10) = hmld(:,:) * tmask(:,:,1) 
    255 #endif 
    256252          fsel(:,:,11) = hmlp(:,:) * tmask(:,:,1) 
    257253          fsel(:,:,12) = fr_i(:,:) * tmask(:,:,1) 
  • trunk/NEMO/OPA_SRC/ZDF/zdfmxl.F90

    r1577 r1585  
    55   !!====================================================================== 
    66   !! History :  1.0  ! 2003-08  (G. Madec)  original code 
    7    !!            3.2  ! 2009-07  (S. Masson, G. Madec)  IOM + 10m ref. 
     7   !!            3.2  ! 2009-07  (S. Masson, G. Madec)  IOM + merge of DO-loop 
    88   !!---------------------------------------------------------------------- 
    9    !!   zdf_mxl       : Compute mixed layer depth and level 
     9   !!   zdf_mxl      : Compute the turbocline and mixed layer depths. 
    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 
    1314   USE in_out_manager  ! I/O manager 
    1415   USE prtctl          ! Print control 
     
    2021   PUBLIC   zdf_mxl    ! called by step.F90 
    2122 
    22    INTEGER , PUBLIC, DIMENSION(jpi,jpj) ::   nmln    !: number of level in the mixed layer 
     23   INTEGER , PUBLIC, DIMENSION(jpi,jpj) ::   nmln    !: number of level in the mixed layer (used by TOP) 
     24   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   hmld    !: mixing layer depth (turbocline)      [m] 
    2325   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   hmlp    !: mixed layer depth  (rho=rho0+zdcrit) [m] 
    2426   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   hmlpt   !: mixed layer depth at t-points        [m] 
     
    3840      !!                  ***  ROUTINE zdfmxl  *** 
    3941      !!                    
    40       !! ** Purpose :   the mixed layer depth with density criteria. 
     42      !! ** Purpose :   Compute the turbocline depth and the mixed layer depth 
     43      !!              with density criteria. 
    4144      !! 
    4245      !! ** Method  :   The mixed layer depth is the shallowest W depth with  
    4346      !!      the density of the corresponding T point (just bellow) bellow a 
    4447      !!      given value defined locally as rho(10m) + zrho_c 
     48      !!               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) 
    4552      !! 
    46       !! ** Action  :   nmln, hmlp, hmlpt 
     53      !! ** Action  :   nmln, hmld, hmlp, hmlpt 
    4754      !!---------------------------------------------------------------------- 
    4855      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index 
    4956      !! 
    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 
     57      INTEGER                     ::   ji, jj, jk          ! dummy loop indices 
     58      INTEGER                     ::   iikn, iiki          ! temporary integer within a do loop 
     59      INTEGER, DIMENSION(jpi,jpj) ::   imld                ! temporary workspace 
     60      REAL(wp)                    ::   zrho_c = 0.01_wp    ! density criterion for mixed layer depth 
     61      REAL(wp)                    ::   zavt_c = 5.e-4_wp   ! Kz criterion for the turbocline depth 
    5362      !!---------------------------------------------------------------------- 
    5463 
     
    6170      ! w-level of the mixing and mixed layers 
    6271      nmln(:,:) = mbathy(:,:)          ! Initialization to the number of w ocean point mbathy 
     72      imld(:,:) = mbathy(:,:) 
    6373      DO jk = jpkm1, nlb10, -1         ! from the bottom to nlb10 
    6474         DO jj = 1, jpj 
    6575            DO ji = 1, jpi 
    6676               IF( rhop(ji,jj,jk) > rhop(ji,jj,nla10) + zrho_c )   nmln(ji,jj) = jk      ! Mixed layer 
     77               IF( avt (ji,jj,jk) < zavt_c                     )   imld(ji,jj) = jk      ! Turbocline  
    6778            END DO 
    6879         END DO 
     
    7182      DO jj = 1, jpj 
    7283         DO ji = 1, jpi 
     84            iiki = imld(ji,jj) 
    7385            iikn = nmln(ji,jj) 
     86            hmld (ji,jj) = fsdepw(ji,jj,iiki  ) * tmask(ji,jj,1)    ! Turbocline depth  
    7487            hmlp (ji,jj) = fsdepw(ji,jj,iikn  ) * tmask(ji,jj,1)    ! Mixed layer depth 
    7588            hmlpt(ji,jj) = fsdept(ji,jj,iikn-1)                     ! depth of the last T-point inside the mixed layer 
    7689         END DO 
    7790      END DO 
    78       CALL iom_put( "mldr10_1" , hmlp )   ! mixed layer depth 
     91      CALL iom_put( "mldr10_1", hmlp )   ! mixed layer depth 
     92      CALL iom_put( "mldkz5"  , hmld )   ! turbocline depth 
    7993       
    8094      IF(ln_ctl)   CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ', ovlap=1 ) 
Note: See TracChangeset for help on using the changeset viewer.