[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 |
---|
[2242] | 8 | use datafile_mod, only: haze_opt_file |
---|
[1947] | 9 | use comcstfi_mod, only: r |
---|
[2133] | 10 | use callkeys_mod, only: continuum,graybody,corrk_recombin, & |
---|
[2050] | 11 | callclouds,callmufi,seashaze,uncoupl_optic_haze |
---|
[1897] | 12 | use tracer_h, only : nmicro,nice |
---|
| 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) |
---|
[135] | 85 | |
---|
[1648] | 86 | integer igas, jgas, ilay |
---|
[253] | 87 | |
---|
[873] | 88 | integer interm |
---|
| 89 | |
---|
[2242] | 90 | ! Variables for haze optics |
---|
| 91 | character(len=200) file_path |
---|
| 92 | logical file_ok |
---|
| 93 | integer dumch |
---|
| 94 | real*8 dumwvl |
---|
| 95 | |
---|
| 96 | real*8 m3as,m3af |
---|
| 97 | real*8 dtauaer_s,dtauaer_f |
---|
| 98 | real*8,save :: rhoaer_s(L_NSPECTI),ssa_s(L_NSPECTI),asf_s(L_NSPECTI) |
---|
| 99 | real*8,save :: rhoaer_f(L_NSPECTI),ssa_f(L_NSPECTI),asf_f(L_NSPECTI) |
---|
| 100 | !$OMP THREADPRIVATE(rhoaer_s,rhoaer_f,ssa_s,ssa_f,asf_s,asf_f) |
---|
| 101 | |
---|
[1897] | 102 | logical,save :: firstcall=.true. |
---|
[2242] | 103 | !$OMP THREADPRIVATE(firstcall) |
---|
[1897] | 104 | |
---|
[2242] | 105 | |
---|
[873] | 106 | !! AS: to save time in computing continuum (see bilinearbig) |
---|
| 107 | IF (.not.ALLOCATED(indi)) THEN |
---|
[878] | 108 | ALLOCATE(indi(L_NSPECTI,ngasmx,ngasmx)) |
---|
[873] | 109 | indi = -9999 ! this initial value means "to be calculated" |
---|
| 110 | ENDIF |
---|
| 111 | |
---|
[2242] | 112 | ! Some initialisation because there's a pb with disr_haze at the limits (nw=1) |
---|
[1792] | 113 | ! I should check this - For now we set vars to zero : better than nans - JVO 2017 |
---|
[2242] | 114 | DHAZE_T(:,:) = 0.0 |
---|
| 115 | SSA_T(:,:) = 0.0 |
---|
| 116 | ASF_T(:,:) = 0.0 |
---|
[1792] | 117 | |
---|
[2242] | 118 | ! Load tabulated haze optical properties if needed. |
---|
| 119 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
| 120 | IF (firstcall .AND. callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN |
---|
| 121 | OPEN(12,file=TRIM(haze_opt_file),form='formatted') ! The file has been inquired in physiq_mod firstcall |
---|
| 122 | READ(12,*) ! dummy header |
---|
| 123 | DO NW=1,L_NSPECTI |
---|
| 124 | READ(12,*) dumch, dumwvl, rhoaer_f(nw), ssa_f(nw), asf_f(nw), rhoaer_s(nw), ssa_s(nw), asf_s(nw) |
---|
| 125 | ENDDO |
---|
| 126 | CLOSE(12) |
---|
| 127 | ENDIF |
---|
| 128 | |
---|
[716] | 129 | !======================================================================= |
---|
| 130 | ! Determine the total gas opacity throughout the column, for each |
---|
| 131 | ! spectral interval, NW, and each Gauss point, NG. |
---|
[135] | 132 | |
---|
[716] | 133 | taugsurf(:,:) = 0.0 |
---|
| 134 | dpr(:) = 0.0 |
---|
| 135 | lkcoef(:,:) = 0.0 |
---|
[135] | 136 | |
---|
[716] | 137 | do K=2,L_LEVELS |
---|
[1947] | 138 | |
---|
[2242] | 139 | ilay = L_NLAYRAD+1 - k/2 ! int. arithmetic => gives the gcm layer index (reversed) |
---|
[1947] | 140 | |
---|
[716] | 141 | DPR(k) = PLEV(K)-PLEV(K-1) |
---|
[135] | 142 | |
---|
[716] | 143 | ! if we have continuum opacities, we need dz |
---|
[1947] | 144 | dz(k) = dpr(k)*R*TMID(K)/(gzlat_ig(ilay)*PMID(K)) |
---|
| 145 | U(k) = Cmk(ilay)*DPR(k) ! only Cmk line in optci.F |
---|
[1648] | 146 | |
---|
| 147 | call tpindex(PMID(K),TMID(K),pfgasref,tgasref,LCOEF,MT(K),MP(K)) |
---|
[135] | 148 | |
---|
[716] | 149 | do LK=1,4 |
---|
| 150 | LKCOEF(K,LK) = LCOEF(LK) |
---|
| 151 | end do |
---|
[918] | 152 | end do ! levels |
---|
[253] | 153 | |
---|
[918] | 154 | do NW=1,L_NSPECTI |
---|
[135] | 155 | |
---|
[918] | 156 | do K=2,L_LEVELS |
---|
[1722] | 157 | |
---|
[2242] | 158 | ilay = L_NLAYRAD+1 - k/2 ! int. arithmetic => gives the gcm layer index (reversed) |
---|
[1648] | 159 | |
---|
[1897] | 160 | IF (callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN |
---|
[2242] | 161 | m3as = pqmo(ilay,2) / 2.0 |
---|
| 162 | m3af = pqmo(ilay,4) / 2.0 |
---|
| 163 | |
---|
| 164 | IF ( ilay .lt. 18 ) THEN |
---|
| 165 | m3as = pqmo(18,2) / 2.0 *dz(k)/dz(18) |
---|
| 166 | m3af = pqmo(18,4) / 2.0 *dz(k)/dz(18) |
---|
| 167 | ENDIF |
---|
[1722] | 168 | |
---|
[2242] | 169 | dtauaer_s = m3as*rhoaer_s(nw) |
---|
| 170 | dtauaer_f = m3af*rhoaer_f(nw) |
---|
| 171 | DHAZE_T(k,nw) = dtauaer_s + dtauaer_f |
---|
| 172 | |
---|
| 173 | IF ( dtauaer_s + dtauaer_f .GT. 1.D-30 ) THEN |
---|
| 174 | SSA_T(k,nw) = ( dtauaer_s*ssa_s(nw) + dtauaer_f*ssa_f(nw) ) / ( dtauaer_s+dtauaer_f ) |
---|
| 175 | ASF_T(k,nw) = ( dtauaer_s*ssa_s(nw)*asf_s(nw) + dtauaer_f*ssa_f(nw)*asf_f(nw) ) & |
---|
| 176 | / ( ssa_s(nw)*dtauaer_s + ssa_f(nw)*dtauaer_f ) |
---|
| 177 | ELSE |
---|
| 178 | DHAZE_T(k,nw) = 0.D0 |
---|
| 179 | SSA_T(k,nw) = 1.0 |
---|
| 180 | ASF_T(k,nw) = 1.0 |
---|
| 181 | ENDIF |
---|
| 182 | |
---|
[1897] | 183 | IF (callclouds.and.firstcall) & |
---|
| 184 | WRITE(*,*) 'WARNING: In optci, optical properties & |
---|
| 185 | &calculations are not implemented yet' |
---|
| 186 | ELSE |
---|
| 187 | ! Call fixed vertical haze profile of extinction - same for all columns |
---|
[2242] | 188 | call disr_haze(dz(k),plev(k),wnoi(nw),DHAZE_T(k,nw),SSA_T(k,nw),ASF_T(k,nw)) |
---|
| 189 | if (seashaze) DHAZE_T(k,nw) = DHAZE_T(k,nw)*seashazefact(k) |
---|
[1897] | 190 | ENDIF |
---|
| 191 | |
---|
[1722] | 192 | DCONT = 0.0d0 ! continuum absorption |
---|
| 193 | |
---|
[873] | 194 | if(continuum.and.(.not.graybody))then |
---|
[716] | 195 | ! include continua if necessary |
---|
| 196 | wn_cont = dble(wnoi(nw)) |
---|
| 197 | T_cont = dble(TMID(k)) |
---|
| 198 | do igas=1,ngasmx |
---|
[135] | 199 | |
---|
[1648] | 200 | p_cont = dble(PMID(k)*scalep*gfrac(igas,ilay)) |
---|
[253] | 201 | |
---|
[961] | 202 | dtemp=0.0d0 |
---|
[716] | 203 | if(igas.eq.igas_N2)then |
---|
[305] | 204 | |
---|
[878] | 205 | interm = indi(nw,igas,igas) |
---|
| 206 | call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
| 207 | indi(nw,igas,igas) = interm |
---|
[253] | 208 | |
---|
[716] | 209 | elseif(igas.eq.igas_H2)then |
---|
[253] | 210 | |
---|
[716] | 211 | ! first do self-induced absorption |
---|
[878] | 212 | interm = indi(nw,igas,igas) |
---|
[873] | 213 | call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
[878] | 214 | indi(nw,igas,igas) = interm |
---|
[253] | 215 | |
---|
[716] | 216 | ! then cross-interactions with other gases |
---|
| 217 | do jgas=1,ngasmx |
---|
[1648] | 218 | p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay)) |
---|
[961] | 219 | dtempc = 0.0d0 |
---|
[716] | 220 | if(jgas.eq.igas_N2)then |
---|
[878] | 221 | interm = indi(nw,igas,jgas) |
---|
| 222 | call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) |
---|
| 223 | indi(nw,igas,jgas) = interm |
---|
[716] | 224 | endif |
---|
| 225 | dtemp = dtemp + dtempc |
---|
| 226 | enddo |
---|
[135] | 227 | |
---|
[1648] | 228 | elseif(igas.eq.igas_CH4)then |
---|
| 229 | |
---|
| 230 | ! first do self-induced absorption |
---|
| 231 | interm = indi(nw,igas,igas) |
---|
| 232 | call interpolateCH4CH4(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
| 233 | indi(nw,igas,igas) = interm |
---|
| 234 | |
---|
| 235 | ! then cross-interactions with other gases |
---|
| 236 | do jgas=1,ngasmx |
---|
| 237 | p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay)) |
---|
| 238 | dtempc = 0.0d0 |
---|
| 239 | if(jgas.eq.igas_N2)then |
---|
| 240 | interm = indi(nw,igas,jgas) |
---|
| 241 | call interpolateN2CH4(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) |
---|
| 242 | indi(nw,igas,jgas) = interm |
---|
| 243 | endif |
---|
| 244 | dtemp = dtemp + dtempc |
---|
| 245 | enddo |
---|
| 246 | |
---|
[716] | 247 | endif |
---|
[135] | 248 | |
---|
[716] | 249 | DCONT = DCONT + dtemp |
---|
[135] | 250 | |
---|
[716] | 251 | enddo |
---|
[135] | 252 | |
---|
[716] | 253 | DCONT = DCONT*dz(k) |
---|
[135] | 254 | |
---|
[716] | 255 | endif |
---|
[135] | 256 | |
---|
[716] | 257 | do ng=1,L_NGAUSS-1 |
---|
[135] | 258 | |
---|
[716] | 259 | ! Now compute TAUGAS |
---|
[135] | 260 | |
---|
[1725] | 261 | ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically |
---|
| 262 | ! the execution time of optci/v -> ~ factor 2 on the whole radiative |
---|
| 263 | ! transfer on the tested simulations ! |
---|
[253] | 264 | |
---|
[2050] | 265 | if (corrk_recombin) then |
---|
| 266 | tmpk = GASI_RECOMB(MT(K):MT(K)+1,MP(K):MP(K)+1,NW,NG) |
---|
| 267 | else |
---|
| 268 | tmpk = GASI(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG) |
---|
| 269 | endif |
---|
[1725] | 270 | |
---|
| 271 | KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASI(MT(K),MP(K),1,NW,NG) |
---|
| 272 | KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASI(MT(K),MP(K)+1,1,NW,NG) |
---|
| 273 | KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASI(MT(K)+1,MP(K)+1,1,NW,NG) |
---|
| 274 | KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASI(MT(K)+1,MP(K),1,NW,NG) |
---|
| 275 | |
---|
| 276 | |
---|
[716] | 277 | ! Interpolate the gaseous k-coefficients to the requested T,P values |
---|
[135] | 278 | |
---|
[716] | 279 | ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) + & |
---|
| 280 | LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4) |
---|
[135] | 281 | |
---|
[716] | 282 | TAUGAS = U(k)*ANS |
---|
[135] | 283 | |
---|
[716] | 284 | TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT |
---|
[918] | 285 | DTAUKI(K,nw,ng) = TAUGAS & |
---|
| 286 | + DCONT & ! For parameterized continuum absorption |
---|
[1648] | 287 | + DHAZE_T(K,NW) ! For Titan haze |
---|
[135] | 288 | |
---|
[716] | 289 | end do |
---|
[135] | 290 | |
---|
[716] | 291 | ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS), |
---|
| 292 | ! which holds continuum opacity only |
---|
[135] | 293 | |
---|
[716] | 294 | NG = L_NGAUSS |
---|
[918] | 295 | DTAUKI(K,nw,ng) = 0.d0 & |
---|
| 296 | + DCONT & ! For parameterized continuum absorption |
---|
[1648] | 297 | + DHAZE_T(K,NW) ! For Titan Haze |
---|
[135] | 298 | |
---|
[716] | 299 | end do |
---|
| 300 | end do |
---|
[135] | 301 | |
---|
[716] | 302 | !======================================================================= |
---|
| 303 | ! Now the full treatment for the layers, where besides the opacity |
---|
| 304 | ! we need to calculate the scattering albedo and asymmetry factors |
---|
[1648] | 305 | ! ====================================================================== |
---|
[135] | 306 | |
---|
[1648] | 307 | ! Haze scattering |
---|
[918] | 308 | DO NW=1,L_NSPECTI |
---|
[1722] | 309 | DO K=2,L_LEVELS |
---|
[1648] | 310 | DHAZES_T(K,NW) = DHAZE_T(K,NW) * SSA_T(K,NW) |
---|
| 311 | ENDDO |
---|
| 312 | ENDDO |
---|
| 313 | |
---|
| 314 | DO NW=1,L_NSPECTI |
---|
| 315 | DO L=1,L_NLAYRAD-1 |
---|
[918] | 316 | K = 2*L+1 |
---|
[1788] | 317 | btemp(L,NW) = DHAZES_T(K,NW) + DHAZES_T(K+1,NW) |
---|
[918] | 318 | END DO ! L vertical loop |
---|
[1648] | 319 | |
---|
[1722] | 320 | ! Last level |
---|
| 321 | L = L_NLAYRAD |
---|
| 322 | K = 2*L+1 |
---|
[1788] | 323 | btemp(L,NW) = DHAZES_T(K,NW) |
---|
[1648] | 324 | |
---|
[918] | 325 | END DO ! NW spectral loop |
---|
| 326 | |
---|
[135] | 327 | |
---|
[716] | 328 | DO NW=1,L_NSPECTI |
---|
| 329 | NG = L_NGAUSS |
---|
[1648] | 330 | DO L=1,L_NLAYRAD-1 |
---|
[135] | 331 | |
---|
[716] | 332 | K = 2*L+1 |
---|
| 333 | DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) + DTAUKI(K+1,NW,NG)! + 1.e-50 |
---|
[135] | 334 | |
---|
[716] | 335 | atemp = 0. |
---|
[961] | 336 | if(DTAUI(L,NW,NG) .GT. 1.0D-9) then |
---|
[1722] | 337 | atemp = atemp + & |
---|
| 338 | ASF_T(K,NW)*DHAZES_T(K,NW) + & |
---|
| 339 | ASF_T(K+1,NW)*DHAZES_T(K+1,NW) |
---|
| 340 | |
---|
[918] | 341 | WBARI(L,nw,ng) = btemp(L,nw) / DTAUI(L,NW,NG) |
---|
[716] | 342 | else |
---|
| 343 | WBARI(L,nw,ng) = 0.0D0 |
---|
[961] | 344 | DTAUI(L,NW,NG) = 1.0D-9 |
---|
[716] | 345 | endif |
---|
[135] | 346 | |
---|
[961] | 347 | if(btemp(L,nw) .GT. 0.0d0) then |
---|
[918] | 348 | cosbi(L,NW,NG) = atemp/btemp(L,nw) |
---|
[716] | 349 | else |
---|
| 350 | cosbi(L,NW,NG) = 0.0D0 |
---|
| 351 | end if |
---|
[135] | 352 | |
---|
[716] | 353 | END DO ! L vertical loop |
---|
[1648] | 354 | |
---|
[1722] | 355 | ! Last level |
---|
[1648] | 356 | |
---|
[1722] | 357 | L = L_NLAYRAD |
---|
| 358 | K = 2*L+1 |
---|
| 359 | DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) ! + 1.e-50 |
---|
[1648] | 360 | |
---|
[1722] | 361 | atemp = 0. |
---|
| 362 | if(DTAUI(L,NW,NG) .GT. 1.0D-9) then |
---|
| 363 | atemp = atemp + ASF_T(K,NW)*DHAZES_T(K,NW) |
---|
| 364 | WBARI(L,nw,ng) = btemp(L,nw) / DTAUI(L,NW,NG) |
---|
| 365 | else |
---|
| 366 | WBARI(L,nw,ng) = 0.0D0 |
---|
| 367 | DTAUI(L,NW,NG) = 1.0D-9 |
---|
| 368 | endif |
---|
[1648] | 369 | |
---|
[1722] | 370 | if(btemp(L,nw) .GT. 0.0d0) then |
---|
| 371 | cosbi(L,NW,NG) = atemp/btemp(L,nw) |
---|
| 372 | else |
---|
| 373 | cosbi(L,NW,NG) = 0.0D0 |
---|
| 374 | end if |
---|
| 375 | |
---|
| 376 | |
---|
[716] | 377 | ! Now the other Gauss points, if needed. |
---|
[135] | 378 | |
---|
[716] | 379 | DO NG=1,L_NGAUSS-1 |
---|
| 380 | IF(TAUGSURF(NW,NG) .gt. TLIMIT) THEN |
---|
[135] | 381 | |
---|
[1648] | 382 | DO L=1,L_NLAYRAD-1 |
---|
[716] | 383 | K = 2*L+1 |
---|
| 384 | DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)+DTAUKI(K+1,NW,NG)! + 1.e-50 |
---|
| 385 | |
---|
[961] | 386 | if(DTAUI(L,NW,NG) .GT. 1.0D-9) then |
---|
[1722] | 387 | |
---|
| 388 | WBARI(L,nw,ng) = btemp(L,nw) / DTAUI(L,NW,NG) |
---|
| 389 | |
---|
[716] | 390 | else |
---|
| 391 | WBARI(L,nw,ng) = 0.0D0 |
---|
[961] | 392 | DTAUI(L,NW,NG) = 1.0D-9 |
---|
[716] | 393 | endif |
---|
| 394 | |
---|
| 395 | cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS) |
---|
| 396 | END DO ! L vertical loop |
---|
[1648] | 397 | |
---|
[1722] | 398 | ! Last level |
---|
| 399 | L = L_NLAYRAD |
---|
| 400 | K = 2*L+1 |
---|
| 401 | DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)! + 1.e-50 |
---|
[1648] | 402 | |
---|
[1722] | 403 | if(DTAUI(L,NW,NG) .GT. 1.0D-9) then |
---|
| 404 | |
---|
[1648] | 405 | WBARI(L,nw,ng) = btemp(L,nw) / DTAUI(L,NW,NG) |
---|
[1722] | 406 | |
---|
| 407 | else |
---|
| 408 | WBARI(L,nw,ng) = 0.0D0 |
---|
| 409 | DTAUI(L,NW,NG) = 1.0D-9 |
---|
| 410 | endif |
---|
| 411 | |
---|
| 412 | cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS) |
---|
[1648] | 413 | |
---|
[716] | 414 | END IF |
---|
| 415 | |
---|
| 416 | END DO ! NG Gauss loop |
---|
| 417 | END DO ! NW spectral loop |
---|
| 418 | |
---|
| 419 | ! Total extinction optical depths |
---|
| 420 | |
---|
[918] | 421 | DO NG=1,L_NGAUSS ! full gauss loop |
---|
| 422 | DO NW=1,L_NSPECTI |
---|
[716] | 423 | TAUCUMI(1,NW,NG)=0.0D0 |
---|
| 424 | DO K=2,L_LEVELS |
---|
| 425 | TAUCUMI(K,NW,NG)=TAUCUMI(K-1,NW,NG)+DTAUKI(K,NW,NG) |
---|
| 426 | END DO |
---|
| 427 | END DO ! end full gauss loop |
---|
| 428 | END DO |
---|
| 429 | |
---|
| 430 | ! be aware when comparing with textbook results |
---|
| 431 | ! (e.g. Pierrehumbert p. 218) that |
---|
| 432 | ! taucumi does not take the <cos theta>=0.5 factor into |
---|
| 433 | ! account. It is the optical depth for a vertically |
---|
| 434 | ! ascending ray with angle theta = 0. |
---|
| 435 | |
---|
[1897] | 436 | if(firstcall) firstcall = .false. |
---|
| 437 | |
---|
[716] | 438 | return |
---|
| 439 | |
---|
| 440 | |
---|
| 441 | end subroutine optci |
---|
| 442 | |
---|
| 443 | |
---|
| 444 | |
---|