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

Last change on this file since 2921 was 2786, checked in by emillour, 3 years ago

Generic GCM:
Very minor fix: remove unecessary test about datadir introduced in r2715.
The real issue was that datadir was not set at this point if initracer
was called by newstart. Should now be fine after newstart update of r2785.
EM

File size: 10.0 KB
RevLine 
[2701]1module generic_cloud_common_h
[2714]2
3    !============================================================================
4    ! Purpose
5    ! -------
6    ! Set up relevant constants and parameters for generic condensable tracers cycles,
7    ! and generic condensable tracers cloud properties
8    !
9    ! Authors
10    ! -------
11    ! Adapted from watercommon_h.F90
12    ! by Lucas Teinturier & Noé Clément (2022)
13    !============================================================================
14
[2690]15    use comcstfi_mod, only: r, cpp, mugaz
16    implicit none
17    real, save :: m  ! molecular mass of the specie (g/mol)
18    real, save :: delta_vapH ! Enthalpy of vaporization (J/mol)
19    real,save :: Tref ! Ref temperature for Clausis-Clapeyron (K)
20    real,save :: Pref ! Reference pressure for Clausius Clapeyron (Pa)
[2725]21    real, save :: epsi_generic ! fractionnal molecular mass (m/mugaz)
22    real,save :: RLVTT_generic !Latent heat of vaporization (J/kg)
[2690]23    real,save :: metallicity_coeff ! Coefficient to take into account the metallicity
[2718]24
25   
[2690]26    contains
27
28    subroutine specie_parameters(specname)
29       
30        implicit none
31
32    !============================================================================
33    !   Load the adequate set of parameters for specname
34    !   We use Clausius Clapeyron.
35    !   Values are based on Vissher (2010) and Morley (2012).
36    !
[2725]37    !   Also set up few others useful parameters, such as epsi_generic=m/mugaz, RLVTT_generic and
[2690]38    !   the metallicity_coeff.
39    !   Authors
40    !   --------
41    !   Lucas Teinturier (2022)
42    !============================================================================
43        ! Inputs
44        character(*), intent(in) :: specname
45        print*,"entree dans specie_parameters avec specname = ",specname
46        if (trim(specname) .eq. "Mg") then
47            print*,"Loading data for Mg"
48            delta_vapH = 521.7E3 !J/mol
49            Tref = 2303.0 !K
50            Pref =   1.0E5 !in Pa
[2697]51            m = 140.6931
[2725]52            RLVTT_generic = delta_vapH/(m/1000.)
[2690]53            metallicity_coeff=1.0*log(10.0)
[2694]54           
[2690]55        else if (trim(specname) .eq. "Na") then
56            print*,"Loading data for Na"
57            delta_vapH = 265.9E3
58            Tref = 1624.0
59            Pref =  1.0E5
[2705]60            m = 78.04454 !Na2S
[2725]61            RLVTT_generic = delta_vapH/(m/1000.)
[2690]62            metallicity_coeff=0.5*log(10.0)
[2694]63           
[2690]64        else if (trim(specname) .eq. "Fe") then
65            print*,"Loading data for Fe"
66            delta_vapH = 401.9E3
67            Tref = 2903.0
68            Pref  = 1.0E5
[2694]69            m = 55.8450
[2725]70            RLVTT_generic = delta_vapH/(m/1000.)
[2690]71            metallicity_coeff=0.0*log(10.0)
[2694]72           
[2690]73        else if (trim(specname) .eq. "Cr") then
74            print*,"Loading data for Cr"
[2702]75            delta_vapH = 394.2E3
[2690]76            Tref = 2749.0
77            Pref  = 1.0E5
[2694]78            m = 51.9961
[2725]79            RLVTT_generic = delta_vapH/(m/1000.)
[2690]80            metallicity_coeff=0.0*log(10.0)
[2694]81           
[2690]82        else if (trim(specname) .eq. "KCl") then
83            print*,"Loading data for KCl"
84            delta_vapH = 217.9E3
85            Tref = 1495.0
86            Pref = 1.0E5
87            metallicity_coeff=0.0*log(10.0)
88            m = 74.5498
[2725]89            RLVTT_generic = delta_vapH/(m/1000.)
[2690]90        else if (trim(specname) .eq. "Mn") then
91            print*,"Loading data for Mn"
92            delta_vapH = 455.8E3
93            Tref = 2064.0
94            Pref = 1.0E5
95            metallicity_coeff=1.0*log(10.0)
[2705]96            m = 87.003049
[2725]97            RLVTT_generic = delta_vapH/(m/1000.)
[2690]98        else if (trim(specname) .eq. "Zn") then
99            print*,"Loading data for Zn"
100            delta_vapH = 303.9E3
101            Tref = 1238.0
102            Pref = 1.0E5
103            metallicity_coeff=1.0*log(10.0)
[2705]104            m = 97.445
[2725]105            RLVTT_generic = delta_vapH/(m/1000.)
[2690]106        else
107            print*,"Unknow species (not in Mg, Fe, Na, KCl, Cr, Mn or Zn)"
108        endif
[2725]109        epsi_generic = m/mugaz
[2690]110    end subroutine specie_parameters
111
[2706]112
113    subroutine specie_parameters_table(specname)
[2710]114        use datafile_mod, only: datadir
[2706]115        implicit none
116    !============================================================================
117    !   Load the adequate set of parameters for specname
118    !   From a table of traceurs
[2707]119    !
[2725]120    !   Also set up few others useful parameters, such as epsi_generic=m/mugaz, RLVTT_generic and
[2707]121    !   the metallicity_coeff.
122    !   Authors
123    !   --------
124    !   Noé Clément (2022)
[2706]125    !============================================================================
126
127        character(*), intent(in) :: specname
[2707]128        integer :: ios
[2706]129        character(len=500):: table_traceurs_line ! table_traceurs_line lines with parameters
[2786]130
[2710]131        open(117,file=trim(datadir)//'/table_tracers_condensable',form='formatted',status='old')
[2706]132
133        read(117,'(A)') table_traceurs_line
134
135        do
[2707]136            read(117,'(A)', iostat=ios) table_traceurs_line
[2706]137           
138            if (index(table_traceurs_line,specname) /= 0) then
139               
[2707]140                if (index(table_traceurs_line,'delta_vapH='   ) /= 0) then
141                    read(table_traceurs_line(index(table_traceurs_line,'delta_vapH=')+len('delta_vapH='):),*) delta_vapH
142                    print*, 'delta_vapH ', delta_vapH
[2706]143                end if
144                if (index(table_traceurs_line,'Tref='   ) /= 0) then
145                    read(table_traceurs_line(index(table_traceurs_line,'Tref=')+len('Tref='):),*) Tref
[2707]146                    print*, 'Tref', Tref
[2706]147                end if
148                if (index(table_traceurs_line,'Pref='   ) /= 0) then
149                    read(table_traceurs_line(index(table_traceurs_line,'Pref=')+len('Pref='):),*) Pref
[2707]150                    print*, 'Pref', Pref
[2706]151                end if
152                if (index(table_traceurs_line,'mass='   ) /= 0) then
153                    read(table_traceurs_line(index(table_traceurs_line,'mass=')+len('mass='):),*) m
[2707]154                    print*, 'mass', m
[2706]155                end if
156                if (index(table_traceurs_line,'metallicity_coeff='   ) /= 0) then
157                   read(table_traceurs_line(index(table_traceurs_line,'metallicity_coeff=')+len('metallicity_coeff='):),*) metallicity_coeff
[2707]158                   print*, 'metallicity_coeff', metallicity_coeff
[2706]159                end if
[2707]160
161                ios=1
[2706]162            end if
[2707]163
164            if (ios /= 0) exit
165           
[2706]166        end do
167
[2708]168        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
169
170            write(*,*), "ERROR : the following traceur is not referenced in the table_tracers_condensable : ", specname
171
172        end if
173       
[2706]174        close(117)
175
[2725]176        RLVTT_generic=delta_vapH/(m/1000.)
[2707]177
[2725]178        write(*,*) 'RLVTT_generic', RLVTT_generic
[2707]179
[2725]180        epsi_generic = m/mugaz
[2707]181
[2706]182    end subroutine specie_parameters_table
183
[2690]184    subroutine Psat_generic(T,p,metallicity,psat,qsat)
185        IMPLICIT NONE
186    !============================================================================
187    !   Clausius-Clapeyron relation :
188    !   d(Ln(psat))/dT = delta_vapH/(RT^2)
189    !   ->Psat = Pref * exp(-delta_vapH/R*(1/T-1/Tref))
190    !
191    !   Authors
192    !   --------
193    !   Lucas Teinturier (2022) adapted from Psat_water of Jeremy Leconte (2012)
194    !============================================================================
195        !Inputs
196        real, intent(in) :: T, p ! Temperature and pressure of the layer (in K and Pas)
197        real, intent(in) :: metallicity ! metallycity (log10)
198        !Outputs
199        real, intent(out) :: psat,qsat !saturation vapor pressure (Pa) and mass mixing ratio at saturation (kg/kg) of the layer
200
201        psat = pref*exp(-delta_vapH/(r*mugaz/1000.)*(1/T-1/Tref)-metallicity_coeff*metallicity) ! in Pa (because pref in Pa)
202
203        if (psat .gt. p) then
[2721]204            qsat = 1. ! is the maximum amount of vapor that a parcel can hold without condensation, in specific concentration.
[2690]205        else
[2725]206            qsat = epsi_generic *psat/(p-(1-epsi_generic)*psat)
[2690]207        endif
208    end subroutine Psat_generic
209
210    subroutine Lcpdqsat_generic(T,p,psat,qsat,dqsat,dlnpsat)
211        implicit none
212
213    !===============================================================================
214    !   Compute dqsat = L/cp* d(qsat)/dT and d(Ln(psat)) = L/cp*d(Ln(psat))/dT
215    !   we have d(ln(psat))/dT =  delta_vapH/R*(1/T^2) for clausius-clapeyron
216    !   r*mugaz/1000. is the perfect gaz constant in the computation of "dummy"
217    !   Authors
218    !   --------
219    !   Lucas Teinturier (2022) adapted from Lcpdqsat_water of Jeremy Leconte (2012)
220    !===============================================================================
221        ! Inputs
222        real, intent(in) :: T ! Temperature (K)
223        real, intent(in) :: p ! Pressure (Pa)
224        real, intent(in) :: psat ! Saturation vapor pressure (Pa)
225        real, intent(in) :: qsat ! Mass mixing ratio at saturation (kg/kg)
226
227        ! Outputs
228        real, intent(out) :: dqsat,dlnpsat
229
230        ! Variables
231        real :: dummy ! used to store d(ln(psat))/dT
232
233        dummy = delta_vapH/((r*mugaz/1000.)*(T**2))
234
235        if (psat .gt. p) then
236            dqsat =0.
237        else
[2725]238            dqsat = (RLVTT_generic/cpp) *qsat*(p/(p-(1-epsi_generic)*psat))*dummy
239            dlnpsat = (RLVTT_generic/cpp) * dummy
[2690]240        endif
241
242    end subroutine Lcpdqsat_generic
243
[2718]244    !==================================================================
245    subroutine Tsat_generic(p,Tsat)
246
247        implicit none
248
249    !==================================================================
250    !     Purpose
251    !     -------
252    !     Compute the saturation temperature
253    !     for a given pressure (Pa)
254    !
255    !     Authors
256    !     -------
257    !     Noé Clément (2022)
258    !
259    !==================================================================
260
261        ! input
262        real p
263   
264        ! output
265        real Tsat
266
267
268        Tsat = 1/(1/Tref - (r*mugaz/1000.)/delta_vapH*Log(p/Pref))
269
270    end subroutine Tsat_generic
271    !==================================================================
272
[2786]273end module generic_cloud_common_h
Note: See TracBrowser for help on using the repository browser.