Changeset 1715 for trunk/LMDZ.GENERIC/libf/phystd/optci.F90
- Timestamp:
- Jun 8, 2017, 1:03:06 AM (7 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/libf/phystd/optci.F90
r1397 r1715 32 32 33 33 real*8 DTAUI(L_NLAYRAD,L_NSPECTI,L_NGAUSS) 34 real*8 DTAUKI(L_LEVELS +1,L_NSPECTI,L_NGAUSS)34 real*8 DTAUKI(L_LEVELS,L_NSPECTI,L_NGAUSS) 35 35 real*8 TAUI(L_NLEVRAD,L_NSPECTI,L_NGAUSS) 36 36 real*8 TAUCUMI(L_LEVELS,L_NSPECTI,L_NGAUSS) … … 42 42 43 43 ! for aerosols 44 real*8 QXIAER(L_LEVELS +1,L_NSPECTI,NAERKIND)45 real*8 QSIAER(L_LEVELS +1,L_NSPECTI,NAERKIND)46 real*8 GIAER(L_LEVELS +1,L_NSPECTI,NAERKIND)47 real*8 TAUAERO(L_LEVELS +1,NAERKIND)48 real*8 TAUAEROLK(L_LEVELS +1,L_NSPECTI,NAERKIND)44 real*8 QXIAER(L_LEVELS,L_NSPECTI,NAERKIND) 45 real*8 QSIAER(L_LEVELS,L_NSPECTI,NAERKIND) 46 real*8 GIAER(L_LEVELS,L_NSPECTI,NAERKIND) 47 real*8 TAUAERO(L_LEVELS,NAERKIND) 48 real*8 TAUAEROLK(L_LEVELS,L_NSPECTI,NAERKIND) 49 49 real*8 TAEROS(L_LEVELS,L_NSPECTI,NAERKIND) 50 50 … … 135 135 end do ! levels 136 136 137 137 ! Spectral dependance of aerosol absorption 138 138 do iaer=1,naerkind 139 139 DO NW=1,L_NSPECTI … … 147 147 148 148 do K=2,L_LEVELS 149 150 ! continuum absorption 151 DCONT = 0.0d0 149 150 DAERO=SUM(TAEROS(K,NW,1:naerkind)) ! aerosol absorption 151 152 DCONT = 0.0d0 ! continuum absorption 152 153 153 154 if(continuum.and.(.not.graybody))then … … 225 226 endif 226 227 227 ! aerosol absorption228 DAERO=SUM(TAEROS(K,NW,1:naerkind))229 230 228 do ng=1,L_NGAUSS-1 231 229 … … 286 284 end do 287 285 288 DTAUKI(L_LEVELS+1,1:L_NSPECTI,1:L_NGAUSS)=0.d0289 286 !======================================================================= 290 287 ! Now the full treatment for the layers, where besides the opacity … … 293 290 do iaer=1,naerkind 294 291 DO NW=1,L_NSPECTI 295 DO K=2,L_LEVELS +1296 TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER)*QSIAER(K,NW,IAER) 292 DO K=2,L_LEVELS 293 TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER)*QSIAER(K,NW,IAER) ! effect of scattering albedo 297 294 ENDDO 298 295 ENDDO … … 300 297 301 298 DO NW=1,L_NSPECTI 302 DO L=1,L_NLAYRAD 299 DO L=1,L_NLAYRAD-1 303 300 K = 2*L+1 304 301 btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) + SUM(TAUAEROLK(K+1,NW,1:naerkind)) 305 302 END DO ! L vertical loop 303 304 ! Last level 305 L = L_NLAYRAD 306 K = 2*L+1 307 btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) 308 306 309 END DO ! NW spectral loop 307 310 … … 309 312 DO NW=1,L_NSPECTI 310 313 NG = L_NGAUSS 311 DO L=1,L_NLAYRAD 314 DO L=1,L_NLAYRAD-1 312 315 313 316 K = 2*L+1 … … 334 337 335 338 END DO ! L vertical loop 339 340 ! Last level 341 342 L = L_NLAYRAD 343 K = 2*L+1 344 DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) ! + 1.e-50 345 346 atemp = 0. 347 if(DTAUI(L,NW,NG) .GT. 1.0D-9) then 348 do iaer=1,naerkind 349 atemp = atemp + GIAER(K,NW,IAER) * TAUAEROLK(K,NW,IAER) 350 end do 351 WBARI(L,nw,ng) = btemp(L,nw) / DTAUI(L,NW,NG) 352 else 353 WBARI(L,nw,ng) = 0.0D0 354 DTAUI(L,NW,NG) = 1.0D-9 355 endif 356 357 if(btemp(L,nw) .GT. 0.0d0) then 358 cosbi(L,NW,NG) = atemp/btemp(L,nw) 359 else 360 cosbi(L,NW,NG) = 0.0D0 361 end if 362 336 363 337 364 ! Now the other Gauss points, if needed. … … 340 367 IF(TAUGSURF(NW,NG) .gt. TLIMIT) THEN 341 368 342 DO L=1,L_NLAYRAD 369 DO L=1,L_NLAYRAD-1 343 370 K = 2*L+1 344 371 DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)+DTAUKI(K+1,NW,NG)! + 1.e-50 … … 355 382 cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS) 356 383 END DO ! L vertical loop 384 385 ! Last level 386 L = L_NLAYRAD 387 K = 2*L+1 388 DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)! + 1.e-50 389 390 if(DTAUI(L,NW,NG) .GT. 1.0D-9) then 391 392 WBARI(L,nw,ng) = btemp(L,nw) / DTAUI(L,NW,NG) 393 394 else 395 WBARI(L,nw,ng) = 0.0D0 396 DTAUI(L,NW,NG) = 1.0D-9 397 endif 398 399 cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS) 400 357 401 END IF 358 402
Note: See TracChangeset
for help on using the changeset viewer.