Ignore:
Timestamp:
Feb 26, 2024, 4:22:34 PM (8 months ago)
Author:
gmilcareck
Message:

Add virtual correction for convective adjustment and turbulent diffusion (turbdiff and vdifc) + correction of allocated tables in physiq_mod for varspec.

Location:
trunk/LMDZ.GENERIC/libf/phystd
Files:
7 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/callkeys_mod.F90

    r3235 r3236  
    4848      logical,save :: generic_condensation
    4949      logical,save :: generic_rain
    50 !$OMP THREADPRIVATE(varactive,varfixed,sedimentation,generic_condensation,generic_rain)
     50      logical,save :: virtual_correction
     51!$OMP THREADPRIVATE(varactive,varfixed,sedimentation,generic_condensation,generic_rain,virtual_correction)
    5152      logical,save :: water ,watercond, waterrain, moistadjustment
    5253!$OMP THREADPRIVATE(water, watercond, waterrain, moistadjustment)
  • trunk/LMDZ.GENERIC/libf/phystd/convadj.F

    r2482 r3236  
    77      USE tracer_h
    88      use comcstfi_mod, only: g
    9       use callkeys_mod, only: tracer,water
     9      use generic_cloud_common_h
     10      use generic_tracer_index_mod, only: generic_tracer_index
     11      use callkeys_mod, only: tracer,water,generic_condensation,
     12     &                        virtual_correction
    1013
    1114      implicit none
     
    5558      REAL sig(nlay+1),sdsig(nlay),dsig(nlay)
    5659      REAL zu(ngrid,nlay),zv(ngrid,nlay)
    57       REAL zh(ngrid,nlay)
     60      REAL zh(ngrid,nlay), zvh(ngrid,nlay)
    5861      REAL zu2(ngrid,nlay),zv2(ngrid,nlay)
    59       REAL zh2(ngrid,nlay), zhc(ngrid,nlay)
     62      REAL zh2(ngrid,nlay),zvh2(ngrid,nlay),zhc(ngrid,nlay)
    6063      REAL zhm,zsm,zdsm,zum,zvm,zalpha,zhmc
    6164
     
    6467      REAL zq(ngrid,nlay,nq), zq2(ngrid,nlay,nq)
    6568      REAL zqm(nq),zqco2m
     69     
     70      integer igcm_generic_vap, igcm_generic_ice
     71      logical call_ice_vap_generic
    6672
    6773      LOGICAL vtest(ngrid),down
     
    7581!     Initialisation
    7682!     --------------
    77 
    78       DO l=1,nlay
    79          DO ig=1,ngrid
    80             zh(ig,l)=ph(ig,l)+pdhfi(ig,l)*ptimestep
    81             zu(ig,l)=pu(ig,l)+pdufi(ig,l)*ptimestep
    82             zv(ig,l)=pv(ig,l)+pdvfi(ig,l)*ptimestep
    83          ENDDO
    84       ENDDO
     83      zh(:,:)=ph(:,:)+pdhfi(:,:)*ptimestep
     84      zu(:,:)=pu(:,:)+pdufi(:,:)*ptimestep
     85      zv(:,:)=pv(:,:)+pdvfi(:,:)*ptimestep
    8586
    8687      if(tracer) then     
    87         DO iq =1, nq
    88          DO l=1,nlay
    89            DO ig=1,ngrid
    90               zq(ig,l,iq)=pq(ig,l,iq)+pdqfi(ig,l,iq)*ptimestep
    91            ENDDO
    92          ENDDO
    93         ENDDO
     88        zq(:,:,:)=pq(:,:,:)+pdqfi(:,:,:)*ptimestep
     89        !zq(:,:,:)=0
    9490      end if
    9591
     
    9995      CALL scopy(ngrid*nlay*nq,zq,1,zq2,1)
    10096
     97
    10198!     -----------------------------
    10299!     Detection of unstable columns
     
    108105      ENDDO
    109106
    110       CALL scopy(ngrid*nlay,zh2,1,zhc,1)
     107      if((generic_condensation) .and. (virtual_correction)) THEN
     108        DO iq=1,nq
     109          call generic_tracer_index(nq,iq,igcm_generic_vap,
     110     &         igcm_generic_ice,call_ice_vap_generic)
     111          if(call_ice_vap_generic) then
     112            zvh(:,:)=zh(:,:)*
     113     &         (1.+zq(:,:,igcm_generic_vap)/epsi_generic)/
     114     &         (1.+zq(:,:,igcm_generic_vap))
     115          endif
     116        ENDDO
     117        CALL scopy(ngrid*nlay,zvh,1,zvh2,1)
     118        CALL scopy(ngrid*nlay,zvh2,1,zhc,1)
     119      else       
     120        CALL scopy(ngrid*nlay,zh2,1,zhc,1)
     121      endif
    111122
    112123!     Find out which grid points are convectively unstable
     
    169180                zdsm = zdsm + dsig(l)
    170181                zhm = zhm + sdsig(l) * (zh2(i, l) - zhm) / zsm
    171                 zhmc = zhm
     182 
     183                if (generic_condensation .and. virtual_correction) then
     184                  zhmc = zhm*
     185     &            (1.+zq(i,l,igcm_generic_vap)/epsi_generic)/
     186     &            (1.+zq(i,l,igcm_generic_vap))
     187                else
     188                  zhmc = zhm
     189                endif
    172190 
    173191!     do we have to extend the column downwards?
     
    327345      ENDDO
    328346
    329       DO l=1,nlay
    330         DO ig=1,ngrid
    331           pdhadj(ig,l)=(zh2(ig,l)-zh(ig,l))/ptimestep
    332           pduadj(ig,l)=(zu2(ig,l)-zu(ig,l))/ptimestep
    333           pdvadj(ig,l)=(zv2(ig,l)-zv(ig,l))/ptimestep
    334         ENDDO
    335       ENDDO
     347      pdhadj(:,:)=(zh2(:,:)-zh(:,:))/ptimestep
     348      pduadj(:,:)=(zu2(:,:)-zu(:,:))/ptimestep
     349      pdvadj(:,:)=(zv2(:,:)-zv(:,:))/ptimestep
    336350
    337351      if(tracer) then
    338         do iq=1, nq
    339           do  l=1,nlay
    340             DO ig=1,ngrid
    341               pdqadj(ig,l,iq)=(zq2(ig,l,iq)-zq(ig,l,iq))/ptimestep
    342             end do
    343           end do
    344         end do
     352        pdqadj(:,:,:)=(zq2(:,:,:)-zq(:,:,:))/ptimestep
    345353      end if
    346354
  • trunk/LMDZ.GENERIC/libf/phystd/inifis_mod.F90

    r3235 r3236  
    10241024     call getin_p("generic_rain",generic_rain)
    10251025     if (is_master) write(*,*)trim(rname)//": generic_rain = ",generic_rain
     1026     
     1027     if (is_master) write(*,*)trim(rname)//": Virtual correction ?"
     1028     virtual_correction=.false. !default value
     1029     call getin_p("virtual_correction",virtual_correction)
     1030     if (is_master) write(*,*)trim(rname)//": virtual_correction = ",virtual_correction
    10261031
    10271032     if (is_master) write(*,*)trim(rname)//": Compute water cycle ?"
  • trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90

    r3235 r3236  
    551551           endif
    552552           close(33)
     553         else
     554           IF (.NOT.ALLOCATED(p_var))  ALLOCATE(p_var(nlayer))
     555           IF (.NOT.ALLOCATED(mu_var))  ALLOCATE(mu_var(nlayer))
     556           IF (.NOT.ALLOCATED(frac_var))  ALLOCATE(frac_var(nlayer,ngasmx))
    553557         endif
    554558
  • trunk/LMDZ.GENERIC/libf/phystd/turbdiff_mod.F90

    r2482 r3236  
    247247!     ------------------------------------------------------
    248248
    249       call vdif_kc(ngrid,nlay,ptimestep,g,pzlev,pzlay,pu,pv,zh,zcdv_true,pq2,zkv,zkh) !JL12 why not call vdif_kc with updated winds and temperature
     249      call vdif_kc(ngrid,nlay,nq,ptimestep,g,pzlev,pzlay,pu,pv,pq,zh,zcdv_true,pq2,zkv,zkh) !JL12 why not call vdif_kc with updated winds and temperature
    250250     
    251251!     Adding eddy mixing to mimic 3D general circulation in 1D
  • trunk/LMDZ.GENERIC/libf/phystd/vdif_kc.F

    r1529 r3236  
    1       SUBROUTINE vdif_kc(ngrid,nlay,dt,g,zlev,zlay,u,v,teta,cd,q2,km,kn)
     1      SUBROUTINE vdif_kc(ngrid,nlay,nq,dt,g,zlev,zlay,u,v,
     2     &                       zq,teta,cd,q2,km,kn)
     3      use generic_cloud_common_h, only: epsi_generic
     4      use generic_tracer_index_mod, only: generic_tracer_index
     5      use callkeys_mod, only: generic_condensation, virtual_correction
    26      IMPLICIT NONE
    37c.......................................................................
     
    2630      INTEGER,INTENT(IN) :: ngrid
    2731      INTEGER,INTENT(IN) :: nlay
     32      INTEGER,INTENT(IN) :: nq
    2833      REAL,INTENT(IN) :: dt,g
    2934      REAL,INTENT(IN) :: zlev(ngrid,nlay+1)
     
    3136      REAL,INTENT(IN) :: u(ngrid,nlay)
    3237      REAL,INTENT(IN) :: v(ngrid,nlay)
     38      REAL,INTENT(IN) :: zq(ngrid,nlay,nq)
    3339      REAL,INTENT(IN) :: teta(ngrid,nlay)
    3440      REAL,INTENT(IN) :: cd(ngrid)
     
    5258      REAL unsdzdec(ngrid,nlay+1)
    5359      REAL q(ngrid,nlay+1)
     60      REAL tetav(ngrid,nlay)
     61      integer iq,igcm_generic_vap, igcm_generic_ice
     62      logical call_ice_vap_generic
    5463c.......................................................................
    5564c
     
    250259        mpre(igrid,1)=m(igrid,1)
    251260      ENDDO
     261
     262c
     263c.......................................................................
     264c  Virtual correction
     265c.......................................................................
     266c
     267
     268      if((generic_condensation) .and. (virtual_correction)) THEN
     269        DO ilev=1,nlay
     270         DO igrid=1,ngrid
     271          DO iq=1,nq
     272            call generic_tracer_index(nq,iq,igcm_generic_vap,
     273     &         igcm_generic_ice,call_ice_vap_generic)
     274            if(call_ice_vap_generic) then
     275               tetav(igrid,ilev) = teta(igrid,ilev)*
     276     &              (1.+zq(igrid,ilev,igcm_generic_vap)/epsi_generic)/
     277     &              (1.+zq(igrid,ilev,igcm_generic_vap))
     278            endif
     279          ENDDO
     280         ENDDO
     281        ENDDO
     282      endif
     283
    252284c
    253285c-----------------------------------------------------------------------
     
    256288c-----------------------------------------------------------------------
    257289c
     290      if((generic_condensation) .and. (virtual_correction)) THEN
    258291        n2(igrid,ilev)=g*unsdzdec(igrid,ilev)
    259      &                   *(teta(igrid,ilev)-teta(igrid,ilev-1))
    260      &                   /(teta(igrid,ilev)+teta(igrid,ilev-1)) *2.E+0
     292     &                 *(tetav(igrid,ilev)-tetav(igrid,ilev-1))
     293     &                 /(tetav(igrid,ilev)+tetav(igrid,ilev-1))*2.E+0
     294      else
     295        n2(igrid,ilev)=g*unsdzdec(igrid,ilev)
     296     &                 *(teta(igrid,ilev)-teta(igrid,ilev-1))
     297     &                 /(teta(igrid,ilev)+teta(igrid,ilev-1)) *2.E+0
     298      endif
    261299c
    262300c --->
  • trunk/LMDZ.GENERIC/libf/phystd/vdifc_mod.F

    r2427 r3236  
    246246!     ------------------------------------------------------
    247247
    248       call vdif_kc(ngrid,nlay,ptimestep,g,pzlev,pzlay
    249      &     ,pu,pv,ph,zcdv_true
     248      call vdif_kc(ngrid,nlay,nq,ptimestep,g,pzlev,pzlay
     249     &     ,pu,pv,pq,ph,zcdv_true
    250250     &     ,pq2,zkv,zkh)
    251251
Note: See TracChangeset for help on using the changeset viewer.