[2050] | 1 | subroutine optci(PQMO,NLAY,PLEV,TLEV,TMID,PMID, & |
---|
[2046] | 2 | DTAUI,TAUCUMI,COSBI,WBARI,TAUGSURF,SEASHAZEFACT) |
---|
[135] | 3 | |
---|
[716] | 4 | use radinc_h |
---|
[2050] | 5 | use radcommon_h, only: gasi,gasi_recomb,tlimit,Cmk,gzlat_ig, & |
---|
[2133] | 6 | tgasref,pfgasref,wnoi,scalep,indi |
---|
[716] | 7 | use gases_h |
---|
[1947] | 8 | use comcstfi_mod, only: r |
---|
[2133] | 9 | use callkeys_mod, only: continuum,graybody,corrk_recombin, & |
---|
[2050] | 10 | callclouds,callmufi,seashaze,uncoupl_optic_haze |
---|
[1897] | 11 | use tracer_h, only : nmicro,nice |
---|
[2133] | 12 | use MMP_OPTICS |
---|
[1897] | 13 | |
---|
[716] | 14 | implicit none |
---|
[135] | 15 | |
---|
[716] | 16 | !================================================================== |
---|
| 17 | ! |
---|
| 18 | ! Purpose |
---|
| 19 | ! ------- |
---|
| 20 | ! Calculates longwave optical constants at each level. For each |
---|
| 21 | ! layer and spectral interval in the IR it calculates WBAR, DTAU |
---|
| 22 | ! and COSBAR. For each level it calculates TAU. |
---|
| 23 | ! |
---|
[1822] | 24 | ! TAUCUMI(L,LW) is the cumulative optical depth at level L (or alternatively |
---|
[716] | 25 | ! at the *bottom* of layer L), LW is the spectral wavelength interval. |
---|
| 26 | ! |
---|
| 27 | ! TLEV(L) - Temperature at the layer boundary (i.e., level) |
---|
| 28 | ! PLEV(L) - Pressure at the layer boundary (i.e., level) |
---|
| 29 | ! |
---|
| 30 | ! Authors |
---|
| 31 | ! ------- |
---|
| 32 | ! Adapted from the NASA Ames code by R. Wordsworth (2009) |
---|
[1822] | 33 | ! Clean and adaptation to Titan by J. Vatant d'Ollone (2016-17) |
---|
[716] | 34 | ! |
---|
| 35 | !================================================================== |
---|
[135] | 36 | |
---|
| 37 | |
---|
[1822] | 38 | !========================================================== |
---|
| 39 | ! Input/Output |
---|
| 40 | !========================================================== |
---|
[2050] | 41 | REAL*8, INTENT(IN) :: PQMO(nlay,nmicro) ! Tracers for microphysics optics (X/m2). |
---|
| 42 | INTEGER, INTENT(IN) :: NLAY ! Number of pressure layers (for pqmo) |
---|
[1822] | 43 | REAL*8, INTENT(IN) :: PLEV(L_LEVELS), TLEV(L_LEVELS) |
---|
| 44 | REAL*8, INTENT(IN) :: TMID(L_LEVELS), PMID(L_LEVELS) |
---|
[2046] | 45 | REAL*8, INTENT(IN) :: SEASHAZEFACT(L_LEVELS) |
---|
[1822] | 46 | |
---|
| 47 | REAL*8, INTENT(OUT) :: DTAUI(L_NLAYRAD,L_NSPECTI,L_NGAUSS) |
---|
| 48 | REAL*8, INTENT(OUT) :: TAUCUMI(L_LEVELS,L_NSPECTI,L_NGAUSS) |
---|
| 49 | REAL*8, INTENT(OUT) :: COSBI(L_NLAYRAD,L_NSPECTI,L_NGAUSS) |
---|
| 50 | REAL*8, INTENT(OUT) :: WBARI(L_NLAYRAD,L_NSPECTI,L_NGAUSS) |
---|
| 51 | REAL*8, INTENT(OUT) :: TAUGSURF(L_NSPECTI,L_NGAUSS-1) |
---|
| 52 | ! ========================================================== |
---|
| 53 | |
---|
[1722] | 54 | real*8 DTAUKI(L_LEVELS,L_NSPECTI,L_NGAUSS) |
---|
[135] | 55 | |
---|
[1648] | 56 | ! Titan customisation |
---|
| 57 | ! J. Vatant d'Ollone (2016) |
---|
[1722] | 58 | real*8 DHAZE_T(L_LEVELS,L_NSPECTI) |
---|
| 59 | real*8 DHAZES_T(L_LEVELS,L_NSPECTI) |
---|
| 60 | real*8 SSA_T(L_LEVELS,L_NSPECTI) |
---|
| 61 | real*8 ASF_T(L_LEVELS,L_NSPECTI) |
---|
[1648] | 62 | ! ========================== |
---|
| 63 | |
---|
[716] | 64 | integer L, NW, NG, K, LK, IAER |
---|
| 65 | integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS) |
---|
| 66 | real*8 ANS, TAUGAS |
---|
| 67 | real*8 DPR(L_LEVELS), U(L_LEVELS) |
---|
| 68 | real*8 LCOEF(4), LKCOEF(L_LEVELS,4) |
---|
[135] | 69 | |
---|
[1788] | 70 | real*8 DCONT |
---|
[716] | 71 | double precision wn_cont, p_cont, p_air, T_cont, dtemp, dtempc |
---|
| 72 | double precision p_cross |
---|
[135] | 73 | |
---|
[716] | 74 | real*8 KCOEF(4) |
---|
[1725] | 75 | |
---|
| 76 | ! temporary variable to reduce memory access time to gasi |
---|
| 77 | real*8 tmpk(2,2) |
---|
[1648] | 78 | |
---|
[716] | 79 | ! temporary variables for multiple aerosol calculation |
---|
[918] | 80 | real*8 atemp |
---|
| 81 | real*8 btemp(L_NLAYRAD,L_NSPECTI) |
---|
[135] | 82 | |
---|
[716] | 83 | ! variables for k in units m^-1 |
---|
[873] | 84 | real*8 dz(L_LEVELS) |
---|
| 85 | !real*8 rho !! see test below |
---|
[135] | 86 | |
---|
[1648] | 87 | integer igas, jgas, ilay |
---|
[253] | 88 | |
---|
[873] | 89 | integer interm |
---|
| 90 | |
---|
[1897] | 91 | real*8 m0as,m3as,m0af,m3af |
---|
| 92 | real*8 ext_s,sca_s,ssa_s,asf_s |
---|
| 93 | real*8 ext_f,sca_f,ssa_f,asf_f |
---|
| 94 | logical,save :: firstcall=.true. |
---|
| 95 | !$OMP THREADPRIVATE(firstcall) |
---|
| 96 | |
---|
[873] | 97 | !! AS: to save time in computing continuum (see bilinearbig) |
---|
| 98 | IF (.not.ALLOCATED(indi)) THEN |
---|
[878] | 99 | ALLOCATE(indi(L_NSPECTI,ngasmx,ngasmx)) |
---|
[873] | 100 | indi = -9999 ! this initial value means "to be calculated" |
---|
| 101 | ENDIF |
---|
| 102 | |
---|
[1792] | 103 | ! Some initialisation beacause there's a pb with disr_haze at the limits (nw=1) |
---|
| 104 | ! I should check this - For now we set vars to zero : better than nans - JVO 2017 |
---|
| 105 | |
---|
| 106 | dhaze_t(:,:) = 0. |
---|
| 107 | ssa_t(:,:) = 0. |
---|
| 108 | asf_t(:,:) = 0. |
---|
| 109 | |
---|
[716] | 110 | !======================================================================= |
---|
| 111 | ! Determine the total gas opacity throughout the column, for each |
---|
| 112 | ! spectral interval, NW, and each Gauss point, NG. |
---|
[135] | 113 | |
---|
[716] | 114 | taugsurf(:,:) = 0.0 |
---|
| 115 | dpr(:) = 0.0 |
---|
| 116 | lkcoef(:,:) = 0.0 |
---|
[135] | 117 | |
---|
[716] | 118 | do K=2,L_LEVELS |
---|
[1947] | 119 | |
---|
| 120 | ilay = k / 2 ! int. arithmetic => gives the gcm layer index |
---|
| 121 | |
---|
[716] | 122 | DPR(k) = PLEV(K)-PLEV(K-1) |
---|
[135] | 123 | |
---|
[716] | 124 | ! if we have continuum opacities, we need dz |
---|
[1947] | 125 | dz(k) = dpr(k)*R*TMID(K)/(gzlat_ig(ilay)*PMID(K)) |
---|
| 126 | U(k) = Cmk(ilay)*DPR(k) ! only Cmk line in optci.F |
---|
[1648] | 127 | |
---|
| 128 | call tpindex(PMID(K),TMID(K),pfgasref,tgasref,LCOEF,MT(K),MP(K)) |
---|
[135] | 129 | |
---|
[716] | 130 | do LK=1,4 |
---|
| 131 | LKCOEF(K,LK) = LCOEF(LK) |
---|
| 132 | end do |
---|
[918] | 133 | end do ! levels |
---|
[253] | 134 | |
---|
[918] | 135 | do NW=1,L_NSPECTI |
---|
[135] | 136 | |
---|
[918] | 137 | do K=2,L_LEVELS |
---|
[1722] | 138 | |
---|
[1648] | 139 | ilay = k / 2 ! int. arithmetic => gives the gcm layer index |
---|
| 140 | |
---|
[1897] | 141 | ! Optical coupling of YAMMS is plugged but inactivated for now |
---|
| 142 | ! as long as the microphysics only isn't fully debugged -- JVO 01/18 |
---|
| 143 | IF (callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN |
---|
[2050] | 144 | m0as = pqmo(ilay,1) |
---|
| 145 | m3as = pqmo(ilay,2) |
---|
| 146 | m0af = pqmo(ilay,3) |
---|
| 147 | m3af = pqmo(ilay,4) |
---|
[1722] | 148 | |
---|
[1897] | 149 | IF (.NOT.mmp_sph_optics_ir(m0as,m3as,nw,ext_s,sca_s,ssa_s,asf_s)) & |
---|
| 150 | CALL abort_gcm("optcv", "Fatal error in mmp_sph_optics_ir", 12) |
---|
| 151 | IF (.NOT.mmp_fra_optics_ir(m0af,m3af,nw,ext_f,sca_f,ssa_f,asf_f)) & |
---|
| 152 | CALL abort_gcm("optcv", "Fatal error in mmp_fra_optics_ir", 12) |
---|
| 153 | dhaze_T(k,nw) = ext_s+ext_f |
---|
| 154 | SSA_T(k,nw) = (sca_s+sca_f)/dhaze_T(k,nw) |
---|
| 155 | ASF_T(k,nw) = (asf_s*sca_s + asf_f*sca_f) /(sca_s+sca_f) |
---|
| 156 | IF (callclouds.and.firstcall) & |
---|
| 157 | WRITE(*,*) 'WARNING: In optci, optical properties & |
---|
| 158 | &calculations are not implemented yet' |
---|
| 159 | ELSE |
---|
| 160 | ! Call fixed vertical haze profile of extinction - same for all columns |
---|
| 161 | call disr_haze(dz(k),plev(k),wnoi(nw),dhaze_T(k,nw),SSA_T(k,nw),ASF_T(k,nw)) |
---|
[2046] | 162 | if (seashaze) dhaze_T(k,nw) = dhaze_T(k,nw)*seashazefact(k) |
---|
[1897] | 163 | ENDIF |
---|
| 164 | |
---|
[1722] | 165 | DCONT = 0.0d0 ! continuum absorption |
---|
| 166 | |
---|
[873] | 167 | if(continuum.and.(.not.graybody))then |
---|
[716] | 168 | ! include continua if necessary |
---|
| 169 | wn_cont = dble(wnoi(nw)) |
---|
| 170 | T_cont = dble(TMID(k)) |
---|
| 171 | do igas=1,ngasmx |
---|
[135] | 172 | |
---|
[1648] | 173 | p_cont = dble(PMID(k)*scalep*gfrac(igas,ilay)) |
---|
[253] | 174 | |
---|
[961] | 175 | dtemp=0.0d0 |
---|
[716] | 176 | if(igas.eq.igas_N2)then |
---|
[305] | 177 | |
---|
[878] | 178 | interm = indi(nw,igas,igas) |
---|
| 179 | call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
| 180 | indi(nw,igas,igas) = interm |
---|
[253] | 181 | |
---|
[716] | 182 | elseif(igas.eq.igas_H2)then |
---|
[253] | 183 | |
---|
[716] | 184 | ! first do self-induced absorption |
---|
[878] | 185 | interm = indi(nw,igas,igas) |
---|
[873] | 186 | call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
[878] | 187 | indi(nw,igas,igas) = interm |
---|
[253] | 188 | |
---|
[716] | 189 | ! then cross-interactions with other gases |
---|
| 190 | do jgas=1,ngasmx |
---|
[1648] | 191 | p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay)) |
---|
[961] | 192 | dtempc = 0.0d0 |
---|
[716] | 193 | if(jgas.eq.igas_N2)then |
---|
[878] | 194 | interm = indi(nw,igas,jgas) |
---|
| 195 | call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) |
---|
| 196 | indi(nw,igas,jgas) = interm |
---|
[716] | 197 | endif |
---|
| 198 | dtemp = dtemp + dtempc |
---|
| 199 | enddo |
---|
[135] | 200 | |
---|
[1648] | 201 | elseif(igas.eq.igas_CH4)then |
---|
| 202 | |
---|
| 203 | ! first do self-induced absorption |
---|
| 204 | interm = indi(nw,igas,igas) |
---|
| 205 | call interpolateCH4CH4(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
| 206 | indi(nw,igas,igas) = interm |
---|
| 207 | |
---|
| 208 | ! then cross-interactions with other gases |
---|
| 209 | do jgas=1,ngasmx |
---|
| 210 | p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay)) |
---|
| 211 | dtempc = 0.0d0 |
---|
| 212 | if(jgas.eq.igas_N2)then |
---|
| 213 | interm = indi(nw,igas,jgas) |
---|
| 214 | call interpolateN2CH4(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) |
---|
| 215 | indi(nw,igas,jgas) = interm |
---|
| 216 | endif |
---|
| 217 | dtemp = dtemp + dtempc |
---|
| 218 | enddo |
---|
| 219 | |
---|
[716] | 220 | endif |
---|
[135] | 221 | |
---|
[716] | 222 | DCONT = DCONT + dtemp |
---|
[135] | 223 | |
---|
[716] | 224 | enddo |
---|
[135] | 225 | |
---|
[716] | 226 | DCONT = DCONT*dz(k) |
---|
[135] | 227 | |
---|
[716] | 228 | endif |
---|
[135] | 229 | |
---|
[716] | 230 | do ng=1,L_NGAUSS-1 |
---|
[135] | 231 | |
---|
[716] | 232 | ! Now compute TAUGAS |
---|
[135] | 233 | |
---|
[1725] | 234 | ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically |
---|
| 235 | ! the execution time of optci/v -> ~ factor 2 on the whole radiative |
---|
| 236 | ! transfer on the tested simulations ! |
---|
[253] | 237 | |
---|
[2050] | 238 | if (corrk_recombin) then |
---|
| 239 | tmpk = GASI_RECOMB(MT(K):MT(K)+1,MP(K):MP(K)+1,NW,NG) |
---|
| 240 | else |
---|
| 241 | tmpk = GASI(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG) |
---|
| 242 | endif |
---|
[1725] | 243 | |
---|
| 244 | KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASI(MT(K),MP(K),1,NW,NG) |
---|
| 245 | KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASI(MT(K),MP(K)+1,1,NW,NG) |
---|
| 246 | KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASI(MT(K)+1,MP(K)+1,1,NW,NG) |
---|
| 247 | KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASI(MT(K)+1,MP(K),1,NW,NG) |
---|
| 248 | |
---|
| 249 | |
---|
[716] | 250 | ! Interpolate the gaseous k-coefficients to the requested T,P values |
---|
[135] | 251 | |
---|
[716] | 252 | ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) + & |
---|
| 253 | LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4) |
---|
[135] | 254 | |
---|
[716] | 255 | TAUGAS = U(k)*ANS |
---|
[135] | 256 | |
---|
[716] | 257 | TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT |
---|
[918] | 258 | DTAUKI(K,nw,ng) = TAUGAS & |
---|
| 259 | + DCONT & ! For parameterized continuum absorption |
---|
[1648] | 260 | + DHAZE_T(K,NW) ! For Titan haze |
---|
[135] | 261 | |
---|
[716] | 262 | end do |
---|
[135] | 263 | |
---|
[716] | 264 | ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS), |
---|
| 265 | ! which holds continuum opacity only |
---|
[135] | 266 | |
---|
[716] | 267 | NG = L_NGAUSS |
---|
[918] | 268 | DTAUKI(K,nw,ng) = 0.d0 & |
---|
| 269 | + DCONT & ! For parameterized continuum absorption |
---|
[1648] | 270 | + DHAZE_T(K,NW) ! For Titan Haze |
---|
[135] | 271 | |
---|
[716] | 272 | end do |
---|
| 273 | end do |
---|
[135] | 274 | |
---|
[716] | 275 | !======================================================================= |
---|
| 276 | ! Now the full treatment for the layers, where besides the opacity |
---|
| 277 | ! we need to calculate the scattering albedo and asymmetry factors |
---|
[1648] | 278 | ! ====================================================================== |
---|
[135] | 279 | |
---|
[1648] | 280 | ! Haze scattering |
---|
[918] | 281 | DO NW=1,L_NSPECTI |
---|
[1722] | 282 | DO K=2,L_LEVELS |
---|
[1648] | 283 | DHAZES_T(K,NW) = DHAZE_T(K,NW) * SSA_T(K,NW) |
---|
| 284 | ENDDO |
---|
| 285 | ENDDO |
---|
| 286 | |
---|
| 287 | DO NW=1,L_NSPECTI |
---|
| 288 | DO L=1,L_NLAYRAD-1 |
---|
[918] | 289 | K = 2*L+1 |
---|
[1788] | 290 | btemp(L,NW) = DHAZES_T(K,NW) + DHAZES_T(K+1,NW) |
---|
[918] | 291 | END DO ! L vertical loop |
---|
[1648] | 292 | |
---|
[1722] | 293 | ! Last level |
---|
| 294 | L = L_NLAYRAD |
---|
| 295 | K = 2*L+1 |
---|
[1788] | 296 | btemp(L,NW) = DHAZES_T(K,NW) |
---|
[1648] | 297 | |
---|
[918] | 298 | END DO ! NW spectral loop |
---|
| 299 | |
---|
[135] | 300 | |
---|
[716] | 301 | DO NW=1,L_NSPECTI |
---|
| 302 | NG = L_NGAUSS |
---|
[1648] | 303 | DO L=1,L_NLAYRAD-1 |
---|
[135] | 304 | |
---|
[716] | 305 | K = 2*L+1 |
---|
| 306 | DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) + DTAUKI(K+1,NW,NG)! + 1.e-50 |
---|
[135] | 307 | |
---|
[716] | 308 | atemp = 0. |
---|
[961] | 309 | if(DTAUI(L,NW,NG) .GT. 1.0D-9) then |
---|
[1722] | 310 | atemp = atemp + & |
---|
| 311 | ASF_T(K,NW)*DHAZES_T(K,NW) + & |
---|
| 312 | ASF_T(K+1,NW)*DHAZES_T(K+1,NW) |
---|
| 313 | |
---|
[918] | 314 | WBARI(L,nw,ng) = btemp(L,nw) / DTAUI(L,NW,NG) |
---|
[716] | 315 | else |
---|
| 316 | WBARI(L,nw,ng) = 0.0D0 |
---|
[961] | 317 | DTAUI(L,NW,NG) = 1.0D-9 |
---|
[716] | 318 | endif |
---|
[135] | 319 | |
---|
[961] | 320 | if(btemp(L,nw) .GT. 0.0d0) then |
---|
[918] | 321 | cosbi(L,NW,NG) = atemp/btemp(L,nw) |
---|
[716] | 322 | else |
---|
| 323 | cosbi(L,NW,NG) = 0.0D0 |
---|
| 324 | end if |
---|
[135] | 325 | |
---|
[716] | 326 | END DO ! L vertical loop |
---|
[1648] | 327 | |
---|
[1722] | 328 | ! Last level |
---|
[1648] | 329 | |
---|
[1722] | 330 | L = L_NLAYRAD |
---|
| 331 | K = 2*L+1 |
---|
| 332 | DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) ! + 1.e-50 |
---|
[1648] | 333 | |
---|
[1722] | 334 | atemp = 0. |
---|
| 335 | if(DTAUI(L,NW,NG) .GT. 1.0D-9) then |
---|
| 336 | atemp = atemp + ASF_T(K,NW)*DHAZES_T(K,NW) |
---|
| 337 | WBARI(L,nw,ng) = btemp(L,nw) / DTAUI(L,NW,NG) |
---|
| 338 | else |
---|
| 339 | WBARI(L,nw,ng) = 0.0D0 |
---|
| 340 | DTAUI(L,NW,NG) = 1.0D-9 |
---|
| 341 | endif |
---|
[1648] | 342 | |
---|
[1722] | 343 | if(btemp(L,nw) .GT. 0.0d0) then |
---|
| 344 | cosbi(L,NW,NG) = atemp/btemp(L,nw) |
---|
| 345 | else |
---|
| 346 | cosbi(L,NW,NG) = 0.0D0 |
---|
| 347 | end if |
---|
| 348 | |
---|
| 349 | |
---|
[716] | 350 | ! Now the other Gauss points, if needed. |
---|
[135] | 351 | |
---|
[716] | 352 | DO NG=1,L_NGAUSS-1 |
---|
| 353 | IF(TAUGSURF(NW,NG) .gt. TLIMIT) THEN |
---|
[135] | 354 | |
---|
[1648] | 355 | DO L=1,L_NLAYRAD-1 |
---|
[716] | 356 | K = 2*L+1 |
---|
| 357 | DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)+DTAUKI(K+1,NW,NG)! + 1.e-50 |
---|
| 358 | |
---|
[961] | 359 | if(DTAUI(L,NW,NG) .GT. 1.0D-9) then |
---|
[1722] | 360 | |
---|
| 361 | WBARI(L,nw,ng) = btemp(L,nw) / DTAUI(L,NW,NG) |
---|
| 362 | |
---|
[716] | 363 | else |
---|
| 364 | WBARI(L,nw,ng) = 0.0D0 |
---|
[961] | 365 | DTAUI(L,NW,NG) = 1.0D-9 |
---|
[716] | 366 | endif |
---|
| 367 | |
---|
| 368 | cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS) |
---|
| 369 | END DO ! L vertical loop |
---|
[1648] | 370 | |
---|
[1722] | 371 | ! Last level |
---|
| 372 | L = L_NLAYRAD |
---|
| 373 | K = 2*L+1 |
---|
| 374 | DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)! + 1.e-50 |
---|
[1648] | 375 | |
---|
[1722] | 376 | if(DTAUI(L,NW,NG) .GT. 1.0D-9) then |
---|
| 377 | |
---|
[1648] | 378 | WBARI(L,nw,ng) = btemp(L,nw) / DTAUI(L,NW,NG) |
---|
[1722] | 379 | |
---|
| 380 | else |
---|
| 381 | WBARI(L,nw,ng) = 0.0D0 |
---|
| 382 | DTAUI(L,NW,NG) = 1.0D-9 |
---|
| 383 | endif |
---|
| 384 | |
---|
| 385 | cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS) |
---|
[1648] | 386 | |
---|
[716] | 387 | END IF |
---|
| 388 | |
---|
| 389 | END DO ! NG Gauss loop |
---|
| 390 | END DO ! NW spectral loop |
---|
| 391 | |
---|
| 392 | ! Total extinction optical depths |
---|
| 393 | |
---|
[918] | 394 | DO NG=1,L_NGAUSS ! full gauss loop |
---|
| 395 | DO NW=1,L_NSPECTI |
---|
[716] | 396 | TAUCUMI(1,NW,NG)=0.0D0 |
---|
| 397 | DO K=2,L_LEVELS |
---|
| 398 | TAUCUMI(K,NW,NG)=TAUCUMI(K-1,NW,NG)+DTAUKI(K,NW,NG) |
---|
| 399 | END DO |
---|
| 400 | END DO ! end full gauss loop |
---|
| 401 | END DO |
---|
| 402 | |
---|
| 403 | ! be aware when comparing with textbook results |
---|
| 404 | ! (e.g. Pierrehumbert p. 218) that |
---|
| 405 | ! taucumi does not take the <cos theta>=0.5 factor into |
---|
| 406 | ! account. It is the optical depth for a vertically |
---|
| 407 | ! ascending ray with angle theta = 0. |
---|
| 408 | |
---|
[1897] | 409 | if(firstcall) firstcall = .false. |
---|
| 410 | |
---|
[716] | 411 | return |
---|
| 412 | |
---|
| 413 | |
---|
| 414 | end subroutine optci |
---|
| 415 | |
---|
| 416 | |
---|
| 417 | |
---|