Changeset 3219 for trunk/LMDZ.MARS


Ignore:
Timestamp:
Feb 16, 2024, 10:23:53 AM (9 months ago)
Author:
llange
Message:

Mars PCM
Small correction for the surface layer scheme: the value of the critical Richardson varys if one used the classic yamada scheme (0.2 in this case) or can be a tunning parameter if one uses the atke scheme
Also add the possibily to write the tke in the diagfi.nc when calling pbl_parameter
LL

Location:
trunk/LMDZ.MARS
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/changelog.txt

    r3217 r3219  
    44964496== 13/02/2024 == EM
    44974497Move "Martian nogcm" from common dynamics to Mars dynamics-physics interface.
     4498
     4499== 16/02/2024 == LL
     4500Small correction for the surface layer scheme: the value of the critical Richardson varys if one used the classic yamada scheme (0.2 in this case) or can be a tunning parameter if one uses the atke scheme
     4501Also add the possibily to write the tke in the diagfi.nc when calling pbl_parameter
  • trunk/LMDZ.MARS/libf/phymars/pbl_parameters_mod.F90

    r3167 r3219  
    1919CONTAINS
    2020
    21       SUBROUTINE pbl_parameters(ngrid,nlay,ps,pplay,pz0,pg,zzlay,zzlev,pu,pv,wstar_in,hfmax,zmax,pts,ph,z_out,n_out, &
     21      SUBROUTINE pbl_parameters(ngrid,nlay,ps,pplay,pz0,pg,zzlay,zzlev,pu,pv,wstar_in,hfmax,zmax,tke,pts,ph,z_out,n_out, &
    2222                                T_out,u_out,ustar,tstar,vhf,vvv)
    2323                               
     
    2525      use write_output_mod, only: write_output
    2626      use turb_mod, only: turb_resolved
    27       use lmdz_atke_turbulence_ini, only : smmin, cinf, cepsilon, pr_slope, pr_asym, pr_neut, ri0, ri1, cn, rpi
     27      use lmdz_atke_turbulence_ini, only : smmin, ric, cinf, cepsilon, pr_slope, pr_asym, pr_neut, ri0, ri1, cn, rpi
    2828
    2929      IMPLICIT NONE
     
    5757!     hfmax(ngrid)        maximum vertical turbulent heat flux in thermals (W/m^2)
    5858!     zmax(ngrid)         height reached by the thermals (pbl height) (m)
     59!     tke(ngrid,nlay+1)   turbulent kinetic energy (J/kg)
    5960!     pts(ngrid)          surface temperature (K)
    6061!     ph(ngrid,nlay)      potential temperature T*(p/ps)^kappa (k)
     
    8889      REAL, INTENT(IN) :: pu(ngrid,nlay),pv(ngrid,nlay)             ! winds
    8990      REAL, INTENT(IN) :: wstar_in(ngrid),hfmax(ngrid),zmax(ngrid)  ! free convection velocity , vertical turbulent heat, pbl height
     91      REAL, INTENT(IN) :: tke(ngrid,nlay+1)                         ! Turbulent kinetic energy
    9092      REAL, INTENT(IN) :: pts(ngrid),ph(ngrid,nlay)                 ! surface temperature, potentiel temperature
    9193      REAL, INTENT(IN) :: z_out(n_out)                              ! altitudes of the interpolation in the surface layer
     
    123125      REAL cdn(ngrid),chn(ngrid)                                    ! neutral momentum and heat drag coefficient (1)
    124126      REAL pz0t                                                     ! initial thermal roughness length z0t for the iterative scheme (m)
    125       REAL ric                                                      ! critical richardson number (1)
     127      REAL ric_colaitis                                             ! critical richardson number (1)
    126128      REAL reynolds(ngrid)                                          ! Reynolds number (1)
    127129      REAL prandtl(ngrid)                                           ! Prandtl number  (1)
     
    139141      REAL sm                                                       ! Stability function computed with the ATKE scheme
    140142      REAL f_ri_cd_min                                              ! minimum of the stability function with the ATKE scheme
    141 
     143      REAL ric_4interp                                              ! critical richardson number which is either the one from Colaitis et al. (2013) or ATKE (1)
    142144!    Code:
    143145!    --------
     
    173175      betah = 5.         
    174176      lambda = (sqrt(bh/bm))/alphah
    175      
    176       ric = betah/(betam**2)
     177      ric_colaitis = betah/(betam**2)
     178     
     179      if(callatke) then
     180         ric_4interp = ric
     181      else
     182         ric_4interp = ric_colaitis
     183      endif
     184     
    177185     
    178186      DO ig=1,ngrid
     
    195203               print*,'warning, infinite Richardson at surface'
    196204               print*,pu(ig,1),pv(ig,1)
    197                rib(ig) = ric
     205               rib(ig) = ric_colaitis
    198206            ENDIF
    199207
     
    221229                   ! STABLE BOUNDARY LAYER :
    222230                   prandtl(ig) = 1.
    223                    IF(rib(ig) .lt. ric) THEN
     231                   IF(rib(ig) .lt. ric_colaitis) THEN
    224232                      ! Assuming alphah=1. and bh=bm for stable conditions :
    225                       fm(ig) = ((ric-rib(ig))/ric)**2
    226                       fh(ig) = ((ric-rib(ig))/ric)**2
     233                      fm(ig) = ((ric_colaitis-rib(ig))/ric_colaitis)**2
     234                      fh(ig) = ((ric_colaitis-rib(ig))/ric_colaitis)**2
    227235                   ELSE
    228236                      ! For Ri>Ric, we consider Ri->Infinity => no turbulent mixing at surface
     
    259267         zout = z_out(n)
    260268         DO ig=1,ngrid
    261             IF (rib(ig) .ge. ric) THEN !stable case, no turbulence
     269            IF (rib(ig) .ge. ric_4interp) THEN !stable case, no turbulence
    262270               ustar(ig) = 0.
    263271               tstar(ig) = 0.
     
    275283                u_out(ig,n) = 0.
    276284            ELSE
    277                IF (rib(ig) .ge. ric) THEN ! ustar=tstar=0  (and fm=fh=0) because no turbulence
     285               IF (rib(ig) .ge. ric_4interp) THEN ! ustar=tstar=0  (and fm=fh=0) because no turbulence
    278286                  u_out(ig,n) = 0
    279287                  Teta_out(ig,n) = pts(ig)
     
    346354
    347355#ifndef MESOSCALE
     356            call write_output('tke_pbl','Tke in the pbl after physic','J/kg',tke(:,:))
    348357            call write_output('rib_pbl','Richardson in pbl parameter','m/s',rib(:))
    349358            call write_output('cdn_pbl','neutral momentum coef','m/s',cdn(:))
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r3207 r3219  
    25812581
    25822582           call pbl_parameters(ngrid,nlayer,ps,zplay,z0,
    2583      & g,zzlay,zzlev,zu,zv,wstar,hfmax_th,zmax_th,tsurf(:,iflat),
     2583     & g,zzlay,zzlev,zu,zv,wstar,hfmax_th,zmax_th,q2,tsurf(:,iflat),
    25842584     & zh,z_out,n_out,T_out,u_out,ustar,tstar,vhf,vvv)
    25852585                   ! pourquoi ustar recalcule ici? fait dans vdifc.
  • trunk/LMDZ.MARS/libf/phymars/vdif_cd_mod.F90

    r3167 r3219  
    1  MODULE vdif_cd_mod
     1MODULE vdif_cd_mod
    22 
    33!======================================================================================================================!
     
    2121   use turb_mod, only: turb_resolved
    2222   use watersat_mod, only: watersat
    23    use lmdz_atke_turbulence_ini, only : smmin, cinf, cepsilon, pr_slope, pr_asym, pr_neut, ri0, ri1, cn, rpi
     23   use lmdz_atke_turbulence_ini, only : smmin, ric, cinf, cepsilon, pr_slope, pr_asym, pr_neut, ri0, ri1, cn, rpi
    2424   
    2525   IMPLICIT NONE
     
    112112
    113113      REAL pz0t                  ! initial thermal roughness length z0t for the iterative scheme (m)
    114       REAL ric                   ! critical richardson number (1)
     114      REAL ric_colaitis          ! critical richardson number (1)
    115115      REAL reynolds(ngrid)       ! Reynolds number (1)
    116116      REAL prandtl(ngrid)        ! Prandtl number  (1)
     
    160160      lambda = (sqrt(bh/bm))/alphah
    161161     
    162       ric = betah/(betam**2)
     162      ric_colaitis = betah/(betam**2)
    163163     
    164164      IF(virtual) then
     
    214214                      print*,'warning, infinite Richardson at surface'
    215215                      print*,pu(ig,1),pv(ig,1)
    216                       rib(ig) = ric         
     216                      rib(ig) = ric_colaitis         
    217217                   ENDIF
    218218       
     
    240240                        ! STABLE BOUNDARY LAYER :
    241241                        prandtl(ig) = 1.
    242                         IF(rib(ig) .lt. ric) THEN
     242                        IF(rib(ig) .lt. ric_colaitis) THEN
    243243               ! Assuming alphah=1. and bh=bm for stable conditions :
    244                            fm(ig) = ((ric-rib(ig))/ric)**2
    245                            fh(ig) = ((ric-rib(ig))/ric)**2
     244                           fm(ig) = ((ric_colaitis-rib(ig))/ric_colaitis)**2
     245                           fh(ig) = ((ric_colaitis-rib(ig))/ric_colaitis)**2
    246246                        ELSE
    247247               ! For Ri>Ric, we consider Ri->Infinity => no turbulent mixing at surface
     
    284284     
    285285     
     286     
     287     
    286288       
    287289
Note: See TracChangeset for help on using the changeset viewer.