source: trunk/LMDZ.GENERIC/libf/phystd/generic_cloud_common_h.F90 @ 2705

Last change on this file since 2705 was 2705, checked in by aslmd, 2 years ago

Fixed the molecular mass of composites generic clouds

File size: 5.8 KB
Line 
1module generic_cloud_common_h
2    use comcstfi_mod, only: r, cpp, mugaz
3    implicit none
4    real, save :: m  ! molecular mass of the specie (g/mol)
5    real, save :: delta_vapH ! Enthalpy of vaporization (J/mol)
6    real,save :: Tref ! Ref temperature for Clausis-Clapeyron (K)
7    real,save :: Pref ! Reference pressure for Clausius Clapeyron (Pa)
8    real, save :: epsi ! fractionnal molecular mass (m/mugaz)
9    real,save :: RLVTT !Latent heat of vaporization (J/kg)
10    real,save :: metallicity_coeff ! Coefficient to take into account the metallicity
11    contains
12
13    subroutine specie_parameters(specname)
14       
15        implicit none
16
17    !============================================================================
18    !   Load the adequate set of parameters for specname
19    !   We use Clausius Clapeyron.
20    !   Values are based on Vissher (2010) and Morley (2012).
21    !
22    !   Also set up few others useful parameters, such as epsi=m/mugaz, RLVTT and
23    !   the metallicity_coeff.
24    !   Authors
25    !   --------
26    !   Lucas Teinturier (2022)
27    !============================================================================
28        ! Inputs
29        character(*), intent(in) :: specname
30        print*,"entree dans specie_parameters avec specname = ",specname
31        if (trim(specname) .eq. "Mg") then
32            print*,"Loading data for Mg"
33            delta_vapH = 521.7E3 !J/mol
34            Tref = 2303.0 !K
35            Pref =   1.0E5 !in Pa
36            m = 140.6931
37            RLVTT = delta_vapH/(m/1000.)
38            metallicity_coeff=1.0*log(10.0)
39           
40        else if (trim(specname) .eq. "Na") then
41            print*,"Loading data for Na"
42            delta_vapH = 265.9E3
43            Tref = 1624.0
44            Pref =  1.0E5
45            m = 78.04454 !Na2S
46            RLVTT = delta_vapH/(m/1000.)
47            metallicity_coeff=0.5*log(10.0)
48           
49        else if (trim(specname) .eq. "Fe") then
50            print*,"Loading data for Fe"
51            delta_vapH = 401.9E3
52            Tref = 2903.0
53            Pref  = 1.0E5
54            m = 55.8450
55            RLVTT = delta_vapH/(m/1000.)
56            metallicity_coeff=0.0*log(10.0)
57           
58        else if (trim(specname) .eq. "Cr") then
59            print*,"Loading data for Cr"
60            delta_vapH = 394.2E3
61            Tref = 2749.0
62            Pref  = 1.0E5
63            m = 51.9961
64            RLVTT = delta_vapH/(m/1000.)
65            metallicity_coeff=0.0*log(10.0)
66           
67        else if (trim(specname) .eq. "KCl") then
68            print*,"Loading data for KCl"
69            delta_vapH = 217.9E3
70            Tref = 1495.0
71            Pref = 1.0E5
72            metallicity_coeff=0.0*log(10.0)
73            m = 74.5498
74            RLVTT = delta_vapH/(m/1000.)
75        else if (trim(specname) .eq. "Mn") then
76            print*,"Loading data for Mn"
77            delta_vapH = 455.8E3
78            Tref = 2064.0
79            Pref = 1.0E5
80            metallicity_coeff=1.0*log(10.0)
81            m = 87.003049
82            RLVTT = delta_vapH/(m/1000.)
83        else if (trim(specname) .eq. "Zn") then
84            print*,"Loading data for Zn"
85            delta_vapH = 303.9E3
86            Tref = 1238.0
87            Pref = 1.0E5
88            metallicity_coeff=1.0*log(10.0)
89            m = 97.445
90            RLVTT = delta_vapH/(m/1000.)
91        else
92            print*,"Unknow species (not in Mg, Fe, Na, KCl, Cr, Mn or Zn)"
93        endif
94        epsi = m/mugaz
95    end subroutine specie_parameters
96
97    subroutine Psat_generic(T,p,metallicity,psat,qsat)
98        IMPLICIT NONE
99    !============================================================================
100    !   Clausius-Clapeyron relation :
101    !   d(Ln(psat))/dT = delta_vapH/(RT^2)
102    !   ->Psat = Pref * exp(-delta_vapH/R*(1/T-1/Tref))
103    !
104    !   Authors
105    !   --------
106    !   Lucas Teinturier (2022) adapted from Psat_water of Jeremy Leconte (2012)
107    !============================================================================
108        !Inputs
109        real, intent(in) :: T, p ! Temperature and pressure of the layer (in K and Pas)
110        real, intent(in) :: metallicity ! metallycity (log10)
111        !Outputs
112        real, intent(out) :: psat,qsat !saturation vapor pressure (Pa) and mass mixing ratio at saturation (kg/kg) of the layer
113
114        psat = pref*exp(-delta_vapH/(r*mugaz/1000.)*(1/T-1/Tref)-metallicity_coeff*metallicity) ! in Pa (because pref in Pa)
115
116        if (psat .gt. p) then
117            qsat = 1.
118        else
119            qsat = epsi *psat/(p-(1-epsi)*psat)
120        endif
121    end subroutine Psat_generic
122
123    subroutine Lcpdqsat_generic(T,p,psat,qsat,dqsat,dlnpsat)
124        implicit none
125
126    !===============================================================================
127    !   Compute dqsat = L/cp* d(qsat)/dT and d(Ln(psat)) = L/cp*d(Ln(psat))/dT
128    !   we have d(ln(psat))/dT =  delta_vapH/R*(1/T^2) for clausius-clapeyron
129    !   r*mugaz/1000. is the perfect gaz constant in the computation of "dummy"
130    !   Authors
131    !   --------
132    !   Lucas Teinturier (2022) adapted from Lcpdqsat_water of Jeremy Leconte (2012)
133    !===============================================================================
134        ! Inputs
135        real, intent(in) :: T ! Temperature (K)
136        real, intent(in) :: p ! Pressure (Pa)
137        real, intent(in) :: psat ! Saturation vapor pressure (Pa)
138        real, intent(in) :: qsat ! Mass mixing ratio at saturation (kg/kg)
139
140        ! Outputs
141        real, intent(out) :: dqsat,dlnpsat
142
143        ! Variables
144        real :: dummy ! used to store d(ln(psat))/dT
145
146        dummy = delta_vapH/((r*mugaz/1000.)*(T**2))
147
148        if (psat .gt. p) then
149            dqsat =0.
150        else
151            dqsat = (RLVTT/cpp) *qsat*(p/(p-(1-epsi)*psat))*dummy
152            dlnpsat = (RLVTT/cpp) * dummy
153        endif
154
155    end subroutine Lcpdqsat_generic
156
157end module generic_cloud_common_h
Note: See TracBrowser for help on using the repository browser.