source: trunk/LMDZ.GENERIC/libf/phystd/condensation_generic_mod.F90 @ 2720

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

new module to get back indexes of cond tracers

File size: 8.0 KB
Line 
1module condensation_generic_mod
2    implicit none
3   
4contains
5   
6    subroutine condensation_generic(ngrid,nlayer,nq,ptimestep, pplev, pplay,    &
7                pt, pq, pdt, pdq, pdtlsc, pdqvaplsc, pdqliqlsc)
8        use ioipsl_getin_p_mod, only: getin_p !-> to get the metallicity
9        use generic_cloud_common_h
10        USE tracer_h
11        USE mod_phys_lmdz_para, only: is_master
12        use generic_tracer_index_mod, only: generic_tracer_index
13        IMPLICIT none
14
15!=======================================================================
16!     
17!     Purpose
18!     -------
19!     Calculates large-scale condensation of generic tracer "tname".
20!     By convention, tname ends with the suffix "_vap", as it represents the
21!     gas phase of the generic tracer
22!     
23!     Authors
24!     -------
25!     Adapted from largescale.F90 by Lucas Teinturier & Noé Clément (2022)
26!     largescale.F90 adapted from the LMDTERRE code by R. Wordsworth (2009)
27!     Original author Z. X. Li (1993)
28!     
29!=========================================================================
30
31        INTEGER, intent(in) :: ngrid,nlayer,nq
32
33!       Arguments
34        REAL, intent(in) :: ptimestep             ! intervalle du temps (s)
35        REAL, intent(in) :: pplev(ngrid,nlayer+1) ! pression a inter-couche
36        REAL, intent(in) :: pplay(ngrid,nlayer)   ! pression au milieu de couche
37        REAL, intent(in) :: pt(ngrid,nlayer)      ! temperature (K)
38        REAL, intent(in) :: pq(ngrid,nlayer,nq)   ! tracer mixing ratio (kg/kg)
39        REAL, intent(in) :: pdt(ngrid,nlayer)     ! physical temperature tendency (K/s)
40        REAL, intent(in) :: pdq(ngrid,nlayer,nq)  ! physical tracer tendency (K/s)
41        ! CHARACTER(*), intent(in) :: tname_vap     ! name of the tracer we consider. BY convention, it ends with _vap !!!
42        REAL, intent(out) :: pdtlsc(ngrid,nlayer)  ! incrementation de la temperature (K)
43        REAL, intent(out) :: pdqvaplsc(ngrid,nlayer,nq) ! incrementation de la vapeur du traceur
44        REAL, intent(out) :: pdqliqlsc(ngrid,nlayer,nq) ! incrementation du traceur liquide
45
46!       Options :
47        real, save :: metallicity !metallicity of planet
48        REAL, SAVE :: qvap_deep   ! deep mixing ratio of water vapor when simulating bottom less planets
49!$OMP THREADPRIVATE(metallicity, qvap_deep)
50
51!       Local variables
52
53        ! to call only one time the ice/vap pair of a tracer
54        logical one_call_ice_vap_generic
55
56        INTEGER i, k , nn, iq
57        INTEGER,PARAMETER :: nitermax=5000
58        INTEGER,PARAMETER :: tau=14400.
59        DOUBLE PRECISION,PARAMETER :: alpha=.1,qthreshold=1.d-8
60        ! JL13: if "careful, T<Tmin in psat water" appears often, you may want to stabilise the model by
61        !                   decreasing alpha and increasing nitermax accordingly
62        DOUBLE PRECISION zq(ngrid)
63        DOUBLE PRECISION zcond(ngrid),zcond_iter
64        DOUBLE PRECISION zqs(ngrid)
65        real zt(ngrid),local_p,psat_tmp,dlnpsat_tmp,Lcp,zqs_temp,zdqs
66        integer igcm_generic_vap, igcm_generic_ice! index of the vap and ice of generic_tracer
67        ! CHARACTER(len=*) :: tname_ice
68        ! evaporation calculations
69        REAL dqevap(ngrid,nlayer),dtevap(ngrid,nlayer)     
70        REAL qevap(ngrid,nlayer,nq)
71        REAL tevap(ngrid,nlayer)
72
73        DOUBLE PRECISION zx_q(ngrid)
74        LOGICAL,SAVE :: firstcall=.true.
75!$OMP THREADPRIVATE(firstcall)
76        IF (firstcall) THEN
77                write(*,*) "value for metallicity? "
78                metallicity=0.0 ! default value
79                call getin_p("metallicity",metallicity)
80                write(*,*) " metallicity = ",metallicity
81
82                write(*,*) "Deep water vapor mixing ratio ? (no effect if negative) "
83                qvap_deep=-1. ! default value
84                call getin_p("qvap_deep",qvap_deep)
85                write(*,*) " qvap_deep = ",qvap_deep   
86
87                firstcall = .false.
88        ENDIF
89!       Initialisation of outputs and local variables
90        pdtlsc(1:ngrid,1:nlayer)  = 0.0
91        pdqvaplsc(1:ngrid,1:nlayer,1:nq)  = 0.0
92        pdqliqlsc(1:ngrid,1:nlayer,1:nq) = 0.0
93        dqevap(1:ngrid,1:nlayer)=0.0
94        dtevap(1:ngrid,1:nlayer)=0.0
95        qevap(1:ngrid,1:nlayer,1:nq)=0.0
96        tevap(1:ngrid,1:nlayer)=0.0
97        ! Let's loop on tracers
98        do iq=1,nq
99                one_call_ice_vap_generic = .false.
100
101                call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,one_call_ice_vap_generic)
102
103                if(one_call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer
104                        m=constants_mass(iq)
105                        delta_vapH=constants_delta_vapH(iq)
106                        Tref=constants_Tref(iq)
107                        Pref=constants_Pref(iq)
108                        epsi=constants_epsi(iq)
109                        RLVTT=constants_RLVTT(iq)
110                        metallicity_coeff=constants_metallicity_coeff(iq)
111
112                        Lcp=RLVTT/cpp ! need to be init here
113
114                        !  Vertical loop (from top to bottom)
115                        DO k = nlayer, 1, -1
116                                zt(1:ngrid)=pt(1:ngrid,k)+pdt(1:ngrid,k)*ptimestep
117
118                                ! Computes Psat and the partial condensation
119                                DO i = 1, ngrid
120                                        local_p=pplay(i,k)
121                                        if(zt(i).le.15.) then
122                                                print*,'in lsc',i,k,zt(i)
123                                        !           zt(i)=15.   ! check too low temperatures
124                                        endif
125                                        zx_q(i) = pq(i,k,igcm_generic_vap)+pdq(i,k,igcm_generic_vap)*ptimestep
126                                        ! iterative process to stabilize the scheme when large water amounts JL12
127                                        zcond(i) = 0.0d0
128                                        Do nn=1,nitermax 
129                                                call Psat_generic(zt(i),local_p,metallicity,psat_tmp,zqs_temp)
130                                                zqs(i)=zqs_temp
131                                                call Lcpdqsat_generic(zt(i),local_p,psat_tmp,zqs_temp,zdqs,dlnpsat_tmp)
132                                                zcond_iter = alpha*(zx_q(i)-zqs(i))/(1.d0+zdqs)
133                                                !zcond can be negative here
134                                                zx_q(i) = zx_q(i) - zcond_iter
135                                                zcond(i) = zcond(i) + zcond_iter
136                                                zt(i) = zt(i) + zcond_iter*Lcp
137                                                if (ABS(zcond_iter/alpha/zqs(i)).lt.qthreshold) exit
138                                                if (nn.eq.nitermax) print*,'itermax in largescale'
139                                        End do ! niter
140                                        zcond(i)=MAX(zcond(i),-(pq(i,k,igcm_generic_ice)+pdq(i,k,igcm_generic_ice)*ptimestep))
141                                        zcond(i) = zcond(i)/ptimestep
142                                ENDDO ! i=1,ngrid
143
144                                !Tendances de t et q
145                                pdqvaplsc(1:ngrid,k,igcm_generic_vap)  = - zcond(1:ngrid)
146                                pdqliqlsc(1:ngrid,k,igcm_generic_ice) = - pdqvaplsc(1:ngrid,k,igcm_generic_vap)
147                                pdtlsc(1:ngrid,k)  = pdtlsc(1:ngrid,k) + pdqliqlsc(1:ngrid,k,igcm_generic_ice)*Lcp
148
149                        Enddo ! k= nlayer, 1, -1
150
151                        if (qvap_deep >= 0.) then
152                                !brings lower generic vapor ratio to a fixed value.
153                                ! tau=3600. seems too fast
154                                pdqvaplsc(1:ngrid,1,igcm_generic_vap) = (qvap_deep - pq(1:ngrid,1,igcm_generic_vap))/tau - pdq(1:ngrid,1,igcm_generic_vap)
155                        endif
156                endif !if(one_call_ice_vap_generic)
157        enddo ! iq=1,nq
158
159    end subroutine condensation_generic
160end module condensation_generic_mod
Note: See TracBrowser for help on using the repository browser.