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 1559 for trunk/NEMO/OPA_SRC/ZDF/zdfmxl.F90 – NEMO

Ignore:
Timestamp:
2009-07-29T16:03:14+02:00 (15 years ago)
Author:
ctlod
Message:

only cosmetic changes

File:
1 edited

Legend:

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

    r1482 r1559  
    44   !! Ocean physics: mixed layer depth  
    55   !!====================================================================== 
    6    !! History : 
    7    !!   9.0  !  03-08  (G. Madec)  OpenMP/autotasking optimization 
     6   !! History :  1.0  ! 2003-08  (G. Madec)  original code 
     7   !!            3.2  ! 2009-07  (S. Masson, G. Madec)  IOM + merge of DO-loop 
    88   !!---------------------------------------------------------------------- 
    99   !!   zdf_mxl      : Compute the turbocline and mixed layer depths. 
    1010   !!---------------------------------------------------------------------- 
    11    !! * Modules used 
    1211   USE oce             ! ocean dynamics and tracers variables 
    1312   USE dom_oce         ! ocean space and time domain variables 
     
    2019   PRIVATE 
    2120 
    22    !! * Routine accessibility 
    23    PUBLIC zdf_mxl           ! called by step.F90 
     21   PUBLIC   zdf_mxl    ! called by step.F90 
    2422 
    25    !! * Shared module variables 
    26    INTEGER, PUBLIC, DIMENSION(jpi,jpj) ::   &   !: 
    27       nmln                  !: number of level in the mixed layer 
    28    REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   &   !: 
    29       hmld ,             &  !: mixing layer depth (turbocline) (m) 
    30       hmlp ,             &  !: mixed layer depth  (rho=rho0+zdcrit) (m) 
    31       hmlpt                 !: mixed layer depth at t-points (m) 
     23   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   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   hmlp    !: mixed layer depth  (rho=rho0+zdcrit) [m] 
     26   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   hmlpt   !: mixed layer depth at t-points        [m] 
    3227 
    33    !! * module variables 
    34    REAL(wp) ::   & 
    35       avt_c = 5.e-4_wp,  &  ! Kz criterion for the turbocline depth 
    36       rho_c = 0.01_wp       ! density criterion for mixed layer depth 
     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 
    3730 
    3831   !! * Substitutions 
    3932#  include "domzgr_substitute.h90" 
    4033   !!---------------------------------------------------------------------- 
    41    !!   OPA 9.0 , LOCEAN-IPSL (2005)  
     34   !! NEMO/OPA 3.2 , LOCEAN-IPSL (2009)  
    4235   !! $Id$  
    43    !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt  
     36   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) 
    4437   !!---------------------------------------------------------------------- 
    4538 
     
    5144      !!                    
    5245      !! ** Purpose :   Compute the turbocline depth and the mixed layer depth 
    53       !!      with density criteria. 
     46      !!              with density criteria. 
    5447      !! 
    5548      !! ** Method  :   The turbocline depth is the depth at which the vertical 
     
    5851      !!      value defined locally (avt_c here taken equal to 5 cm/s2) 
    5952      !! 
    60       !! ** Action  : 
     53      !! ** Action  :   nmln, hmld, hmlp, hmlpt 
     54      !!---------------------------------------------------------------------- 
     55      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index 
    6156      !! 
    62       !! History : 
    63       !!        !  94-11  (M. Imbard)  Original code 
    64       !!   8.0  !  96-01  (E. Guilyardi)  sub mixed layer temp. 
    65       !!   8.1  !  97-07  (G. Madec)  optimization 
    66       !!   8.5  !  02-06  (G. Madec)  F90: Free form and module 
    67       !!---------------------------------------------------------------------- 
    68       !! * Arguments 
    69       INTEGER, INTENT( in ) ::   kt         ! ocean time-step index 
    70  
    71       !! * Local declarations 
    7257      INTEGER ::   ji, jj, jk     ! dummy loop indices 
    73       INTEGER ::   ik             ! temporary integer 
    74       INTEGER, DIMENSION(jpi,jpj) ::   & 
    75          imld                     ! temporary workspace 
     58      INTEGER ::   iki, ikn       ! temporary integer 
     59      INTEGER, DIMENSION(jpi,jpj) ::   imld   ! temporary workspace 
    7660      !!---------------------------------------------------------------------- 
    7761 
     
    8266      ENDIF 
    8367 
    84  
    85       ! 1. Turbocline depth 
    86       ! ------------------- 
    87       ! last w-level at which avt<avt_c (starting from the bottom jk=jpk) 
    88       ! (since avt(.,.,jpk)=0, we have jpk=< imld =< 2 ) 
    89       DO jk = jpk, 2, -1 
     68      ! 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 
    9072         DO jj = 1, jpj 
    9173            DO ji = 1, jpi 
    92                IF( avt(ji,jj,jk) < avt_c ) imld(ji,jj) = jk  
     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 
    9376            END DO 
    9477         END DO 
    9578      END DO 
    96  
    97       ! Turbocline depth and sub-turbocline temperature 
     79      ! depth of the mixing and mixed layers 
    9880      DO jj = 1, jpj 
    9981         DO ji = 1, jpi 
    100             ik = imld(ji,jj) 
    101             hmld (ji,jj) = fsdepw(ji,jj,ik) * tmask(ji,jj,1) 
     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 
    10287         END DO 
    10388      END DO 
    10489      CALL iom_put( "mldturb", hmld )   ! turbocline depth 
    105  
    106 !!gm idea 
    107 !!    
    108 !!gm  DO jk = jpk, 2, -1 
    109 !!gm     DO jj = 1, jpj 
    110 !!gm        DO ji = 1, jpi 
    111 !!gm           IF( avt(ji,jj,jk) < avt_c ) hmld(ji,jj) = fsdepw(ji,jj,jk) * tmask(ji,jj,1) 
    112 !!gm        END DO 
    113 !!gm     END DO 
    114 !!gm  END DO 
    115 !!gm 
    116  
    117       ! 2. Mixed layer depth 
    118       ! -------------------- 
    119       ! Initialization to the number of w ocean point mbathy 
    120       nmln(:,:) = mbathy(:,:) 
    121  
    122       ! Last w-level at which rhop>=rho surf+rho_c (starting from jpk-1) 
    123       ! (rhop defined at t-point, thus jk-1 for w-level just above) 
    124       DO jk = jpkm1, 2, -1 
    125          DO jj = 1, jpj 
    126             DO ji = 1, jpi 
    127                IF( rhop(ji,jj,jk) > rhop(ji,jj,1) + rho_c )   nmln(ji,jj) = jk 
    128             END DO 
    129          END DO 
    130       END DO 
    131  
    132       ! Mixed layer depth 
    133       DO jj = 1, jpj 
    134          DO ji = 1, jpi 
    135             ik = nmln(ji,jj) 
    136             hmlp (ji,jj) = fsdepw(ji,jj,ik) * tmask(ji,jj,1) 
    137             hmlpt(ji,jj) = fsdept(ji,jj,ik-1) 
    138          END DO 
    139       END DO 
    140       CALL iom_put( "mld010", hmlp )   ! mixed layer depth 
     90      CALL iom_put( "mld010" , hmlp )   ! mixed layer depth 
    14191       
    14292      IF(ln_ctl)   CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmld, clinfo2=' hmld : ', ovlap=1 ) 
    143  
     93      ! 
    14494   END SUBROUTINE zdf_mxl 
    14595 
Note: See TracChangeset for help on using the changeset viewer.