Changeset 863 for trunk/LMDZ.GENERIC/libf/phystd
- Timestamp:
- Jan 23, 2013, 2:09:25 PM (12 years ago)
- Location:
- trunk/LMDZ.GENERIC/libf/phystd
- Files:
-
- 7 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/libf/phystd/ave_stelspec.F90
r786 r863 98 98 Case(6) 99 99 file_id='/stellar_spectra/BD_Teff-1600K.txt' 100 tstellar=1600. ! Segura et al. (2003)100 tstellar=1600. 101 101 file_id_lam='/stellar_spectra/lamBD.txt' 102 102 Nfine=5000 103 103 Case(7) 104 104 file_id='/stellar_spectra/BD_Teff-1000K.txt' 105 tstellar=1000. ! Segura et al. (2003)105 tstellar=1000. 106 106 file_id_lam='/stellar_spectra/lamBD.txt' 107 107 Nfine=5000 108 Case(8) 109 file_id='/stellar_spectra/Flux_K5_Teff4700_logg4.5_Met-0.5_BTsettle.dat' 110 tstellar=4700. 111 file_id_lam='/stellar_spectra/lambda_K5_Teff4700_logg4.5_Met-0.5_BTsettle.dat' 112 Nfine=3986 108 113 Case Default 109 114 print*,'Error: unknown star type chosen' -
trunk/LMDZ.GENERIC/libf/phystd/dimphys.h
r787 r863 9 9 ! nsoilmx : number of subterranean layers 10 10 !integer, parameter :: nsoilmx = 4 ! for a test 11 integer, parameter :: nsoilmx = 18 ! for z1=0.0002 cm, depth = 18 m => mars case12 !integer, parameter :: nsoilmx = 13 ! for z1=0.03 cm, depth = 104.8 m => earth case11 integer, parameter :: nsoilmx = 18 ! for z1=0.0002 m, depth = 18 m => mars case 12 !integer, parameter :: nsoilmx = 13 ! for z1=0.03 m, depth = 104.8 m => earth case 13 13 !----------------------------------------------------------------------- 14 14 -
trunk/LMDZ.GENERIC/libf/phystd/hydrol.F90
r787 r863 91 91 real ztsurf(ngrid) 92 92 93 integer ivap, iliq, iice 94 save ivap, iliq, iice 95 96 logical firstcall 97 save firstcall 93 integer, save :: ivap, iliq, iice 94 95 logical, save :: firstcall 98 96 99 97 data firstcall /.true./ … … 333 331 enddo 334 332 335 endif 333 endif 334 336 335 337 336 ! Re-add the albedo effects of CO2 ice if necessary -
trunk/LMDZ.GENERIC/libf/phystd/largescale.F90
r787 r863 138 138 zcond(i) = zcond(i) + zcond_iter 139 139 if (ABS(zcond_iter/alpha).lt.qthreshold) exit 140 zt(i) = zt(i) + zcond_iter*RLVTT/RCPD 140 zt(i) = zt(i) + zcond_iter*RLVTT/RCPD*rneb(i,k) 141 141 call Psat_water(zt(i),pplay(i,k),psat_tmp,zqs(i)) 142 142 End do ! niter -
trunk/LMDZ.GENERIC/libf/phystd/physiq.F90
r862 r863 240 240 241 241 real zdmassmr(ngrid,nlayermx),zdpsrfmr(ngrid) 242 real zdmassmr_col(ngridmx) 242 243 243 244 real zdqdif(ngrid,nlayermx,nq), zdqsdif(ngrid,nq) … … 1242 1243 + dqmoist(1:ngrid,1:nlayermx,igcm_h2o_vap) & 1243 1244 + dqvaplscale(1:ngrid,1:nlayermx) ) 1244 1245 1246 do ig = 1, ngrid 1247 zdmassmr_col(ig)=SUM(zdmassmr(ig,1:nlayermx)) 1248 enddo 1245 1249 1246 1250 call writediagfi(ngrid,"mass_evap","mass gain"," ",3,zdmassmr) 1251 call writediagfi(ngrid,"mass_evap_col","mass gain col"," ",2,zdmassmr_col) 1247 1252 call writediagfi(ngrid,"mass","mass"," ",3,mass) 1248 1253 … … 1702 1707 ! Temporary inclusions for heating diagnostics 1703 1708 ! call writediagfi(ngrid,"zdtdyn","Dyn. heating","T s-1",3,zdtdyn) 1704 !call writediagfi(ngrid,"zdtsw","SW heating","T s-1",3,zdtsw)1705 !call writediagfi(ngrid,"zdtlw","LW heating","T s-1",3,zdtlw)1706 !call writediagfi(ngrid,"dtrad","radiative heating","K s-1",3,dtrad)1709 call writediagfi(ngrid,"zdtsw","SW heating","T s-1",3,zdtsw) 1710 call writediagfi(ngrid,"zdtlw","LW heating","T s-1",3,zdtlw) 1711 call writediagfi(ngrid,"dtrad","radiative heating","K s-1",3,dtrad) 1707 1712 1708 1713 ! debugging -
trunk/LMDZ.GENERIC/libf/phystd/radii_mod.F90
r858 r863 119 119 120 120 !================================================================== 121 ! subroutine h2o_reffrad(ngrid,nq,pq,pt,reffrad,nueffrad)122 121 subroutine h2o_reffrad(ngrid,pq,pt,reffrad,nueffrad) 123 122 !================================================================== … … 131 130 ! 132 131 !================================================================== 133 ! use radinc_h, only: naerkind134 132 use watercommon_h, Only: T_h2O_ice_liq,T_h2O_ice_clouds,rhowater,rhowaterice 135 ! use aerosol_mod, only : iaero_h2o136 ! USE tracer_h, only: igcm_h2o_ice137 133 Implicit none 138 134 … … 143 139 144 140 integer,intent(in) :: ngrid 145 ! intent,integer(in) :: nq 146 147 ! real, intent(in) :: pq(ngrid,nlayermx,nq) !tracer mixing ratios (kg/kg) 141 148 142 real, intent(in) :: pq(ngrid,nlayermx) !water ice mixing ratios (kg/kg) 149 143 real, intent(in) :: pt(ngrid,nlayermx) !temperature (K) 150 144 real, intent(out) :: reffrad(ngrid,nlayermx) !aerosol radii 151 ! real, intent(out) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii152 145 real, intent(out) :: nueffrad(ngrid,nlayermx) ! dispersion 153 ! real, intent(out) :: nueffrad(ngrid,nlayermx,naerkind) ! dispersion154 146 155 147 integer :: ig,l … … 163 155 zfice = 1.0 - (pt(ig,l)-T_h2O_ice_clouds) / (T_h2O_ice_liq-T_h2O_ice_clouds) 164 156 zfice = MIN(MAX(zfice,0.0),1.0) 165 ! reffrad(ig,l,iaero_h2o)= rad_h2o * (1.-zfice) + rad_h2o_ice * zfice166 157 reffrad(ig,l)= rad_h2o * (1.-zfice) + rad_h2o_ice * zfice 167 ! nueffrad(ig,l,iaero_h2o) = coef_chaud * (1.-zfice) + coef_froid * zfice168 158 nueffrad(ig,l) = coef_chaud * (1.-zfice) + coef_froid * zfice 169 159 enddo … … 174 164 zfice = 1.0 - (pt(ig,l)-T_h2O_ice_clouds) / (T_h2O_ice_liq-T_h2O_ice_clouds) 175 165 zfice = MIN(MAX(zfice,0.0),1.0) 176 ! zrad_liq = CBRT( 3*pq(ig,l,igcm_h2o_ice)/(4*Nmix_h2o*pi*rhowater) )177 166 zrad_liq = CBRT( 3*pq(ig,l)/(4*Nmix_h2o*pi*rhowater) ) 178 ! zrad_ice = CBRT( 3*pq(ig,l,igcm_h2o_ice)/(4*Nmix_h2o_ice*pi*rhowaterice) )179 167 zrad_ice = CBRT( 3*pq(ig,l)/(4*Nmix_h2o_ice*pi*rhowaterice) ) 180 ! nueffrad(ig,l,iaero_h2o) = coef_chaud * (1.-zfice) + coef_froid * zfice181 168 nueffrad(ig,l) = coef_chaud * (1.-zfice) + coef_froid * zfice 182 169 zrad = zrad_liq * (1.-zfice) + zrad_ice * zfice 183 ! reffrad(ig,l,iaero_h2o) = min(max(zrad,1.e-6),100.e-6) 184 reffrad(ig,l) = min(max(zrad,1.e-6),100 .e-6)170 171 reffrad(ig,l) = min(max(zrad,1.e-6),1000.e-6) 185 172 enddo 186 173 enddo … … 204 191 !================================================================== 205 192 use watercommon_h, Only: rhowater,rhowaterice 206 ! USE tracer_h207 193 Implicit none 208 194 … … 225 211 else 226 212 reffliq(1:ngrid,1:nlayermx) = CBRT( 3*pql(1:ngrid,1:nlayermx)/(4*Nmix_h2o*pi*rhowater) ) 227 reffliq(1:ngrid,1:nlayermx) = min(max(reffliq(1:ngrid,1:nlayermx),1.e-6),100 .e-6)213 reffliq(1:ngrid,1:nlayermx) = min(max(reffliq(1:ngrid,1:nlayermx),1.e-6),1000.e-6) 228 214 reffice(1:ngrid,1:nlayermx) = CBRT( 3*pql(1:ngrid,1:nlayermx)/(4*Nmix_h2o_ice*pi*rhowaterice) ) 229 reffice(1:ngrid,1:nlayermx) = min(max(reffice(1:ngrid,1:nlayermx),1.e-6),100 .e-6)215 reffice(1:ngrid,1:nlayermx) = min(max(reffice(1:ngrid,1:nlayermx),1.e-6),1000.e-6) 230 216 end if 231 217 … … 247 233 ! 248 234 !================================================================== 249 ! use radinc_h, only: naerkind250 ! use aerosol_mod, only : iaero_co2251 235 USE tracer_h, only:igcm_co2_ice,rho_co2 252 236 Implicit none … … 260 244 261 245 real, intent(in) :: pq(ngrid,nlayermx,nq) !tracer mixing ratios (kg/kg) 262 ! real, intent(out) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 263 real, intent(out) :: reffrad(ngrid,nlayermx) !aerosols radii (K) 246 real, intent(out) :: reffrad(ngrid,nlayermx) !co2 ice particles radii (K) 264 247 265 248 integer :: ig,l … … 270 253 271 254 if (radfixed) then 272 ! reffrad(1:ngrid,1:nlayermx,iaero_co2) = 5.e-5 ! CO2 ice273 255 reffrad(1:ngrid,1:nlayermx) = 5.e-5 ! CO2 ice 274 256 else … … 276 258 do ig=1,ngrid 277 259 zrad = CBRT( 3*pq(ig,l,igcm_co2_ice)/(4*Nmix_co2*pi*rho_co2) ) 278 ! reffrad(ig,l,iaero_co2) = min(max(zrad,1.e-6),100.e-6)279 260 reffrad(ig,l) = min(max(zrad,1.e-6),100.e-6) 280 261 enddo … … 299 280 ! 300 281 !================================================================== 301 ! use radinc_h, only: naerkind 302 ! use aerosol_mod, only : iaero_dust 303 Implicit none 304 305 #include "callkeys.h" 306 #include "dimensions.h" 307 #include "dimphys.h" 308 309 integer,intent(in) :: ngrid 310 311 ! real, intent(out) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 312 real, intent(out) :: reffrad(ngrid,nlayermx) !aerosols radii (K) 282 Implicit none 283 284 #include "callkeys.h" 285 #include "dimensions.h" 286 #include "dimphys.h" 287 288 integer,intent(in) :: ngrid 289 290 real, intent(out) :: reffrad(ngrid,nlayermx) !dust particles radii (K) 313 291 314 ! reffrad(1:ngrid,1:nlayermx,iaero_dust) = 2.e-6 ! dust315 292 reffrad(1:ngrid,1:nlayermx) = 2.e-6 ! dust 316 293 … … 331 308 ! 332 309 !================================================================== 333 ! use radinc_h, only: naerkind 334 ! use aerosol_mod, only : iaero_h2so4 335 Implicit none 336 337 #include "callkeys.h" 338 #include "dimensions.h" 339 #include "dimphys.h" 340 341 integer,intent(in) :: ngrid 342 343 ! real, intent(out) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 344 real, intent(out) :: reffrad(ngrid,nlayermx) !aerosols radii (K) 310 Implicit none 311 312 #include "callkeys.h" 313 #include "dimensions.h" 314 #include "dimphys.h" 315 316 integer,intent(in) :: ngrid 317 318 real, intent(out) :: reffrad(ngrid,nlayermx) !h2so4 particle radii (K) 345 319 346 ! reffrad(1:ngrid,1:nlayermx,iaero_h2so4) = 1.e-6 ! h2so4347 320 reffrad(1:ngrid,1:nlayermx) = 1.e-6 ! h2so4 348 321 -
trunk/LMDZ.GENERIC/libf/phystd/rain.F90
r858 r863 209 209 210 210 if(zt(i,k).gt.Tsat(i,k))then 211 ! treat the case where all liquid water should boil212 zqev=MIN((zt(i,k)-Tsat(i,k))*RCPD*l2c(i,k)/RLVTT ,zrfl(i))211 !! treat the case where all liquid water should boil 212 zqev=MIN((zt(i,k)-Tsat(i,k))*RCPD*l2c(i,k)/RLVTT/ptimestep,zrfl(i)) 213 213 zrfl(i)=MAX(zrfl(i)-zqev,0.) 214 d_q(i,k)=zqev/l2c(i,k) 214 d_q(i,k)=zqev/l2c(i,k)*ptimestep 215 215 d_t(i,k) = - d_q(i,k) * RLVTT/RCPD 216 216 else
Note: See TracChangeset
for help on using the changeset viewer.