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

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

error messages when setting a traceur.def with wrong names , or with vap/ice tracer missing for the same specie , or if the specie is not in table_tracers_condensable

File size: 8.5 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
98    subroutine specie_parameters_table(specname)
99       
100        implicit none
101    !============================================================================
102    !   Load the adequate set of parameters for specname
103    !   From a table of traceurs
104    !
105    !   Also set up few others useful parameters, such as epsi=m/mugaz, RLVTT and
106    !   the metallicity_coeff.
107    !   Authors
108    !   --------
109    !   Noé Clément (2022)
110    !============================================================================
111
112        character(*), intent(in) :: specname
113        integer :: ios
114        character(len=500):: table_traceurs_line ! table_traceurs_line lines with parameters
115
116        open(117,file='table_tracers_condensable',form='formatted',status='old')
117
118        read(117,'(A)') table_traceurs_line
119
120        do
121            read(117,'(A)', iostat=ios) table_traceurs_line
122           
123            if (index(table_traceurs_line,specname) /= 0) then
124               
125                if (index(table_traceurs_line,'delta_vapH='   ) /= 0) then
126                    read(table_traceurs_line(index(table_traceurs_line,'delta_vapH=')+len('delta_vapH='):),*) delta_vapH
127                    print*, 'delta_vapH ', delta_vapH
128                end if
129                if (index(table_traceurs_line,'Tref='   ) /= 0) then
130                    read(table_traceurs_line(index(table_traceurs_line,'Tref=')+len('Tref='):),*) Tref
131                    print*, 'Tref', Tref
132                end if
133                if (index(table_traceurs_line,'Pref='   ) /= 0) then
134                    read(table_traceurs_line(index(table_traceurs_line,'Pref=')+len('Pref='):),*) Pref
135                    print*, 'Pref', Pref
136                end if
137                if (index(table_traceurs_line,'mass='   ) /= 0) then
138                    read(table_traceurs_line(index(table_traceurs_line,'mass=')+len('mass='):),*) m
139                    print*, 'mass', m
140                end if
141                if (index(table_traceurs_line,'metallicity_coeff='   ) /= 0) then
142                   read(table_traceurs_line(index(table_traceurs_line,'metallicity_coeff=')+len('metallicity_coeff='):),*) metallicity_coeff
143                   print*, 'metallicity_coeff', metallicity_coeff
144                end if
145
146                ios=1
147            end if
148
149            if (ios /= 0) exit
150           
151        end do
152
153        if (ios .eq. -1) then ! If the file has been read but there is no data on the speci, then we have ios equal to -1
154
155            write(*,*), "ERROR : the following traceur is not referenced in the table_tracers_condensable : ", specname
156
157        end if
158       
159        close(117)
160
161        RLVTT=delta_vapH/(m/1000.)
162
163        write(*,*) 'RLVTT', RLVTT
164
165        epsi = m/mugaz
166
167    end subroutine specie_parameters_table
168
169    subroutine Psat_generic(T,p,metallicity,psat,qsat)
170        IMPLICIT NONE
171    !============================================================================
172    !   Clausius-Clapeyron relation :
173    !   d(Ln(psat))/dT = delta_vapH/(RT^2)
174    !   ->Psat = Pref * exp(-delta_vapH/R*(1/T-1/Tref))
175    !
176    !   Authors
177    !   --------
178    !   Lucas Teinturier (2022) adapted from Psat_water of Jeremy Leconte (2012)
179    !============================================================================
180        !Inputs
181        real, intent(in) :: T, p ! Temperature and pressure of the layer (in K and Pas)
182        real, intent(in) :: metallicity ! metallycity (log10)
183        !Outputs
184        real, intent(out) :: psat,qsat !saturation vapor pressure (Pa) and mass mixing ratio at saturation (kg/kg) of the layer
185
186        psat = pref*exp(-delta_vapH/(r*mugaz/1000.)*(1/T-1/Tref)-metallicity_coeff*metallicity) ! in Pa (because pref in Pa)
187
188        if (psat .gt. p) then
189            qsat = 1.
190        else
191            qsat = epsi *psat/(p-(1-epsi)*psat)
192        endif
193    end subroutine Psat_generic
194
195    subroutine Lcpdqsat_generic(T,p,psat,qsat,dqsat,dlnpsat)
196        implicit none
197
198    !===============================================================================
199    !   Compute dqsat = L/cp* d(qsat)/dT and d(Ln(psat)) = L/cp*d(Ln(psat))/dT
200    !   we have d(ln(psat))/dT =  delta_vapH/R*(1/T^2) for clausius-clapeyron
201    !   r*mugaz/1000. is the perfect gaz constant in the computation of "dummy"
202    !   Authors
203    !   --------
204    !   Lucas Teinturier (2022) adapted from Lcpdqsat_water of Jeremy Leconte (2012)
205    !===============================================================================
206        ! Inputs
207        real, intent(in) :: T ! Temperature (K)
208        real, intent(in) :: p ! Pressure (Pa)
209        real, intent(in) :: psat ! Saturation vapor pressure (Pa)
210        real, intent(in) :: qsat ! Mass mixing ratio at saturation (kg/kg)
211
212        ! Outputs
213        real, intent(out) :: dqsat,dlnpsat
214
215        ! Variables
216        real :: dummy ! used to store d(ln(psat))/dT
217
218        dummy = delta_vapH/((r*mugaz/1000.)*(T**2))
219
220        if (psat .gt. p) then
221            dqsat =0.
222        else
223            dqsat = (RLVTT/cpp) *qsat*(p/(p-(1-epsi)*psat))*dummy
224            dlnpsat = (RLVTT/cpp) * dummy
225        endif
226
227    end subroutine Lcpdqsat_generic
228
229end module generic_cloud_common_h
Note: See TracBrowser for help on using the repository browser.