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 10897 – NEMO

Changeset 10897


Ignore:
Timestamp:
2019-04-26T12:17:40+02:00 (5 years ago)
Author:
davestorkey
Message:

UKMO/NEMO_4.0_GO8_package branch: Copy over changes from previous version of the branch.

  1. Implement datestamping of restarts.
  2. Local diagnostics, particularly vertically-interpolated MLD diagnostics.
  3. Restrict write-out in tra_bbl_init with IF(lwp).
Location:
NEMO/branches/UKMO/NEMO_4.0_GO8_package
Files:
8 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/UKMO/NEMO_4.0_GO8_package/cfgs/SHARED/field_def_nemo-oce.xml

    r10823 r10897  
    6262        <field id="mldr10_1max"  long_name="Max of Mixed Layer Depth (dsigma = 0.01 wrt 10m)"   field_ref="mldr10_1"   operation="maximum"                                                                          /> 
    6363        <field id="mldr10_1min"  long_name="Min of Mixed Layer Depth (dsigma = 0.01 wrt 10m)"   field_ref="mldr10_1"   operation="minimum"                                                                          /> 
     64         <field id="mldzint_1"    long_name="Mixed Layer Depth interpolated"                     standard_name="ocean_mixed_layer_thickness"                                                       unit="m"          /> 
     65         <field id="mldzint_2"    long_name="Mixed Layer Depth interpolated"                     standard_name="ocean_mixed_layer_thickness"                                                       unit="m"          /> 
     66         <field id="mldzint_3"    long_name="Mixed Layer Depth interpolated"                     standard_name="ocean_mixed_layer_thickness"                                                       unit="m"          /> 
     67         <field id="mldzint_4"    long_name="Mixed Layer Depth interpolated"                     standard_name="ocean_mixed_layer_thickness"                                                       unit="m"          /> 
     68         <field id="mldzint_5"    long_name="Mixed Layer Depth interpolated"                     standard_name="ocean_mixed_layer_thickness"                                                       unit="m"          /> 
     69         <field id="mldhtc_1"     long_name="Mixed Layer Depth integrated heat content"          standard_name="integral_of_sea_water_potential_temperature_wrt_depth_expressed_as_heat_content"   unit="J/m2"       /> 
     70         <field id="mldhtc_2"     long_name="Mixed Layer Depth integrated heat content"          standard_name="integral_of_sea_water_potential_temperature_wrt_depth_expressed_as_heat_content"   unit="J/m2"       /> 
     71         <field id="mldhtc_3"     long_name="Mixed Layer Depth integrated heat content"          standard_name="integral_of_sea_water_potential_temperature_wrt_depth_expressed_as_heat_content"   unit="J/m2"       /> 
     72         <field id="mldhtc_4"     long_name="Mixed Layer Depth integrated heat content"          standard_name="integral_of_sea_water_potential_temperature_wrt_depth_expressed_as_heat_content"   unit="J/m2"       /> 
     73         <field id="mldhtc_5"     long_name="Mixed Layer Depth integrated heat content"          standard_name="integral_of_sea_water_potential_temperature_wrt_depth_expressed_as_heat_content"   unit="J/m2"       /> 
    6474        <field id="heatc"        long_name="Heat content vertically integrated"                 standard_name="integral_of_sea_water_potential_temperature_wrt_depth_expressed_as_heat_content"   unit="J/m2"       /> 
    6575        <field id="saltc"        long_name="Salt content vertically integrated"                                                                                                                   unit="1e-3*kg/m2" /> 
     
    357367        <field id="uoce"         long_name="ocean current along i-axis"                             standard_name="sea_water_x_velocity"        unit="m/s"        grid_ref="grid_U_3D" /> 
    358368        <field id="uoce_e3u"     long_name="ocean current along i-axis  (thickness weighted)"                                                   unit="m/s"        grid_ref="grid_U_3D"  > uoce * e3u </field> 
     369        <field id="uoce2_e3u"    long_name="ocean current along i-axis squared (thickness weighted)"                                            unit="m3/s2"      grid_ref="grid_U_3D"  > uoce * uoce * e3u </field> 
    359370        <field id="ssu"          long_name="ocean surface current along i-axis"                                                                 unit="m/s"                             /> 
    360371        <field id="sbu"          long_name="ocean bottom current along i-axis"                                                                  unit="m/s"                             /> 
     
    411422        <field id="voce"         long_name="ocean current along j-axis"                             standard_name="sea_water_y_velocity"        unit="m/s"        grid_ref="grid_V_3D" /> 
    412423        <field id="voce_e3v"     long_name="ocean current along j-axis  (thickness weighted)"                                                   unit="m/s"        grid_ref="grid_V_3D"  > voce * e3v </field> 
     424        <field id="voce2_e3v"    long_name="ocean current along j-axis squared (thickness weighted)"                                            unit="m3/s2"      grid_ref="grid_V_3D"  > voce * voce * e3v </field> 
    413425        <field id="ssv"          long_name="ocean surface current along j-axis"                                                                 unit="m/s"                             /> 
    414426        <field id="sbv"          long_name="ocean bottom current along j-axis"                                                                  unit="m/s"                             /> 
  • NEMO/branches/UKMO/NEMO_4.0_GO8_package/src/ICE/icerst.F90

    r10888 r10897  
    2121   USE in_out_manager ! I/O manager 
    2222   USE iom            ! I/O manager library 
     23   USE ioipsl, ONLY : ju2ymds    ! for calendar 
    2324   USE lib_mpp        ! MPP library 
    2425   USE lib_fortran    ! fortran utilities (glob_sum + no signed zero) 
     
    4647      INTEGER, INTENT(in) ::   kt       ! number of iteration 
    4748      ! 
     49      INTEGER             ::   iyear, imonth, iday 
     50      REAL (wp)           ::   zsec 
     51      REAL (wp)           ::   zfjulday 
    4852      CHARACTER(len=20)   ::   clkt     ! ocean time-step define as a character 
    4953      CHARACTER(len=50)   ::   clname   ! ice output restart file name 
     
    6064         IF( nitrst <= nitend .AND. nitrst > 0 ) THEN 
    6165            ! beware of the format used to write kt (default is i8.8, that should be large enough...) 
    62             IF( nitrst > 99999999 ) THEN   ;   WRITE(clkt, *       ) nitrst 
    63             ELSE                           ;   WRITE(clkt, '(i8.8)') nitrst 
     66            IF ( ln_rstdate ) THEN 
     67               zfjulday = fjulday + (2*nn_fsbc+1)*rdt / rday 
     68               IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday )   zfjulday = REAL(NINT(zfjulday),wp)   ! avoid truncation error 
     69               CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec )            
     70               WRITE(clkt, '(i4.4,2i2.2)') iyear, imonth, iday 
     71            ELSE 
     72               IF( nitrst > 99999999 ) THEN   ;   WRITE(clkt, *       ) nitrst 
     73               ELSE                           ;   WRITE(clkt, '(i8.8)') nitrst 
     74               ENDIF 
    6475            ENDIF 
    6576            ! create the file 
  • NEMO/branches/UKMO/NEMO_4.0_GO8_package/src/OCE/DOM/domain.F90

    r10888 r10897  
    292292         &             nn_it000, nn_itend , nn_date0    , nn_time0     , nn_leapy  , nn_istate ,     & 
    293293         &             nn_stock, nn_write , ln_mskland  , ln_clobber   , nn_chunksz, nn_euler  ,     & 
    294          &             ln_cfmeta, ln_iscpl, ln_xios_read, nn_wxios 
     294         &             ln_cfmeta, ln_iscpl, ln_xios_read, nn_wxios, ln_rstdate 
    295295      NAMELIST/namdom/ ln_linssh, rn_isfhmin, rn_rdt, rn_atfp, ln_crs, ln_meshmask 
    296296#if defined key_netcdf4 
     
    338338         WRITE(numout,*) '      frequency of output file        nn_write        = ', nn_write 
    339339         WRITE(numout,*) '      mask land points                ln_mskland      = ', ln_mskland 
     340         WRITE(numout,*) '      date-stamp restart files        ln_rstdate = ', ln_rstdate 
    340341         WRITE(numout,*) '      additional CF standard metadata ln_cfmeta       = ', ln_cfmeta 
    341342         WRITE(numout,*) '      overwrite an existing file      ln_clobber      = ', ln_clobber 
  • NEMO/branches/UKMO/NEMO_4.0_GO8_package/src/OCE/ICB/icbrst.F90

    r10888 r10897  
    2525   USE netcdf         ! netcdf routines for IO 
    2626   USE iom 
     27   USE ioipsl, ONLY : ju2ymds    ! for calendar 
    2728   USE icb_oce        ! define iceberg arrays 
    2829   USE icbutl         ! iceberg utility routines 
     
    190191      INTEGER ::   jn   ! dummy loop index 
    191192      INTEGER ::   ix_dim, iy_dim, ik_dim, in_dim 
    192       CHARACTER(len=256)     :: cl_path 
    193       CHARACTER(len=256)     :: cl_filename 
     193      INTEGER             ::   iyear, imonth, iday 
     194      REAL (wp)           ::   zsec 
     195      REAL (wp)           ::   zfjulday 
     196      CHARACTER(len=256)  :: cl_path 
     197      CHARACTER(len=256)  :: cl_filename 
     198      CHARACTER(LEN=20)   ::   clkt     ! ocean time-step deine as a character 
    194199      TYPE(iceberg), POINTER :: this 
    195200      TYPE(point)  , POINTER :: pt 
     
    206211         cl_path = TRIM(cn_ocerst_outdir) 
    207212         IF( cl_path(LEN_TRIM(cl_path):) /= '/' ) cl_path = TRIM(cl_path) // '/' 
     213         IF ( ln_rstdate ) THEN 
     214            zfjulday = fjulday + rdt / rday 
     215            IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday )   zfjulday = REAL(NINT(zfjulday),wp)   ! avoid truncation error 
     216            CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec )            
     217            WRITE(clkt, '(i4.4,2i2.2)') iyear, imonth, iday 
     218         ELSE 
     219            IF( kt > 999999999 ) THEN   ;   WRITE(clkt, *       ) kt 
     220            ELSE                        ;   WRITE(clkt, '(i8.8)') kt 
     221            ENDIF 
     222         ENDIF 
    208223         IF( lk_mpp ) THEN 
    209             WRITE(cl_filename,'(A,"_icebergs_",I8.8,"_restart_",I4.4,".nc")') TRIM(cexper), kt, narea-1 
     224            WRITE(cl_filename,'(A,"_icebergs_",A,"_restart_",I4.4,".nc")') TRIM(cexper), TRIM(ADJUSTL(clkt)), narea-1 
    210225         ELSE 
    211             WRITE(cl_filename,'(A,"_icebergs_",I8.8,"_restart.nc")') TRIM(cexper), kt 
     226            WRITE(cl_filename,'(A,"_icebergs_",A,"_restart.nc")') TRIM(cexper), TRIM(ADJUSTL(clkt)) 
    212227         ENDIF 
    213228         IF ( lwp .AND. nn_verbose_level >= 0) WRITE(numout,'(2a)') 'icebergs, write_restart: creating ',  & 
  • NEMO/branches/UKMO/NEMO_4.0_GO8_package/src/OCE/IOM/in_out_manager.F90

    r10888 r10897  
    4040   INTEGER, DIMENSION(10) :: nn_stocklist  !: restart dump times 
    4141   LOGICAL       ::   ln_mskland       !: mask land points in NetCDF outputs (costly: + ~15%) 
     42   LOGICAL       ::   ln_rstdate       !: T=> stamp output restart files with date instead of timestep 
    4243   LOGICAL       ::   ln_cfmeta        !: output additional data to netCDF files required for compliance with the CF metadata standard 
    4344   LOGICAL       ::   ln_clobber       !: clobber (overwrite) an existing file 
  • NEMO/branches/UKMO/NEMO_4.0_GO8_package/src/OCE/IOM/restart.F90

    r10888 r10897  
    2727   USE in_out_manager  ! I/O manager 
    2828   USE iom             ! I/O module 
     29   USE ioipsl, ONLY : ju2ymds    ! for calendar 
    2930   USE diurnal_bulk 
    3031   USE lib_mpp         ! distribued memory computing library 
     
    5960      INTEGER, INTENT(in) ::   kt     ! ocean time-step 
    6061      !! 
     62      INTEGER             ::   iyear, imonth, iday 
     63      REAL (wp)           ::   zsec 
     64      REAL (wp)           ::   zfjulday 
    6165      CHARACTER(LEN=20)   ::   clkt     ! ocean time-step deine as a character 
    6266      CHARACTER(LEN=50)   ::   clname   ! ocean output restart file name 
     
    8892         IF( nitrst <= nitend .AND. nitrst > 0 ) THEN  
    8993            ! beware of the format used to write kt (default is i8.8, that should be large enough...) 
    90             IF( nitrst > 999999999 ) THEN   ;   WRITE(clkt, *       ) nitrst 
    91             ELSE                            ;   WRITE(clkt, '(i8.8)') nitrst 
     94            IF ( ln_rstdate ) THEN 
     95               zfjulday = fjulday + rdt / rday 
     96               IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday )   zfjulday = REAL(NINT(zfjulday),wp)   ! avoid truncation error 
     97               CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec )            
     98               WRITE(clkt, '(i4.4,2i2.2)') iyear, imonth, iday 
     99            ELSE 
     100               IF( nitrst > 999999999 ) THEN   ;   WRITE(clkt, *       ) nitrst 
     101               ELSE                            ;   WRITE(clkt, '(i8.8)') nitrst 
     102               ENDIF 
    92103            ENDIF 
    93104            ! create the file 
  • NEMO/branches/UKMO/NEMO_4.0_GO8_package/src/OCE/TRA/trabbl.F90

    r10888 r10897  
    513513      IF( tra_bbl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'tra_bbl_init : unable to allocate arrays' ) 
    514514      ! 
    515       IF( nn_bbl_adv == 1 )    WRITE(numout,*) '       * Advective BBL using upper velocity' 
    516       IF( nn_bbl_adv == 2 )    WRITE(numout,*) '       * Advective BBL using velocity = F( delta rho)' 
     515      IF(lwp) THEN 
     516         IF( nn_bbl_adv == 1 )    WRITE(numout,*) '       * Advective BBL using upper velocity' 
     517         IF( nn_bbl_adv == 2 )    WRITE(numout,*) '       * Advective BBL using velocity = F( delta rho)' 
     518      ENDIF 
    517519      ! 
    518520      !                             !* vertical index of  "deep" bottom u- and v-points 
  • NEMO/branches/UKMO/NEMO_4.0_GO8_package/src/OCE/ZDF/zdfmxl.F90

    r10888 r10897  
    1515   USE trc_oce  , ONLY: l_offline         ! ocean space and time domain variables 
    1616   USE zdf_oce        ! ocean vertical physics 
     17   USE eosbn2         ! for zdf_mxl_zint 
    1718   ! 
    1819   USE in_out_manager ! I/O manager 
     
    3132   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmlp    !: mixed layer depth  (rho=rho0+zdcrit) [m]   (used by LDF) 
    3233   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   hmlpt   !: depth of the last T-point inside the mixed layer [m] (used by LDF) 
     34   REAL(wp), PUBLIC, ALLOCATABLE,       DIMENSION(:,:) ::   hmld_zint  !: vertically-interpolated mixed layer depth   [m] 
     35   REAL(wp), PUBLIC, ALLOCATABLE,       DIMENSION(:,:) ::   htc_mld    ! Heat content of hmld_zint 
     36   LOGICAL, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)    :: ll_found   ! Is T_b to be found by interpolation ? 
     37   LOGICAL, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)  :: ll_belowml ! Flag points below mixed layer when ll_found=F 
    3338 
    3439   REAL(wp), PUBLIC ::   rho_c = 0.01_wp    !: density criterion for mixed layer depth 
    3540   REAL(wp), PUBLIC ::   avt_c = 5.e-4_wp   ! Kz criterion for the turbocline depth 
     41 
     42   TYPE, PUBLIC :: MXL_ZINT   !: Structure for MLD defs 
     43      INTEGER   :: mld_type   ! mixed layer type      
     44      REAL(wp)  :: zref       ! depth of initial T_ref 
     45      REAL(wp)  :: dT_crit    ! Critical temp diff 
     46      REAL(wp)  :: iso_frac   ! Fraction of rn_dT_crit  
     47   END TYPE MXL_ZINT 
    3648 
    3749   !!---------------------------------------------------------------------- 
     
    4860      zdf_mxl_alloc = 0      ! set to zero if no array to be allocated 
    4961      IF( .NOT. ALLOCATED( nmln ) ) THEN 
    50          ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), STAT= zdf_mxl_alloc ) 
     62         ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), hmld_zint(jpi,jpj),     & 
     63   &          htc_mld(jpi,jpj), ll_found(jpi,jpj), ll_belowml(jpi,jpj,jpk), STAT= zdf_mxl_alloc )          
    5164         ! 
    5265         CALL mpp_sum ( 'zdfmxl', zdf_mxl_alloc ) 
     
    137150      ENDIF 
    138151      ! 
     152      ! Vertically-interpolated mixed-layer depth diagnostic 
     153      CALL zdf_mxl_zint( kt ) 
     154      ! 
    139155      IF(ln_ctl)   CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ' ) 
    140156      ! 
    141157   END SUBROUTINE zdf_mxl 
     158 
     159   SUBROUTINE zdf_mxl_zint_mld( sf )  
     160      !!----------------------------------------------------------------------------------  
     161      !!                    ***  ROUTINE zdf_mxl_zint_mld  ***  
     162      !                                                                         
     163      !   Calculate vertically-interpolated mixed layer depth diagnostic.  
     164      !             
     165      !   This routine can calculate the mixed layer depth diagnostic suggested by 
     166      !   Kara et al, 2000, JGR, 105, 16803, but is more general and can calculate 
     167      !   vertically-interpolated mixed-layer depth diagnostics with other parameter 
     168      !   settings set in the namzdf_mldzint namelist.   
     169      !  
     170      !   If mld_type=1 the mixed layer depth is calculated as the depth at which the   
     171      !   density has increased by an amount equivalent to a temperature difference of   
     172      !   0.8C at the surface.  
     173      !  
     174      !   For other values of mld_type the mixed layer is calculated as the depth at   
     175      !   which the temperature differs by 0.8C from the surface temperature.   
     176      !                                                                         
     177      !   David Acreman, Daley Calvert                                       
     178      !  
     179      !!-----------------------------------------------------------------------------------  
     180 
     181      TYPE(MXL_ZINT), INTENT(in)  :: sf 
     182 
     183      ! Diagnostic criteria 
     184      INTEGER   :: nn_mld_type   ! mixed layer type      
     185      REAL(wp)  :: rn_zref       ! depth of initial T_ref 
     186      REAL(wp)  :: rn_dT_crit    ! Critical temp diff 
     187      REAL(wp)  :: rn_iso_frac   ! Fraction of rn_dT_crit used 
     188 
     189      ! Local variables 
     190      REAL(wp), PARAMETER :: zepsilon = 1.e-30          ! local small value 
     191      INTEGER, DIMENSION(jpi,jpj) :: ikmt          ! number of active tracer levels  
     192      INTEGER, DIMENSION(jpi,jpj) :: ik_ref        ! index of reference level  
     193      INTEGER, DIMENSION(jpi,jpj) :: ik_iso        ! index of last uniform temp level  
     194      REAL, DIMENSION(jpi,jpj,jpk)  :: zT            ! Temperature or density  
     195      REAL, DIMENSION(jpi,jpj)    :: ppzdep        ! depth for use in calculating d(rho)  
     196      REAL, DIMENSION(jpi,jpj)    :: zT_ref        ! reference temperature  
     197      REAL    :: zT_b                                   ! base temperature  
     198      REAL, DIMENSION(jpi,jpj,jpk)  :: zdTdz         ! gradient of zT  
     199      REAL, DIMENSION(jpi,jpj,jpk)  :: zmoddT        ! Absolute temperature difference  
     200      REAL    :: zdz                                    ! depth difference  
     201      REAL    :: zdT                                    ! temperature difference  
     202      REAL, DIMENSION(jpi,jpj)    :: zdelta_T      ! difference critereon  
     203      REAL, DIMENSION(jpi,jpj)    :: zRHO1, zRHO2  ! Densities  
     204      INTEGER :: ji, jj, jk                             ! loop counter  
     205 
     206      !!-------------------------------------------------------------------------------------  
     207      !   
     208      ! Unpack structure 
     209      nn_mld_type = sf%mld_type 
     210      rn_zref     = sf%zref 
     211      rn_dT_crit  = sf%dT_crit 
     212      rn_iso_frac = sf%iso_frac 
     213 
     214      ! Set the mixed layer depth criterion at each grid point  
     215      IF( nn_mld_type == 0 ) THEN 
     216         zdelta_T(:,:) = rn_dT_crit 
     217         zT(:,:,:) = rhop(:,:,:) 
     218      ELSE IF( nn_mld_type == 1 ) THEN 
     219         ppzdep(:,:)=0.0  
     220         call eos ( tsn(:,:,1,:), ppzdep(:,:), zRHO1(:,:) )  
     221! Use zT temporarily as a copy of tsn with rn_dT_crit added to SST  
     222! [assumes number of tracers less than number of vertical levels]  
     223         zT(:,:,1:jpts)=tsn(:,:,1,1:jpts)  
     224         zT(:,:,jp_tem)=zT(:,:,1)+rn_dT_crit  
     225         CALL eos( zT(:,:,1:jpts), ppzdep(:,:), zRHO2(:,:) )  
     226         zdelta_T(:,:) = abs( zRHO1(:,:) - zRHO2(:,:) ) * rau0  
     227         ! RHO from eos (2d version) doesn't calculate north or east halo:  
     228         CALL lbc_lnk( 'zdfmxl', zdelta_T, 'T', 1. )  
     229         zT(:,:,:) = rhop(:,:,:)  
     230      ELSE  
     231         zdelta_T(:,:) = rn_dT_crit                       
     232         zT(:,:,:) = tsn(:,:,:,jp_tem)                            
     233      END IF  
     234 
     235      ! Calculate the gradient of zT and absolute difference for use later  
     236      DO jk = 1 ,jpk-2  
     237         zdTdz(:,:,jk)  =    ( zT(:,:,jk+1) - zT(:,:,jk) ) / e3w_n(:,:,jk+1)  
     238         zmoddT(:,:,jk) = abs( zT(:,:,jk+1) - zT(:,:,jk) )  
     239      END DO  
     240 
     241      ! Find density/temperature at the reference level (Kara et al use 10m).           
     242      ! ik_ref is the index of the box centre immediately above or at the reference level  
     243      ! Find rn_zref in the array of model level depths and find the ref     
     244      ! density/temperature by linear interpolation.                                    
     245      DO jk = jpkm1, 2, -1  
     246         WHERE ( gdept_n(:,:,jk) > rn_zref )  
     247           ik_ref(:,:) = jk - 1  
     248           zT_ref(:,:) = zT(:,:,jk-1) + zdTdz(:,:,jk-1) * ( rn_zref - gdept_n(:,:,jk-1) )  
     249         END WHERE  
     250      END DO  
     251 
     252      ! If the first grid box centre is below the reference level then use the  
     253      ! top model level to get zT_ref  
     254      WHERE ( gdept_n(:,:,1) > rn_zref )   
     255         zT_ref = zT(:,:,1)  
     256         ik_ref = 1  
     257      END WHERE  
     258 
     259      ! The number of active tracer levels is 1 less than the number of active w levels  
     260      ikmt(:,:) = mbkt(:,:) - 1  
     261 
     262      ! Initialize / reset 
     263      ll_found(:,:) = .false. 
     264 
     265      IF ( rn_iso_frac - zepsilon > 0. ) THEN 
     266         ! Search for a uniform density/temperature region where adjacent levels           
     267         ! differ by less than rn_iso_frac * deltaT.                                       
     268         ! ik_iso is the index of the last level in the uniform layer   
     269         ! ll_found indicates whether the mixed layer depth can be found by interpolation  
     270         ik_iso(:,:)   = ik_ref(:,:)  
     271         DO jj = 1, nlcj  
     272            DO ji = 1, nlci  
     273!CDIR NOVECTOR  
     274               DO jk = ik_ref(ji,jj), ikmt(ji,jj)-1  
     275                  IF ( zmoddT(ji,jj,jk) > ( rn_iso_frac * zdelta_T(ji,jj) ) ) THEN  
     276                     ik_iso(ji,jj)   = jk  
     277                     ll_found(ji,jj) = ( zmoddT(ji,jj,jk) > zdelta_T(ji,jj) )  
     278                     EXIT  
     279                  END IF  
     280               END DO  
     281            END DO  
     282         END DO  
     283 
     284         ! Use linear interpolation to find depth of mixed layer base where possible  
     285         hmld_zint(:,:) = rn_zref  
     286         DO jj = 1, jpj  
     287            DO ji = 1, jpi  
     288               IF (ll_found(ji,jj) .and. tmask(ji,jj,1) == 1.0) THEN  
     289                  zdz =  abs( zdelta_T(ji,jj) / zdTdz(ji,jj,ik_iso(ji,jj)) )  
     290                  hmld_zint(ji,jj) = gdept_n(ji,jj,ik_iso(ji,jj)) + zdz  
     291               END IF  
     292            END DO  
     293         END DO  
     294      END IF 
     295 
     296      ! If ll_found = .false. then calculate MLD using difference of zdelta_T     
     297      ! from the reference density/temperature  
     298  
     299! Prevent this section from working on land points  
     300      WHERE ( tmask(:,:,1) /= 1.0 )  
     301         ll_found = .true.  
     302      END WHERE  
     303  
     304      DO jk=1, jpk  
     305         ll_belowml(:,:,jk) = abs( zT(:,:,jk) - zT_ref(:,:) ) >= zdelta_T(:,:)   
     306      END DO  
     307  
     308! Set default value where interpolation cannot be used (ll_found=false)   
     309      DO jj = 1, jpj  
     310         DO ji = 1, jpi  
     311            IF ( .not. ll_found(ji,jj) )  hmld_zint(ji,jj) = gdept_n(ji,jj,ikmt(ji,jj))  
     312         END DO  
     313      END DO  
     314 
     315      DO jj = 1, jpj  
     316         DO ji = 1, jpi  
     317!CDIR NOVECTOR  
     318            DO jk = ik_ref(ji,jj)+1, ikmt(ji,jj)  
     319               IF ( ll_found(ji,jj) ) EXIT  
     320               IF ( ll_belowml(ji,jj,jk) ) THEN                 
     321                  zT_b = zT_ref(ji,jj) + zdelta_T(ji,jj) * SIGN(1.0, zdTdz(ji,jj,jk-1) )  
     322                  zdT  = zT_b - zT(ji,jj,jk-1)                                       
     323                  zdz  = zdT / zdTdz(ji,jj,jk-1)                                        
     324                  hmld_zint(ji,jj) = gdept_n(ji,jj,jk-1) + zdz  
     325                  EXIT                                                    
     326               END IF  
     327            END DO  
     328         END DO  
     329      END DO  
     330 
     331      hmld_zint(:,:) = hmld_zint(:,:)*tmask(:,:,1)  
     332      !   
     333   END SUBROUTINE zdf_mxl_zint_mld 
     334 
     335   SUBROUTINE zdf_mxl_zint_htc( kt ) 
     336      !!---------------------------------------------------------------------- 
     337      !!                  ***  ROUTINE zdf_mxl_zint_htc  *** 
     338      !!  
     339      !! ** Purpose :    
     340      !! 
     341      !! ** Method  :    
     342      !!---------------------------------------------------------------------- 
     343 
     344      INTEGER, INTENT(in) ::   kt   ! ocean time-step index 
     345 
     346      INTEGER :: ji, jj, jk 
     347      INTEGER :: ikmax 
     348      REAL(wp) :: zc, zcoef 
     349      ! 
     350      INTEGER,  ALLOCATABLE, DIMENSION(:,:) ::   ilevel 
     351      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zthick_0, zthick 
     352 
     353      !!---------------------------------------------------------------------- 
     354 
     355      IF( .NOT. ALLOCATED(ilevel) ) THEN 
     356         ALLOCATE( ilevel(jpi,jpj), zthick_0(jpi,jpj), & 
     357         &         zthick(jpi,jpj), STAT=ji ) 
     358         IF( lk_mpp  )   CALL mpp_sum( 'zdfmxl', ji ) 
     359         IF( ji /= 0 )   CALL ctl_stop( 'STOP', 'zdf_mxl_zint_htc : unable to allocate arrays' ) 
     360      ENDIF 
     361 
     362      ! Find last whole model T level above the MLD 
     363      ilevel(:,:)   = 0 
     364      zthick_0(:,:) = 0._wp 
     365 
     366      DO jk = 1, jpkm1   
     367         DO jj = 1, jpj 
     368            DO ji = 1, jpi                     
     369               zthick_0(ji,jj) = zthick_0(ji,jj) + e3t_n(ji,jj,jk) 
     370               IF( zthick_0(ji,jj) < hmld_zint(ji,jj) )   ilevel(ji,jj) = jk 
     371            END DO 
     372         END DO 
     373         WRITE(numout,*) 'zthick_0(jk =',jk,') =',zthick_0(2,2) 
     374         WRITE(numout,*) 'gdepw_n(jk+1 =',jk+1,') =',gdepw_n(2,2,jk+1) 
     375      END DO 
     376 
     377      ! Surface boundary condition 
     378      IF( ln_linssh ) THEN  ;   zthick(:,:) = sshn(:,:)   ;   htc_mld(:,:) = tsn(:,:,1,jp_tem) * sshn(:,:) * tmask(:,:,1)    
     379      ELSE                  ;   zthick(:,:) = 0._wp       ;   htc_mld(:,:) = 0._wp                                    
     380      ENDIF 
     381 
     382      ! Deepest whole T level above the MLD 
     383      ikmax = MIN( MAXVAL( ilevel(:,:) ), jpkm1 ) 
     384 
     385      ! Integration down to last whole model T level 
     386      DO jk = 1, ikmax 
     387         DO jj = 1, jpj 
     388            DO ji = 1, jpi 
     389               zc = e3t_n(ji,jj,jk) * REAL( MIN( MAX( 0, ilevel(ji,jj) - jk + 1 ) , 1  )  )    ! 0 below ilevel 
     390               zthick(ji,jj) = zthick(ji,jj) + zc 
     391               htc_mld(ji,jj) = htc_mld(ji,jj) + zc * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) 
     392            END DO 
     393         END DO 
     394      END DO 
     395 
     396      ! Subsequent partial T level 
     397      zthick(:,:) = hmld_zint(:,:) - zthick(:,:)   !   remaining thickness to reach MLD 
     398 
     399      DO jj = 1, jpj 
     400         DO ji = 1, jpi 
     401            htc_mld(ji,jj) = htc_mld(ji,jj) + tsn(ji,jj,ilevel(ji,jj)+1,jp_tem)  &  
     402      &                      * MIN( e3t_n(ji,jj,ilevel(ji,jj)+1), zthick(ji,jj) ) * tmask(ji,jj,ilevel(ji,jj)+1) 
     403         END DO 
     404      END DO 
     405 
     406      WRITE(numout,*) 'htc_mld(after) =',htc_mld(2,2) 
     407 
     408      ! Convert to heat content 
     409      zcoef = rau0 * rcp 
     410      htc_mld(:,:) = zcoef * htc_mld(:,:) 
     411 
     412   END SUBROUTINE zdf_mxl_zint_htc 
     413 
     414   SUBROUTINE zdf_mxl_zint( kt ) 
     415      !!---------------------------------------------------------------------- 
     416      !!                  ***  ROUTINE zdf_mxl_zint  *** 
     417      !!  
     418      !! ** Purpose :    
     419      !! 
     420      !! ** Method  :    
     421      !!---------------------------------------------------------------------- 
     422 
     423      INTEGER, INTENT(in) ::   kt   ! ocean time-step index 
     424 
     425      INTEGER :: ios 
     426      INTEGER :: jn 
     427 
     428      INTEGER :: nn_mld_diag = 0    ! number of diagnostics 
     429 
     430      CHARACTER(len=1) :: cmld 
     431 
     432      TYPE(MXL_ZINT) :: sn_mld1, sn_mld2, sn_mld3, sn_mld4, sn_mld5 
     433      TYPE(MXL_ZINT), SAVE, DIMENSION(5) ::   mld_diags 
     434 
     435      NAMELIST/namzdf_mldzint/ nn_mld_diag, sn_mld1, sn_mld2, sn_mld3, sn_mld4, sn_mld5 
     436 
     437      !!---------------------------------------------------------------------- 
     438       
     439      IF( kt == nit000 ) THEN 
     440         REWIND( numnam_ref )              ! Namelist namzdf_mldzint in reference namelist  
     441         READ  ( numnam_ref, namzdf_mldzint, IOSTAT = ios, ERR = 901) 
     442901      IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_mldzint in reference namelist', lwp ) 
     443 
     444         REWIND( numnam_cfg )              ! Namelist namzdf_mldzint in configuration namelist  
     445         READ  ( numnam_cfg, namzdf_mldzint, IOSTAT = ios, ERR = 902 ) 
     446902      IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_mldzint in configuration namelist', lwp ) 
     447         IF(lwm) WRITE ( numond, namzdf_mldzint ) 
     448 
     449         IF( nn_mld_diag > 5 )   CALL ctl_stop( 'STOP', 'zdf_mxl_ini: Specify no more than 5 MLD definitions' ) 
     450 
     451         mld_diags(1) = sn_mld1 
     452         mld_diags(2) = sn_mld2 
     453         mld_diags(3) = sn_mld3 
     454         mld_diags(4) = sn_mld4 
     455         mld_diags(5) = sn_mld5 
     456 
     457         IF( lwp .AND. (nn_mld_diag > 0) ) THEN 
     458            WRITE(numout,*) '=============== Vertically-interpolated mixed layer ================' 
     459            WRITE(numout,*) '(Diagnostic number, nn_mld_type, rn_zref, rn_dT_crit, rn_iso_frac)' 
     460            DO jn = 1, nn_mld_diag 
     461               WRITE(numout,*) 'MLD criterion',jn,':' 
     462               WRITE(numout,*) '    nn_mld_type =', mld_diags(jn)%mld_type 
     463               WRITE(numout,*) '    rn_zref ='    , mld_diags(jn)%zref 
     464               WRITE(numout,*) '    rn_dT_crit =' , mld_diags(jn)%dT_crit 
     465               WRITE(numout,*) '    rn_iso_frac =', mld_diags(jn)%iso_frac 
     466            END DO 
     467            WRITE(numout,*) '====================================================================' 
     468         ENDIF 
     469      ENDIF 
     470 
     471      IF( nn_mld_diag > 0 ) THEN 
     472         DO jn = 1, nn_mld_diag 
     473            WRITE(cmld,'(I1)') jn 
     474            IF( iom_use( "mldzint_"//cmld ) .OR. iom_use( "mldhtc_"//cmld ) ) THEN 
     475               CALL zdf_mxl_zint_mld( mld_diags(jn) ) 
     476 
     477               IF( iom_use( "mldzint_"//cmld ) ) THEN 
     478                  CALL iom_put( "mldzint_"//cmld, hmld_zint(:,:) ) 
     479               ENDIF 
     480 
     481               IF( iom_use( "mldhtc_"//cmld ) )  THEN 
     482                  CALL zdf_mxl_zint_htc( kt ) 
     483                  CALL iom_put( "mldhtc_"//cmld , htc_mld(:,:)   ) 
     484               ENDIF 
     485            ENDIF 
     486         END DO 
     487      ENDIF 
     488 
     489   END SUBROUTINE zdf_mxl_zint 
    142490 
    143491   !!====================================================================== 
Note: See TracChangeset for help on using the changeset viewer.