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