subroutine optci(PLEV,TLEV,DTAUI,TAUCUMI, & QXIAER,QSIAER,GIAER,COSBI,WBARI,TAUAERO, & TMID,PMID,TAUGSURF,QVAR) use radinc_h use radcommon_h, only: gasi, tlimit, wrefVAR, Cmk,tgasref,pfgasref,wnoi implicit none !================================================================== ! ! Purpose ! ------- ! Calculates longwave optical constants at each level. For each ! layer and spectral interval in the IR it calculates WBAR, DTAU ! and COSBAR. For each level it calculates TAU. ! ! TAUI(L,LW) is the cumulative optical depth at level L (or alternatively ! at the *bottom* of layer L), LW is the spectral wavelength interval. ! ! TLEV(L) - Temperature at the layer boundary (i.e., level) ! PLEV(L) - Pressure at the layer boundary (i.e., level) ! ! Authors ! ------- ! Adapted from the NASA Ames code by R. Wordsworth (2009) ! !================================================================== #include "comcstfi.h" #include "callkeys.h" real*8 DTAUI(L_NLAYRAD,L_NSPECTI,L_NGAUSS) real*8 DTAUKI(L_LEVELS+1,L_NSPECTI,L_NGAUSS) real*8 TAUI(L_NLEVRAD,L_NSPECTI,L_NGAUSS) real*8 TAUCUMI(L_LEVELS,L_NSPECTI,L_NGAUSS) real*8 PLEV(L_LEVELS) real*8 TLEV(L_LEVELS) real*8 TMID(L_LEVELS), PMID(L_LEVELS) real*8 COSBI(L_NLAYRAD,L_NSPECTI,L_NGAUSS) real*8 WBARI(L_NLAYRAD,L_NSPECTI,L_NGAUSS) ! For aerosols real*8 QXIAER(L_LEVELS+1,L_NSPECTI,NAERKIND) real*8 QSIAER(L_LEVELS+1,L_NSPECTI,NAERKIND) real*8 GIAER(L_LEVELS+1,L_NSPECTI,NAERKIND) real*8 TAUAERO(L_LEVELS+1,NAERKIND) real*8 TAUAEROLK(L_LEVELS+1,L_NSPECTI,NAERKIND) real*8 TAEROS(L_LEVELS,L_NSPECTI,NAERKIND) integer L, NW, NG, K, LK, IAER integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS) real*8 ANS, TAUGAS real*8 DPR(L_LEVELS), U(L_LEVELS) real*8 LCOEF(4), LKCOEF(L_LEVELS,4) real*8 taugsurf(L_NSPECTI,L_NGAUSS-1) real*8 dco2 ! mixing ratio variables real*8 QVAR(L_LEVELS), WRATIO(L_LEVELS) real*8 KCOEF(4) integer NVAR(L_LEVELS) ! temporary variables for multiple aerosol calculation real*8 atemp real*8 btemp(L_NLAYRAD,L_NSPECTI) ! variables for k in units m^-1 real*8 rho, dz !======================================================================= ! Determine the total gas opacity throughout the column, for each ! spectral interval, NW, and each Gauss point, NG. ! write(*,*)'L_LEVELS',L_LEVELS ! write(*,*)'L_NSPECTI',L_NSPECTI DTAUI(:,:,:)=0. DTAUKI(:,:,:)=0. DO NG=1,L_NGAUSS-1 do NW=1,L_NSPECTI TAUGSURF(NW,NG) = 0.0D0 end do end do do K=2,L_LEVELS DPR(k) = PLEV(K)-PLEV(K-1) ! rho = PLEV(K)/(R*TMID(K)) rho = PMID(K)/(R*TMID(K)) dz = -DPR(k)/(g*rho) !print*,'rho=',rho !print*,'dz=',dz U(k) = Cmk*DPR(k) ! only Cmk line in optci.F ! soon to be replaced by m^-1 !!! call tpindex(PMID(K),TMID(K),QVAR(K),pfgasref,tgasref,WREFVAR, & LCOEF,MT(K),MP(K),NVAR(K),WRATIO(K)) do LK=1,4 LKCOEF(K,LK) = LCOEF(LK) end do DO NW=1,L_NSPECTI do iaer=1,naerkind TAEROS(K,NW,IAER) = TAUAERO(K,IAER) * QXIAER(K,NW,IAER) ! write(22,*) 'TB17 Taero IR:',K,NW,IAER,TAEROS(K,NW,IAER) end do END DO end do ! levels do K=2,L_LEVELS do nw=1,L_NSPECTI DCO2 = 0.0 ! continuum absorption (no longer used) do ng=1,L_NGAUSS-1 ! Now compute TAUGAS ! Interpolate between mixing ratios ! WRATIO = 0.0 if the requested amount is equal to, or outside the ! the data range if (L_REFVAR.eq.1)then ! added by RW for special no variable case KCOEF(1) = GASI(MT(K),MP(K),1,NW,NG) KCOEF(2) = GASI(MT(K),MP(K)+1,1,NW,NG) KCOEF(3) = GASI(MT(K)+1,MP(K)+1,1,NW,NG) KCOEF(4) = GASI(MT(K)+1,MP(K),1,NW,NG) else KCOEF(1) = GASI(MT(K),MP(K),NVAR(K),NW,NG) + WRATIO(K)* & (GASI(MT(K),MP(K),NVAR(K)+1,NW,NG) - & GASI(MT(K),MP(K),NVAR(K),NW,NG)) KCOEF(2) = GASI(MT(K),MP(K)+1,NVAR(K),NW,NG)+ WRATIO(K)* & (GASI(MT(K),MP(K)+1,NVAR(K)+1,NW,NG) - & GASI(MT(K),MP(K)+1,NVAR(K),NW,NG)) KCOEF(3)=GASI(MT(K)+1,MP(K)+1,NVAR(K),NW,NG)+WRATIO(K)* & (GASI(MT(K)+1,MP(K)+1,NVAR(K)+1,NW,NG) - & GASI(MT(K)+1,MP(K)+1,NVAR(K),NW,NG)) KCOEF(4) =GASI(MT(K)+1,MP(K),NVAR(K),NW,NG) + WRATIO(K)* & (GASI(MT(K)+1,MP(K),NVAR(K)+1,NW,NG) - & GASI(MT(K)+1,MP(K),NVAR(K),NW,NG)) endif ! Interpolate the gaseous k-coefficients to the requested T,P values ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) + & LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4) TAUGAS = U(k)*ANS TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS DTAUKI(K,nw,ng) = TAUGAS do iaer=1,naerkind DTAUKI(K,nw,ng) = DTAUKI(K,nw,ng) + TAEROS(K,NW,IAER) & + DCO2 ! For Kasting CIA end do end do ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS), ! which holds continuum opacity only NG = L_NGAUSS DTAUKI(K,nw,ng) = 0.0 do iaer=1,naerkind DTAUKI(K,nw,ng) = DTAUKI(K,nw,ng) + TAEROS(K,NW,IAER) & + DCO2 ! For parameterized continuum absorption end do ! a bug was found here!! end do end do !======================================================================= ! Now the full treatment for the layers, where besides the opacity ! we need to calculate the scattering albedo and asymmetry factors DO NW=1,L_NSPECTI DO K=2,L_LEVELS do iaer=1,naerkind TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER)*QSIAER(K,NW,IAER) end do ENDDO ENDDO DO NW=1,L_NSPECTI NG = L_NGAUSS DO L=1,L_NLAYRAD-1 K = 2*L+1 DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) + DTAUKI(K+1,NW,NG)! + 1.e-50 atemp = 0. btemp(L,NW) = 0. do iaer=1,naerkind atemp = atemp + & GIAER(K,NW,IAER) * TAUAEROLK(K,NW,IAER) + & GIAER(K+1,NW,IAER) * TAUAEROLK(K+1,NW,IAER) btemp(L,NW) = btemp(L,NW) + TAUAEROLK(K,NW,IAER) + TAUAEROLK(K+1,NW,IAER) ! * + 1.e-10 end do if(DTAUI(L,NW,NG) .GT. 1.0E-9) then WBARI(L,nw,ng) = btemp(L,NW) / DTAUI(L,NW,NG) else WBARI(L,nw,ng) = 0.0D0 DTAUI(L,NW,NG) = 1.0E-9 endif if(btemp(L,NW) .GT. 0.0) then cosbi(L,NW,NG) = atemp/btemp(L,NW) else cosbi(L,NW,NG) = 0.0D0 end if END DO ! L vertical loop ! Last level L = L_NLAYRAD K = 2*L+1 DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) ! + 1.e-50 btemp(L,NW) = 0 do iaer=1,naerkind btemp(L,NW) = btemp(L,NW) + TAUAEROLK(K,NW,IAER) enddo atemp = 0. if(DTAUI(L,NW,NG) .GT. 1.0D-9) then do iaer=1,naerkind atemp = atemp + GIAER(K,NW,IAER) * TAUAEROLK(K,NW,IAER) end do WBARI(L,nw,ng) = btemp(L,NW) / DTAUI(L,NW,NG) else WBARI(L,nw,ng) = 0.0D0 DTAUI(L,NW,NG) = 1.0D-9 endif if(btemp(L,NW) .GT. 0.0d0) then cosbi(L,NW,NG) = atemp/btemp(L,NW) else cosbi(L,NW,NG) = 0.0D0 end if ! Now the other Gauss points, if needed. DO NG=1,L_NGAUSS-1 IF(TAUGSURF(NW,NG) .gt. TLIMIT) THEN DO L=1,L_NLAYRAD K = 2*L+1 DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)+DTAUKI(K+1,NW,NG)! + 1.e-50 if(DTAUI(L,NW,NG) .GT. 1.0E-9) then WBARI(L,nw,ng) = btemp(L,NW) / DTAUI(L,NW,NG) else WBARI(L,nw,ng) = 0.0D0 DTAUI(L,NW,NG) = 1.0E-9 endif cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS) END DO ! L vertical loop END IF END DO ! NG Gauss loop END DO ! NW spectral loop ! Total extinction optical depths DO NW=1,L_NSPECTI DO NG=1,L_NGAUSS ! full gauss loop TAUI(1,NW,NG)=0.0D0 DO L=1,L_NLAYRAD TAUI(L+1,NW,NG)=TAUI(L,NW,NG)+DTAUI(L,NW,NG) END DO TAUCUMI(1,NW,NG)=0.0D0 DO K=2,L_LEVELS TAUCUMI(K,NW,NG)=TAUCUMI(K-1,NW,NG)+DTAUKI(K,NW,NG) END DO END DO ! end full gauss loop END DO return end subroutine optci