Ignore:
Timestamp:
Aug 7, 2013, 4:21:01 PM (12 years ago)
Author:
jleconte
Message:

07/08/2013 == JL

  • Water cycle in double precision (largescale+moistadj)
  • Improved wate rayleigh.
  • First step for rayleigh with variable species. Now, just need to change optcv.
  • changed some interpolation indices in callcorrk to limit dependency of OLR on the number of layers
File:
1 edited

Legend:

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

    r989 r1016  
    175175      end subroutine Tsat_water
    176176
     177!==================================================================
     178      subroutine Psat_waterDP(T,p,psat,qsat)
     179
     180         implicit none
     181
     182!==================================================================
     183!     Purpose
     184!     -------
     185!     Compute the saturation vapor pressure and mass mixing ratio at saturation (kg/kg)
     186!     for a given pressure (Pa) and temperature (K)
     187!     DOUBLE PRECISION
     188!     Based on the Tetens formula from L.Li physical parametrization manual
     189!
     190!     Authors
     191!     -------
     192!     Jeremy Leconte (2012)
     193!
     194!==================================================================
     195
     196!        input
     197         double precision, intent(in) :: T, p
     198 
     199!        output
     200         double precision psat,qsat
     201
     202! JL12 variables for tetens formula
     203         double precision,parameter :: Pref_solid_liquid=611.14d0
     204         double precision,parameter :: Trefvaporization=35.86D0
     205         double precision,parameter :: Trefsublimation=7.66d0
     206         double precision,parameter :: Tmin=8.d0
     207         double precision,parameter :: r3vaporization=17.269d0
     208         double precision,parameter :: r3sublimation=21.875d0
     209
     210! checked vs. old watersat data 14/05/2012 by JL.
     211
     212         if (T.gt.T_h2O_ice_liq) then
     213            psat = Pref_solid_liquid*Exp(r3vaporization*(T-T_h2O_ice_liq)/(T-Trefvaporization)) ! liquid / vapour
     214         else if (T.lt.Tmin) then
     215            print*, "careful, T<Tmin in psat water"
     216            psat = Pref_solid_liquid*Exp(r3sublimation*(Tmin-T_h2O_ice_liq)/(Tmin-Trefsublimation)) ! min psat 
     217         else                 
     218            psat = Pref_solid_liquid*Exp(r3sublimation*(T-T_h2O_ice_liq)/(T-Trefsublimation)) ! solid / vapour
     219         endif
     220         if(psat.gt.p) then
     221            qsat=1.d0
     222         else
     223            qsat=epsi*psat/(p-(1.d0-epsi)*psat)
     224         endif
     225         return
     226      end subroutine Psat_waterDP
     227
     228
     229
     230
     231!==================================================================
     232      subroutine Lcpdqsat_waterDP(T,p,psat,qsat,dqsat,dlnpsat)
     233
     234         implicit none
     235
     236!==================================================================
     237!     Purpose
     238!     -------
     239!     Compute dqsat=L/cp*d (q_sat)/d T and dlnpsat=L/cp d(ln Psat)/d T
     240!     for a given temperature (K)!
     241!     Based on the Tetens formula from L.Li physical parametrization manual
     242!
     243!     Authors
     244!     -------
     245!     Jeremy Leconte (2012)
     246!
     247!==================================================================
     248
     249!        input
     250         double precision T, p, psat, qsat
     251 
     252!        output
     253         double precision dqsat,dlnpsat
     254
     255! JL12 variables for tetens formula
     256         double precision,parameter :: Pref_solid_liquid=611.14d0
     257         double precision,parameter :: Trefvaporization=35.86d0
     258         double precision,parameter :: Tmin=8.d0
     259         double precision,parameter :: Trefsublimation=7.66d0
     260         double precision,parameter :: r3vaporization=17.269d0
     261         double precision,parameter :: r3sublimation=21.875d0
     262
     263         double precision :: dummy
     264
     265         if (psat.gt.p) then
     266            dqsat=0.d0
     267            return
     268         endif
     269
     270         if (T.gt.T_h2O_ice_liq) then
     271            dummy = r3vaporization*(T_h2O_ice_liq-Trefvaporization)/(T-Trefvaporization)**2  ! liquid / vapour
     272         else if (T.lt.Tmin) then
     273            print*, "careful, T<Tmin in Lcp psat water"
     274            dummy = r3sublimation*(T_h2O_ice_liq-Trefsublimation)/(Tmin-Trefsublimation)**2  ! solid / vapour
     275         else               
     276            dummy = r3sublimation*(T_h2O_ice_liq-Trefsublimation)/(T-Trefsublimation)**2  ! solid / vapour
     277         endif
     278
     279         dqsat=RLVTT/RCPD*qsat*(p/(p-(1.d0-epsi)*psat))*dummy
     280         dlnpsat=RLVTT/RCPD*dummy
     281         return
     282      end subroutine Lcpdqsat_waterDP
     283
    177284
    178285end module watercommon_h
Note: See TracChangeset for help on using the changeset viewer.