| 1 | module condensation_generic_mod |
|---|
| 2 | implicit none |
|---|
| 3 | |
|---|
| 4 | contains |
|---|
| 5 | |
|---|
| 6 | subroutine condensation_generic(ngrid,nlayer,nq,ptimestep, pplev, pplay, & |
|---|
| 7 | pt, pq, pdt, pdq, pdtlsc, pdqvaplsc, pdqliqlsc, rneb) |
|---|
| 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 | REAL, intent(out) :: rneb(ngrid,nlayer,nq) ! fraction nuageuse |
|---|
| 46 | |
|---|
| 47 | ! Options : |
|---|
| 48 | real, save :: metallicity !metallicity of planet |
|---|
| 49 | REAL, SAVE :: qvap_deep ! deep mixing ratio of water vapor when simulating bottom less planets |
|---|
| 50 | !$OMP THREADPRIVATE(metallicity, qvap_deep) |
|---|
| 51 | |
|---|
| 52 | ! Local variables |
|---|
| 53 | |
|---|
| 54 | ! to call only one time the ice/vap pair of a tracer |
|---|
| 55 | logical call_ice_vap_generic |
|---|
| 56 | |
|---|
| 57 | INTEGER i, k , nn, iq |
|---|
| 58 | INTEGER,PARAMETER :: nitermax=5000 |
|---|
| 59 | INTEGER,PARAMETER :: tau=86400 ! tau is in seconds and must not be lower than the physical step time. |
|---|
| 60 | DOUBLE PRECISION,PARAMETER :: alpha=.1,qthreshold=1.d-8 |
|---|
| 61 | ! JL13: if "careful, T<Tmin in psat water" appears often, you may want to stabilise the model by |
|---|
| 62 | ! decreasing alpha and increasing nitermax accordingly |
|---|
| 63 | DOUBLE PRECISION zq(ngrid) |
|---|
| 64 | DOUBLE PRECISION zcond(ngrid),zcond_iter |
|---|
| 65 | DOUBLE PRECISION zqs(ngrid) |
|---|
| 66 | real zt(ngrid),local_p,psat_tmp,dlnpsat_tmp,Lcp,zqs_temp,zdqs |
|---|
| 67 | integer igcm_generic_vap, igcm_generic_ice! index of the vap and ice of generic_tracer |
|---|
| 68 | ! CHARACTER(len=*) :: tname_ice |
|---|
| 69 | ! evaporation calculations |
|---|
| 70 | REAL dqevap(ngrid,nlayer),dtevap(ngrid,nlayer) |
|---|
| 71 | REAL qevap(ngrid,nlayer,nq) |
|---|
| 72 | REAL tevap(ngrid,nlayer) |
|---|
| 73 | |
|---|
| 74 | DOUBLE PRECISION zx_q(ngrid) |
|---|
| 75 | DOUBLE PRECISION zy_q(ngrid) |
|---|
| 76 | LOGICAL,SAVE :: firstcall=.true. |
|---|
| 77 | !$OMP THREADPRIVATE(firstcall) |
|---|
| 78 | IF (firstcall) THEN |
|---|
| 79 | write(*,*) "value for metallicity? " |
|---|
| 80 | metallicity=0.0 ! default value |
|---|
| 81 | call getin_p("metallicity",metallicity) |
|---|
| 82 | write(*,*) " metallicity = ",metallicity |
|---|
| 83 | |
|---|
| 84 | write(*,*) "Deep water vapor mixing ratio ? (no effect if negative) " |
|---|
| 85 | qvap_deep=-1. ! default value |
|---|
| 86 | call getin_p("qvap_deep",qvap_deep) |
|---|
| 87 | write(*,*) " qvap_deep = ",qvap_deep |
|---|
| 88 | |
|---|
| 89 | firstcall = .false. |
|---|
| 90 | ENDIF |
|---|
| 91 | ! Initialisation of outputs and local variables |
|---|
| 92 | pdtlsc(1:ngrid,1:nlayer) = 0.0 |
|---|
| 93 | pdqvaplsc(1:ngrid,1:nlayer,1:nq) = 0.0 |
|---|
| 94 | pdqliqlsc(1:ngrid,1:nlayer,1:nq) = 0.0 |
|---|
| 95 | dqevap(1:ngrid,1:nlayer)=0.0 |
|---|
| 96 | dtevap(1:ngrid,1:nlayer)=0.0 |
|---|
| 97 | qevap(1:ngrid,1:nlayer,1:nq)=0.0 |
|---|
| 98 | tevap(1:ngrid,1:nlayer)=0.0 |
|---|
| 99 | rneb(1:ngrid,1:nlayer,1:nq) = 0.0 |
|---|
| 100 | ! Let's loop on tracers |
|---|
| 101 | do iq=1,nq |
|---|
| 102 | |
|---|
| 103 | call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic) |
|---|
| 104 | |
|---|
| 105 | if(call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer |
|---|
| 106 | m=constants_mass(iq) |
|---|
| 107 | delta_vapH=constants_delta_vapH(iq) |
|---|
| 108 | Tref=constants_Tref(iq) |
|---|
| 109 | Pref=constants_Pref(iq) |
|---|
| 110 | epsi_generic=constants_epsi_generic(iq) |
|---|
| 111 | RLVTT_generic=constants_RLVTT_generic(iq) |
|---|
| 112 | metallicity_coeff=constants_metallicity_coeff(iq) |
|---|
| 113 | |
|---|
| 114 | Lcp=RLVTT_generic/cpp ! need to be init here |
|---|
| 115 | |
|---|
| 116 | ! Vertical loop (from top to bottom) |
|---|
| 117 | DO k = nlayer, 1, -1 |
|---|
| 118 | zt(1:ngrid)=pt(1:ngrid,k)+pdt(1:ngrid,k)*ptimestep |
|---|
| 119 | |
|---|
| 120 | ! Computes Psat and the partial condensation |
|---|
| 121 | DO i = 1, ngrid |
|---|
| 122 | local_p=pplay(i,k) |
|---|
| 123 | if(zt(i).le.15.) then |
|---|
| 124 | print*,'in lsc',i,k,zt(i) |
|---|
| 125 | ! zt(i)=15. ! check too low temperatures |
|---|
| 126 | endif |
|---|
| 127 | zx_q(i) = pq(i,k,igcm_generic_vap)+pdq(i,k,igcm_generic_vap)*ptimestep |
|---|
| 128 | |
|---|
| 129 | call Psat_generic(zt(i),local_p,metallicity,psat_tmp,zqs_temp) |
|---|
| 130 | zy_q(i) = pq(i,k,igcm_generic_ice)+pdq(i,k,igcm_generic_ice)*ptimestep |
|---|
| 131 | |
|---|
| 132 | if ((zx_q(i) .le. zqs_temp) .and. (zy_q(i) .eq. 0.)) then |
|---|
| 133 | ! if we are are not saturated and if there is no ice |
|---|
| 134 | ! then no change |
|---|
| 135 | |
|---|
| 136 | zcond(i) = 0.0d0 |
|---|
| 137 | |
|---|
| 138 | else ! if we are saturated : we must evaporate |
|---|
| 139 | ! if there is ice : we must check if we can condensate |
|---|
| 140 | |
|---|
| 141 | ! iterative process to stabilize the scheme when large water amounts JL12 |
|---|
| 142 | zcond(i) = 0.0d0 |
|---|
| 143 | Do nn=1,nitermax |
|---|
| 144 | call Psat_generic(zt(i),local_p,metallicity,psat_tmp,zqs_temp) |
|---|
| 145 | zqs(i)=zqs_temp |
|---|
| 146 | call Lcpdqsat_generic(zt(i),local_p,psat_tmp,zqs_temp,zdqs,dlnpsat_tmp) |
|---|
| 147 | zcond_iter = alpha*(zx_q(i)-zqs(i))/(1.d0+zdqs) |
|---|
| 148 | !zcond can be negative here |
|---|
| 149 | zx_q(i) = zx_q(i) - zcond_iter |
|---|
| 150 | zcond(i) = zcond(i) + zcond_iter |
|---|
| 151 | zt(i) = zt(i) + zcond_iter*Lcp |
|---|
| 152 | if (ABS(zcond_iter/alpha/zqs(i)).lt.qthreshold) exit |
|---|
| 153 | if (nn.eq.nitermax) print*,'itermax in largescale' |
|---|
| 154 | End do ! niter |
|---|
| 155 | |
|---|
| 156 | ! if zcond(i) > 0, zcond(i) is the amount of vapor that we can condensate |
|---|
| 157 | ! because we are saturated. zcond(i) will not change below |
|---|
| 158 | ! if zcond(i) < 0, zcond(i) is the amount of ice that we can evaporate. |
|---|
| 159 | ! We can not evaporate more than the existing ice, |
|---|
| 160 | ! so the line below is to check how much we can evaporate. |
|---|
| 161 | ! If there is no ice available, zcond(i) will be 0. (first condidition of "if") |
|---|
| 162 | |
|---|
| 163 | zcond(i)=MAX(zcond(i),-(pq(i,k,igcm_generic_ice)+pdq(i,k,igcm_generic_ice)*ptimestep)) |
|---|
| 164 | |
|---|
| 165 | endif |
|---|
| 166 | |
|---|
| 167 | if (zcond(i) .gt. 0.) then |
|---|
| 168 | rneb(i,k,iq)=1 |
|---|
| 169 | else |
|---|
| 170 | rneb(i,k,iq)=0. |
|---|
| 171 | endif |
|---|
| 172 | |
|---|
| 173 | zcond(i) = zcond(i)/ptimestep |
|---|
| 174 | ENDDO ! i=1,ngrid |
|---|
| 175 | |
|---|
| 176 | !Tendances de t et q |
|---|
| 177 | pdqvaplsc(1:ngrid,k,igcm_generic_vap) = - zcond(1:ngrid) |
|---|
| 178 | pdqliqlsc(1:ngrid,k,igcm_generic_ice) = - pdqvaplsc(1:ngrid,k,igcm_generic_vap) |
|---|
| 179 | pdtlsc(1:ngrid,k) = pdtlsc(1:ngrid,k) + pdqliqlsc(1:ngrid,k,igcm_generic_ice)*Lcp |
|---|
| 180 | |
|---|
| 181 | Enddo ! k= nlayer, 1, -1 |
|---|
| 182 | |
|---|
| 183 | if (qvap_deep >= 0.) then |
|---|
| 184 | !brings lower generic vapor ratio to a fixed value. |
|---|
| 185 | ! tau is in seconds and must not be lower than the physical step time. |
|---|
| 186 | pdqvaplsc(1:ngrid,1,igcm_generic_vap) = (qvap_deep - pq(1:ngrid,1,igcm_generic_vap))/tau - pdq(1:ngrid,1,igcm_generic_vap) |
|---|
| 187 | endif |
|---|
| 188 | endif !if(call_ice_vap_generic) |
|---|
| 189 | enddo ! iq=1,nq |
|---|
| 190 | |
|---|
| 191 | end subroutine condensation_generic |
|---|
| 192 | end module condensation_generic_mod |
|---|