| 1 | ! |
|---|
| 2 | ! $Header$ |
|---|
| 3 | ! |
|---|
| 4 | subroutine calltherm(dtime & |
|---|
| 5 | & ,pplay,paprs,pphi,weak_inversion & |
|---|
| 6 | & ,u_seri,v_seri,t_seri,q_seri,zqsat,debut & |
|---|
| 7 | & ,d_u_ajs,d_v_ajs,d_t_ajs,d_q_ajs & |
|---|
| 8 | & ,fm_therm,entr_therm,zqasc,clwcon0,lmax,ratqscth, & |
|---|
| 9 | & ratqsdiff,zqsatth,Ale_bl,Alp_bl,lalim_conv,wght_th, & |
|---|
| 10 | & zmax0,f0) |
|---|
| 11 | |
|---|
| 12 | implicit none |
|---|
| 13 | #include "dimensions.h" |
|---|
| 14 | #include "dimphy.h" |
|---|
| 15 | #include "thermcell.h" |
|---|
| 16 | #include "iniprint.h" |
|---|
| 17 | |
|---|
| 18 | ! A inclure eventuellement dans les fichiers de configuration |
|---|
| 19 | data r_aspect_thermals,l_mix_thermals,tho_thermals/2.,30.,0./ |
|---|
| 20 | data w2di_thermals/0/ |
|---|
| 21 | |
|---|
| 22 | REAL dtime |
|---|
| 23 | LOGICAL debut |
|---|
| 24 | REAL u_seri(klon,klev),v_seri(klon,klev) |
|---|
| 25 | REAL t_seri(klon,klev),q_seri(klon,klev),qmemoire(klon,klev) |
|---|
| 26 | REAL weak_inversion(klon) |
|---|
| 27 | REAL paprs(klon,klev+1) |
|---|
| 28 | REAL pplay(klon,klev) |
|---|
| 29 | REAL pphi(klon,klev) |
|---|
| 30 | real zlev(klon,klev+1) |
|---|
| 31 | !test: on sort lentr et a* pour alimenter KE |
|---|
| 32 | REAL wght_th(klon,klev) |
|---|
| 33 | INTEGER lalim_conv(klon) |
|---|
| 34 | |
|---|
| 35 | !FH Update Thermiques |
|---|
| 36 | REAL d_t_ajs(klon,klev), d_q_ajs(klon,klev) |
|---|
| 37 | REAL d_u_ajs(klon,klev),d_v_ajs(klon,klev) |
|---|
| 38 | real fm_therm(klon,klev+1),entr_therm(klon,klev) |
|---|
| 39 | |
|---|
| 40 | !******************************************************** |
|---|
| 41 | ! declarations |
|---|
| 42 | real fmc_therm(klon,klev+1),zqasc(klon,klev) |
|---|
| 43 | real zqla(klon,klev) |
|---|
| 44 | real wmax_sec(klon) |
|---|
| 45 | real zmax_sec(klon) |
|---|
| 46 | real f_sec(klon) |
|---|
| 47 | real detrc_therm(klon,klev) |
|---|
| 48 | save fmc_therm, detrc_therm |
|---|
| 49 | real clwcon0(klon,klev) |
|---|
| 50 | real zqsat(klon,klev) |
|---|
| 51 | real zw_sec(klon,klev+1) |
|---|
| 52 | integer lmix_sec(klon) |
|---|
| 53 | integer lmax(klon) |
|---|
| 54 | real ratqscth(klon,klev) |
|---|
| 55 | real ratqsdiff(klon,klev) |
|---|
| 56 | real zqsatth(klon,klev) |
|---|
| 57 | !nouvelles variables pour la convection |
|---|
| 58 | real Ale_bl(klon) |
|---|
| 59 | real Alp_bl(klon) |
|---|
| 60 | real Ale(klon) |
|---|
| 61 | real Alp(klon) |
|---|
| 62 | !RC |
|---|
| 63 | !on garde le zmax du pas de temps precedent |
|---|
| 64 | real zmax0(klon), f0(klon) |
|---|
| 65 | !******************************************************** |
|---|
| 66 | |
|---|
| 67 | |
|---|
| 68 | ! variables locales |
|---|
| 69 | REAL d_t_the(klon,klev), d_q_the(klon,klev) |
|---|
| 70 | REAL d_u_the(klon,klev),d_v_the(klon,klev) |
|---|
| 71 | ! |
|---|
| 72 | real zfm_therm(klon,klev+1),zentr_therm(klon,klev),zdt |
|---|
| 73 | save zentr_therm,zfm_therm |
|---|
| 74 | |
|---|
| 75 | integer i,k |
|---|
| 76 | !******************************************************** |
|---|
| 77 | |
|---|
| 78 | ! Modele du thermique |
|---|
| 79 | ! =================== |
|---|
| 80 | ! print*,'thermiques: WARNING on passe t au lieu de t_seri' |
|---|
| 81 | |
|---|
| 82 | |
|---|
| 83 | fm_therm(:,:)=0. |
|---|
| 84 | entr_therm(:,:)=0. |
|---|
| 85 | Ale_bl(:)=0. |
|---|
| 86 | Alp_bl(:)=0. |
|---|
| 87 | if (prt_level.ge.10) then |
|---|
| 88 | print*,'thermV4 nsplit: ',nsplit_thermals,' weak_inversion' |
|---|
| 89 | endif |
|---|
| 90 | |
|---|
| 91 | ! tests sur les valeurs negatives de l'eau |
|---|
| 92 | do k=1,klev |
|---|
| 93 | do i=1,klon |
|---|
| 94 | if (.not.q_seri(i,k).ge.0.) then |
|---|
| 95 | if (prt_level.ge.10) then |
|---|
| 96 | print*,'WARN eau<0 avant therm i=',i,' k=',k & |
|---|
| 97 | & ,' dq,q',d_q_the(i,k),q_seri(i,k) |
|---|
| 98 | endif |
|---|
| 99 | q_seri(i,k)=1.e-15 |
|---|
| 100 | endif |
|---|
| 101 | enddo |
|---|
| 102 | enddo |
|---|
| 103 | |
|---|
| 104 | |
|---|
| 105 | zdt=dtime/float(nsplit_thermals) |
|---|
| 106 | do isplit=1,nsplit_thermals |
|---|
| 107 | |
|---|
| 108 | if (iflag_thermals.eq.1) then |
|---|
| 109 | CALL thermcell_2002(klon,klev,zdt & |
|---|
| 110 | & ,pplay,paprs,pphi & |
|---|
| 111 | & ,u_seri,v_seri,t_seri,q_seri & |
|---|
| 112 | & ,d_u_the,d_v_the,d_t_the,d_q_the & |
|---|
| 113 | & ,zfm_therm,zentr_therm & |
|---|
| 114 | & ,r_aspect_thermals,30.,w2di_thermals & |
|---|
| 115 | & ,tho_thermals,3) |
|---|
| 116 | else if (iflag_thermals.eq.2) then |
|---|
| 117 | CALL thermcell_sec(klon,klev,zdt & |
|---|
| 118 | & ,pplay,paprs,pphi,zlev & |
|---|
| 119 | & ,u_seri,v_seri,t_seri,q_seri & |
|---|
| 120 | & ,d_u_the,d_v_the,d_t_the,d_q_the & |
|---|
| 121 | & ,zfm_therm,zentr_therm & |
|---|
| 122 | & ,r_aspect_thermals,30.,w2di_thermals & |
|---|
| 123 | & ,tho_thermals,3) |
|---|
| 124 | else if (iflag_thermals.eq.3) then |
|---|
| 125 | CALL thermcell(klon,klev,zdt & |
|---|
| 126 | & ,pplay,paprs,pphi & |
|---|
| 127 | & ,u_seri,v_seri,t_seri,q_seri & |
|---|
| 128 | & ,d_u_the,d_v_the,d_t_the,d_q_the & |
|---|
| 129 | & ,zfm_therm,zentr_therm & |
|---|
| 130 | & ,r_aspect_thermals,l_mix_thermals,w2di_thermals & |
|---|
| 131 | & ,tho_thermals,3) |
|---|
| 132 | else if (iflag_thermals.eq.10) then |
|---|
| 133 | CALL thermcell_eau(klon,klev,zdt & |
|---|
| 134 | & ,pplay,paprs,pphi & |
|---|
| 135 | & ,u_seri,v_seri,t_seri,q_seri & |
|---|
| 136 | & ,d_u_the,d_v_the,d_t_the,d_q_the & |
|---|
| 137 | & ,zfm_therm,zentr_therm & |
|---|
| 138 | & ,r_aspect_thermals,l_mix_thermals,w2di_thermals & |
|---|
| 139 | & ,tho_thermals,3) |
|---|
| 140 | else if (iflag_thermals.eq.11) then |
|---|
| 141 | stop'cas non prevu dans calltherm' |
|---|
| 142 | ! CALL thermcell_pluie(klon,klev,zdt & |
|---|
| 143 | ! & ,pplay,paprs,pphi,zlev & |
|---|
| 144 | ! & ,u_seri,v_seri,t_seri,q_seri & |
|---|
| 145 | ! & ,d_u_the,d_v_the,d_t_the,d_q_the & |
|---|
| 146 | ! & ,zfm_therm,zentr_therm,zqla & |
|---|
| 147 | ! & ,r_aspect_thermals,l_mix_thermals,w2di_thermals & |
|---|
| 148 | ! & ,tho_thermals,3) |
|---|
| 149 | else if (iflag_thermals.eq.12) then |
|---|
| 150 | CALL calcul_sec(klon,klev,zdt & |
|---|
| 151 | & ,pplay,paprs,pphi,zlev & |
|---|
| 152 | & ,u_seri,v_seri,t_seri,q_seri & |
|---|
| 153 | & ,zmax_sec,wmax_sec,zw_sec,lmix_sec & |
|---|
| 154 | & ,r_aspect_thermals,l_mix_thermals,w2di_thermals & |
|---|
| 155 | & ,tho_thermals) |
|---|
| 156 | ! CALL calcul_sec_entr(klon,klev,zdt |
|---|
| 157 | ! s ,pplay,paprs,pphi,zlev,debut |
|---|
| 158 | ! s ,u_seri,v_seri,t_seri,q_seri |
|---|
| 159 | ! s ,zmax_sec,wmax_sec,zw_sec,lmix_sec |
|---|
| 160 | ! s ,r_aspect_thermals,l_mix_thermals,w2di_thermals |
|---|
| 161 | ! s ,tho_thermals,3) |
|---|
| 162 | ! CALL thermcell_pluie_detr(klon,klev,zdt & |
|---|
| 163 | ! & ,pplay,paprs,pphi,zlev,debut & |
|---|
| 164 | ! & ,u_seri,v_seri,t_seri,q_seri & |
|---|
| 165 | ! & ,d_u_the,d_v_the,d_t_the,d_q_the & |
|---|
| 166 | ! & ,zfm_therm,zentr_therm,zqla,lmax & |
|---|
| 167 | ! & ,zmax_sec,wmax_sec,zw_sec,lmix_sec & |
|---|
| 168 | ! & ,ratqscth,ratqsdiff,zqsatth & |
|---|
| 169 | ! & ,r_aspect_thermals,l_mix_thermals,w2di_thermals & |
|---|
| 170 | ! & ,tho_thermals) |
|---|
| 171 | else if (iflag_thermals.ge.13) then |
|---|
| 172 | CALL thermcell_main(klon,klev,zdt & |
|---|
| 173 | & ,pplay,paprs,pphi,debut & |
|---|
| 174 | & ,u_seri,v_seri,t_seri,q_seri & |
|---|
| 175 | & ,d_u_the,d_v_the,d_t_the,d_q_the & |
|---|
| 176 | & ,zfm_therm,zentr_therm,zqla,lmax & |
|---|
| 177 | & ,ratqscth,ratqsdiff,zqsatth & |
|---|
| 178 | & ,r_aspect_thermals,l_mix_thermals,w2di_thermals & |
|---|
| 179 | & ,tho_thermals,Ale,Alp,lalim_conv,wght_th & |
|---|
| 180 | & ,zmax0,f0) |
|---|
| 181 | endif |
|---|
| 182 | |
|---|
| 183 | |
|---|
| 184 | DO i=1,klon |
|---|
| 185 | DO k=1,klev |
|---|
| 186 | IF(iflag_thermals.lt.14.or.weak_inversion(i).gt.0.5) THEN |
|---|
| 187 | |
|---|
| 188 | ! transformation de la derivee en tendance |
|---|
| 189 | d_t_the(i,k)=d_t_the(i,k)*dtime/float(nsplit_thermals) |
|---|
| 190 | d_u_the(i,k)=d_u_the(i,k)*dtime/float(nsplit_thermals) |
|---|
| 191 | d_v_the(i,k)=d_v_the(i,k)*dtime/float(nsplit_thermals) |
|---|
| 192 | d_q_the(i,k)=d_q_the(i,k)*dtime/float(nsplit_thermals) |
|---|
| 193 | fm_therm(i,k)=fm_therm(i,k) & |
|---|
| 194 | & +zfm_therm(i,k)/float(nsplit_thermals) |
|---|
| 195 | entr_therm(i,k)=entr_therm(i,k) & |
|---|
| 196 | & +zentr_therm(i,k)/float(nsplit_thermals) |
|---|
| 197 | fm_therm(:,klev+1)=0. |
|---|
| 198 | |
|---|
| 199 | |
|---|
| 200 | |
|---|
| 201 | ! accumulation de la tendance |
|---|
| 202 | d_t_ajs(i,k)=d_t_ajs(i,k)+d_t_the(i,k) |
|---|
| 203 | d_u_ajs(i,k)=d_u_ajs(i,k)+d_u_the(i,k) |
|---|
| 204 | d_v_ajs(i,k)=d_v_ajs(i,k)+d_v_the(i,k) |
|---|
| 205 | d_q_ajs(i,k)=d_q_ajs(i,k)+d_q_the(i,k) |
|---|
| 206 | |
|---|
| 207 | ! incrementation des variables meteo |
|---|
| 208 | t_seri(i,k) = t_seri(i,k) + d_t_the(i,k) |
|---|
| 209 | u_seri(i,k) = u_seri(i,k) + d_u_the(i,k) |
|---|
| 210 | v_seri(i,k) = v_seri(i,k) + d_v_the(i,k) |
|---|
| 211 | qmemoire(i,k)=q_seri(i,k) |
|---|
| 212 | q_seri(i,k) = q_seri(i,k) + d_q_the(i,k) |
|---|
| 213 | ENDIF |
|---|
| 214 | ENDDO |
|---|
| 215 | ENDDO |
|---|
| 216 | |
|---|
| 217 | DO i=1,klon |
|---|
| 218 | fm_therm(i,klev+1)=0. |
|---|
| 219 | Ale_bl(i)=Ale_bl(i)+Ale(i)/float(nsplit_thermals) |
|---|
| 220 | ! write(22,*)'ALE CALLTHERM',Ale_bl(i),Ale(i) |
|---|
| 221 | Alp_bl(i)=Alp_bl(i)+Alp(i)/float(nsplit_thermals) |
|---|
| 222 | ! write(23,*)'ALP CALLTHERM',Alp_bl(i),Alp(i) |
|---|
| 223 | ENDDO |
|---|
| 224 | |
|---|
| 225 | ! tests sur les valeurs negatives de l'eau |
|---|
| 226 | DO k = 1, klev |
|---|
| 227 | DO i = 1, klon |
|---|
| 228 | if (.not.q_seri(i,k).ge.0.) then |
|---|
| 229 | if (prt_level.ge.10) then |
|---|
| 230 | print*,'WARN eau<0 apres therm i=',i,' k=',k & |
|---|
| 231 | & ,' dq,q',d_q_the(i,k),q_seri(i,k), & |
|---|
| 232 | & 'fm=',zfm_therm(i,k),'entr=',entr_therm(i,k) |
|---|
| 233 | endif |
|---|
| 234 | q_seri(i,k)=1.e-15 |
|---|
| 235 | ! stop |
|---|
| 236 | endif |
|---|
| 237 | ENDDO |
|---|
| 238 | ENDDO |
|---|
| 239 | ! tests sur les valeurs de la temperature |
|---|
| 240 | DO k = 1, klev |
|---|
| 241 | DO i = 1, klon |
|---|
| 242 | if ((t_seri(i,k).lt.50.) .or. & |
|---|
| 243 | & (t_seri(i,k).gt.370.)) then |
|---|
| 244 | print*,'WARN temp apres therm i=',i,' k=',k & |
|---|
| 245 | & ,' t_seri',t_seri(i,k) |
|---|
| 246 | ! CALL abort |
|---|
| 247 | endif |
|---|
| 248 | ENDDO |
|---|
| 249 | ENDDO |
|---|
| 250 | |
|---|
| 251 | enddo ! isplit |
|---|
| 252 | |
|---|
| 253 | ! |
|---|
| 254 | !*************************************************************** |
|---|
| 255 | ! calcul du flux ascencant conservatif |
|---|
| 256 | ! print*,'<<<<calcul flux ascendant conservatif' |
|---|
| 257 | |
|---|
| 258 | fmc_therm=0. |
|---|
| 259 | do k=1,klev |
|---|
| 260 | do i=1,klon |
|---|
| 261 | if (entr_therm(i,k).gt.0.) then |
|---|
| 262 | fmc_therm(i,k+1)=fmc_therm(i,k)+entr_therm(i,k) |
|---|
| 263 | else |
|---|
| 264 | fmc_therm(i,k+1)=fmc_therm(i,k) |
|---|
| 265 | endif |
|---|
| 266 | detrc_therm(i,k)=(fmc_therm(i,k+1)-fm_therm(i,k+1)) & |
|---|
| 267 | & -(fmc_therm(i,k)-fm_therm(i,k)) |
|---|
| 268 | enddo |
|---|
| 269 | enddo |
|---|
| 270 | |
|---|
| 271 | |
|---|
| 272 | !**************************************************************** |
|---|
| 273 | ! calcul de l'humidite dans l'ascendance |
|---|
| 274 | ! print*,'<<<<calcul de lhumidite dans thermique' |
|---|
| 275 | !CR:on ne le calcule que pour le cas sec |
|---|
| 276 | if (iflag_thermals.le.11) then |
|---|
| 277 | do i=1,klon |
|---|
| 278 | zqasc(i,1)=q_seri(i,1) |
|---|
| 279 | do k=2,klev |
|---|
| 280 | if (fmc_therm(i,k+1).gt.1.e-6) then |
|---|
| 281 | zqasc(i,k)=(fmc_therm(i,k)*zqasc(i,k-1) & |
|---|
| 282 | & +entr_therm(i,k)*q_seri(i,k))/fmc_therm(i,k+1) |
|---|
| 283 | !CR:test on asseche le thermique |
|---|
| 284 | ! zqasc(i,k)=zqasc(i,k)/2. |
|---|
| 285 | ! else |
|---|
| 286 | ! zqasc(i,k)=q_seri(i,k) |
|---|
| 287 | endif |
|---|
| 288 | enddo |
|---|
| 289 | enddo |
|---|
| 290 | |
|---|
| 291 | |
|---|
| 292 | ! calcul de l'eau condensee dans l'ascendance |
|---|
| 293 | ! print*,'<<<<calcul de leau condensee dans thermique' |
|---|
| 294 | do i=1,klon |
|---|
| 295 | do k=1,klev |
|---|
| 296 | clwcon0(i,k)=zqasc(i,k)-zqsat(i,k) |
|---|
| 297 | if (clwcon0(i,k).lt.0. .or. & |
|---|
| 298 | & (fm_therm(i,k+1)+detrc_therm(i,k)).lt.1.e-6) then |
|---|
| 299 | clwcon0(i,k)=0. |
|---|
| 300 | endif |
|---|
| 301 | enddo |
|---|
| 302 | enddo |
|---|
| 303 | else |
|---|
| 304 | do i=1,klon |
|---|
| 305 | do k=1,klev |
|---|
| 306 | clwcon0(i,k)=zqla(i,k) |
|---|
| 307 | if (clwcon0(i,k).lt.0. .or. & |
|---|
| 308 | & (fm_therm(i,k+1)+detrc_therm(i,k)).lt.1.e-6) then |
|---|
| 309 | clwcon0(i,k)=0. |
|---|
| 310 | endif |
|---|
| 311 | enddo |
|---|
| 312 | enddo |
|---|
| 313 | endif |
|---|
| 314 | !******************************************************************* |
|---|
| 315 | |
|---|
| 316 | |
|---|
| 317 | return |
|---|
| 318 | |
|---|
| 319 | end |
|---|