| 1 | SUBROUTINE OPTCV(PQMO,NLAY,ZLEV,PLEV,TMID,PMID, & |
|---|
| 2 | DTAUV,TAUV,TAUCUMV,WBARV,COSBV,TAURAY,TAUGSURF,SEASHAZEFACT,& |
|---|
| 3 | DIAG_OPTH,DIAG_OPTT,CDCOLUMN) |
|---|
| 4 | |
|---|
| 5 | use radinc_h |
|---|
| 6 | use radcommon_h, only: gasv,gasv_recomb,tlimit,Cmk,gzlat_ig, & |
|---|
| 7 | tgasref,pfgasref,wnov,scalep,indv |
|---|
| 8 | use gases_h |
|---|
| 9 | use datafile_mod, only: haze_opt_file |
|---|
| 10 | use comcstfi_mod, only: pi,r |
|---|
| 11 | use callkeys_mod, only: continuum,graybody,callgasvis,corrk_recombin, & |
|---|
| 12 | callclouds,callmufi,seashaze,uncoupl_optic_haze,& |
|---|
| 13 | opt4clouds,FHVIS,FCVIS |
|---|
| 14 | use tracer_h, only: nmicro,nice,ices_indx |
|---|
| 15 | |
|---|
| 16 | implicit none |
|---|
| 17 | |
|---|
| 18 | !================================================================== |
|---|
| 19 | ! |
|---|
| 20 | ! Purpose |
|---|
| 21 | ! ------- |
|---|
| 22 | ! Calculates shortwave optical constants at each level. |
|---|
| 23 | ! |
|---|
| 24 | ! Authors |
|---|
| 25 | ! ------- |
|---|
| 26 | ! Adapted from the NASA Ames code by R. Wordsworth (2009) |
|---|
| 27 | ! |
|---|
| 28 | ! Modified |
|---|
| 29 | ! -------- |
|---|
| 30 | ! J. Vatant d'Ollone (2016-17) |
|---|
| 31 | ! --> Clean and adaptation to Titan |
|---|
| 32 | ! B. de Batz de Trenquelléon (2022-2023) |
|---|
| 33 | ! --> Clean and correction |
|---|
| 34 | ! --> New optics added for Titan's clouds |
|---|
| 35 | ! |
|---|
| 36 | !================================================================== |
|---|
| 37 | ! |
|---|
| 38 | ! THIS SUBROUTINE SETS THE OPTICAL CONSTANTS IN THE VISUAL |
|---|
| 39 | ! IT CALCULATES FOR EACH LAYER, FOR EACH SPECTRAL INTERVAL IN THE VISUAL |
|---|
| 40 | ! LAYER: WBAR, DTAU, COSBAR |
|---|
| 41 | ! LEVEL: TAU |
|---|
| 42 | ! |
|---|
| 43 | ! TAUV(L,NW,NG) is the cumulative optical depth at the top of radiation code |
|---|
| 44 | ! layer L. NW is spectral wavelength interval, ng the Gauss point index. |
|---|
| 45 | ! |
|---|
| 46 | ! TLEV(L) - Temperature at the layer boundary |
|---|
| 47 | ! PLEV(L) - Pressure at the layer boundary (i.e. level) |
|---|
| 48 | ! GASV(NT,NPS,NW,NG) - Visible k-coefficients |
|---|
| 49 | ! |
|---|
| 50 | !------------------------------------------------------------------- |
|---|
| 51 | |
|---|
| 52 | |
|---|
| 53 | !========================================================== |
|---|
| 54 | ! Input/Output |
|---|
| 55 | !========================================================== |
|---|
| 56 | REAL*8, INTENT(IN) :: PQMO(nlay,nmicro) ! Tracers for microphysics optics (X/m2). |
|---|
| 57 | INTEGER, INTENT(IN) :: NLAY ! Number of pressure layers (for pqmo) |
|---|
| 58 | REAL*8, INTENT(IN) :: ZLEV(NLAY+1) |
|---|
| 59 | REAL*8, INTENT(IN) :: PLEV(L_LEVELS) |
|---|
| 60 | REAL*8, INTENT(IN) :: TMID(L_LEVELS), PMID(L_LEVELS) |
|---|
| 61 | REAL*8, INTENT(IN) :: TAURAY(L_NSPECTV) |
|---|
| 62 | REAL*8, INTENT(IN) :: SEASHAZEFACT(L_LEVELS) |
|---|
| 63 | INTEGER, INTENT(IN) :: CDCOLUMN |
|---|
| 64 | |
|---|
| 65 | REAL*8, INTENT(OUT) :: DTAUV(L_NLAYRAD,L_NSPECTV,L_NGAUSS) |
|---|
| 66 | REAL*8, INTENT(OUT) :: TAUV(L_NLEVRAD,L_NSPECTV,L_NGAUSS) |
|---|
| 67 | REAL*8, INTENT(OUT) :: TAUCUMV(L_LEVELS,L_NSPECTV,L_NGAUSS) |
|---|
| 68 | REAL*8, INTENT(OUT) :: COSBV(L_NLAYRAD,L_NSPECTV,L_NGAUSS) |
|---|
| 69 | REAL*8, INTENT(OUT) :: WBARV(L_NLAYRAD,L_NSPECTV,L_NGAUSS) |
|---|
| 70 | REAL*8, INTENT(OUT) :: TAUGSURF(L_NSPECTV,L_NGAUSS-1) |
|---|
| 71 | REAL*8, INTENT(OUT) :: DIAG_OPTH(L_LEVELS,L_NSPECTV,6) |
|---|
| 72 | REAL*8, INTENT(OUT) :: DIAG_OPTT(L_LEVELS,L_NSPECTV,6) |
|---|
| 73 | ! ========================================================== |
|---|
| 74 | |
|---|
| 75 | real*8 DTAUKV(L_LEVELS,L_NSPECTV,L_NGAUSS) |
|---|
| 76 | |
|---|
| 77 | ! Titan customisation |
|---|
| 78 | ! J. Vatant d'Ollone (2016) |
|---|
| 79 | real*8 DHAZE_T(L_LEVELS,L_NSPECTV) |
|---|
| 80 | real*8 DHAZES_T(L_LEVELS,L_NSPECTV) |
|---|
| 81 | real*8 SSA_T(L_LEVELS,L_NSPECTV) |
|---|
| 82 | real*8 ASF_T(L_LEVELS,L_NSPECTV) |
|---|
| 83 | ! ========================== |
|---|
| 84 | |
|---|
| 85 | integer L, NW, NG, K, LK, IAER |
|---|
| 86 | integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS) |
|---|
| 87 | real*8 ANS, TAUGAS |
|---|
| 88 | real*8 TRAY(L_LEVELS,L_NSPECTV) |
|---|
| 89 | real*8 DPR(L_LEVELS), U(L_LEVELS) |
|---|
| 90 | real*8 LCOEF(4), LKCOEF(L_LEVELS,4) |
|---|
| 91 | |
|---|
| 92 | real*8 DCONT |
|---|
| 93 | real*8 DRAYAER |
|---|
| 94 | double precision wn_cont, p_cont, p_air, T_cont, dtemp, dtempc |
|---|
| 95 | double precision p_cross |
|---|
| 96 | |
|---|
| 97 | real*8 KCOEF(4) |
|---|
| 98 | |
|---|
| 99 | ! temporary variable to reduce memory access time to gasv |
|---|
| 100 | real*8 tmpk(2,2) |
|---|
| 101 | |
|---|
| 102 | ! temporary variables for multiple aerosol calculation |
|---|
| 103 | real*8 atemp(L_NLAYRAD,L_NSPECTV) |
|---|
| 104 | real*8 btemp(L_NLAYRAD,L_NSPECTV) |
|---|
| 105 | real*8 ctemp(L_NLAYRAD,L_NSPECTV) |
|---|
| 106 | |
|---|
| 107 | ! variables for k in units m^-1 |
|---|
| 108 | real*8 dz(L_LEVELS) |
|---|
| 109 | |
|---|
| 110 | integer igas, jgas, ilay |
|---|
| 111 | |
|---|
| 112 | integer interm |
|---|
| 113 | |
|---|
| 114 | ! Variables for haze optics |
|---|
| 115 | character(len=200) file_path |
|---|
| 116 | logical file_ok |
|---|
| 117 | integer dumch |
|---|
| 118 | real*8 dumwvl |
|---|
| 119 | integer ilev_cutoff ! pressure level index of cutoff |
|---|
| 120 | real*8 corr_haze |
|---|
| 121 | real*8 corr_alb |
|---|
| 122 | |
|---|
| 123 | ! Variables for new optics |
|---|
| 124 | integer iq, iw, FTYPE, CTYPE |
|---|
| 125 | real*8 m0as,m0af,m0ccn,m3as,m3af,m3ccn,m3cld |
|---|
| 126 | real*8 dtauaer_s,dtauaer_f,dtau_ccn,dtau_cld |
|---|
| 127 | real*8,save :: rhoaer_s(L_NSPECTV),ssa_s(L_NSPECTV),asf_s(L_NSPECTV) |
|---|
| 128 | real*8,save :: rhoaer_f(L_NSPECTV),ssa_f(L_NSPECTV),asf_f(L_NSPECTV) |
|---|
| 129 | real*8,save :: ssa_ccn(L_NSPECTV),asf_ccn(L_NSPECTV) |
|---|
| 130 | real*8,save :: ssa_cld(L_NSPECTV),asf_cld(L_NSPECTV) |
|---|
| 131 | !$OMP THREADPRIVATE(rhoaer_s,rhoaer_f,ssa_s,ssa_f,ssa_cld,asf_s,asf_f,asf_cld) |
|---|
| 132 | |
|---|
| 133 | logical,save :: firstcall=.true. |
|---|
| 134 | !$OMP THREADPRIVATE(firstcall) |
|---|
| 135 | |
|---|
| 136 | |
|---|
| 137 | !! AS: to save time in computing continuum (see bilinearbig) |
|---|
| 138 | IF (.not.ALLOCATED(indv)) THEN |
|---|
| 139 | ALLOCATE(indv(L_NSPECTV,ngasmx,ngasmx)) |
|---|
| 140 | indv = -9999 ! this initial value means "to be calculated" |
|---|
| 141 | ENDIF |
|---|
| 142 | |
|---|
| 143 | ! Some initialisation beacause there's a pb with disr_haze at the limits (nw=1) |
|---|
| 144 | ! I should check this - For now we set vars to zero : better than nans - JVO 2017 |
|---|
| 145 | DHAZE_T(:,:) = 0.0 |
|---|
| 146 | SSA_T(:,:) = 0.0 |
|---|
| 147 | ASF_T(:,:) = 0.0 |
|---|
| 148 | |
|---|
| 149 | ! Load tabulated haze optical properties if needed. |
|---|
| 150 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 151 | IF (firstcall .AND. callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN |
|---|
| 152 | OPEN(12,file=TRIM(haze_opt_file),form='formatted') ! The file has been inquired in physiq_mod firstcall |
|---|
| 153 | READ(12,*) ! dummy header |
|---|
| 154 | DO NW=1,L_NSPECTI |
|---|
| 155 | READ(12,*) ! there's IR 1st |
|---|
| 156 | ENDDO |
|---|
| 157 | DO NW=1,L_NSPECTV |
|---|
| 158 | READ(12,*) dumch, dumwvl, rhoaer_f(nw), ssa_f(nw), asf_f(nw), rhoaer_s(nw), ssa_s(nw), asf_s(nw) |
|---|
| 159 | ENDDO |
|---|
| 160 | CLOSE(12) |
|---|
| 161 | ENDIF |
|---|
| 162 | |
|---|
| 163 | !======================================================================= |
|---|
| 164 | ! Determine the total gas opacity throughout the column, for each |
|---|
| 165 | ! spectral interval, NW, and each Gauss point, NG. |
|---|
| 166 | ! Calculate the continuum opacities, i.e., those that do not depend on |
|---|
| 167 | ! NG, the Gauss index. |
|---|
| 168 | |
|---|
| 169 | taugsurf(:,:) = 0.0 |
|---|
| 170 | dpr(:) = 0.0 |
|---|
| 171 | lkcoef(:,:) = 0.0 |
|---|
| 172 | |
|---|
| 173 | ! Level of cutoff |
|---|
| 174 | !~~~~~~~~~~~~~~~~ |
|---|
| 175 | ilev_cutoff = 26 |
|---|
| 176 | |
|---|
| 177 | do K=2,L_LEVELS |
|---|
| 178 | ilay = L_NLAYRAD+1 - k/2 ! int. arithmetic => gives the gcm layer index (reversed) |
|---|
| 179 | DPR(k) = PLEV(K)-PLEV(K-1) |
|---|
| 180 | |
|---|
| 181 | ! if we have continuum opacities, we need dz |
|---|
| 182 | dz(k) = dpr(k)*R*TMID(K)/(gzlat_ig(ilay)*PMID(K)) |
|---|
| 183 | U(k) = Cmk(ilay)*DPR(k) ! only Cmk line in optcv.F |
|---|
| 184 | |
|---|
| 185 | call tpindex(PMID(K),TMID(K),pfgasref,tgasref,LCOEF,MT(K),MP(K)) |
|---|
| 186 | |
|---|
| 187 | do LK=1,4 |
|---|
| 188 | LKCOEF(K,LK) = LCOEF(LK) |
|---|
| 189 | end do |
|---|
| 190 | end do ! L_LEVELS |
|---|
| 191 | |
|---|
| 192 | ! Rayleigh scattering |
|---|
| 193 | do NW=1,L_NSPECTV |
|---|
| 194 | TRAY(1:4,NW) = 1.d-30 |
|---|
| 195 | do K=5,L_LEVELS |
|---|
| 196 | TRAY(K,NW) = TAURAY(NW) * DPR(K) |
|---|
| 197 | end do ! L_LEVELS |
|---|
| 198 | end do |
|---|
| 199 | |
|---|
| 200 | DIAG_OPTH(:,:,:) = 0.D0 |
|---|
| 201 | DIAG_OPTT(:,:,:) = 0.D0 |
|---|
| 202 | |
|---|
| 203 | do NW=1,L_NSPECTV |
|---|
| 204 | ! We ignore K=1... |
|---|
| 205 | do K=2,L_LEVELS |
|---|
| 206 | ! int. arithmetic => gives the gcm layer index (reversed) |
|---|
| 207 | ilay = L_NLAYRAD+1 - k/2 |
|---|
| 208 | |
|---|
| 209 | ! Optics coupled with the microphysics : |
|---|
| 210 | IF (callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN |
|---|
| 211 | |
|---|
| 212 | !========================================================================================== |
|---|
| 213 | ! Old optics (must have callclouds = .false.): |
|---|
| 214 | !========================================================================================== |
|---|
| 215 | IF (.NOT. opt4clouds) THEN |
|---|
| 216 | m3as = pqmo(ilay,2) / 2.0 |
|---|
| 217 | m3af = pqmo(ilay,4) / 2.0 |
|---|
| 218 | ! Cut-off (here for p = 2.7e3Pa / alt = 70km) |
|---|
| 219 | IF (ilay .lt. ilev_cutoff) THEN |
|---|
| 220 | m3as = pqmo(ilev_cutoff,2) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) |
|---|
| 221 | m3af = pqmo(ilev_cutoff,4) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) |
|---|
| 222 | ENDIF |
|---|
| 223 | |
|---|
| 224 | dtauaer_s = m3as*rhoaer_s(nw) |
|---|
| 225 | dtauaer_f = m3af*rhoaer_f(nw) |
|---|
| 226 | |
|---|
| 227 | !========================================================================================== |
|---|
| 228 | ! New optics : |
|---|
| 229 | !========================================================================================== |
|---|
| 230 | ELSE |
|---|
| 231 | iw = (L_NSPECTV + 1) - NW ! Visible first and return |
|---|
| 232 | !----------------------------- |
|---|
| 233 | ! HAZE (Spherical + Fractal) : |
|---|
| 234 | !----------------------------- |
|---|
| 235 | FTYPE = 1 |
|---|
| 236 | |
|---|
| 237 | ! Spherical aerosols : |
|---|
| 238 | !--------------------- |
|---|
| 239 | CTYPE = 5 |
|---|
| 240 | m0as = pqmo(ilay,1) / 2.0 |
|---|
| 241 | m3as = pqmo(ilay,2) / 2.0 |
|---|
| 242 | ! If not callclouds : must have a cut-off (here for p = 2.7e3Pa / alt = 70km) |
|---|
| 243 | IF (.NOT. callclouds) THEN |
|---|
| 244 | IF (ilay .lt. ilev_cutoff) THEN |
|---|
| 245 | m0as = pqmo(ilev_cutoff,1) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) |
|---|
| 246 | m3as = pqmo(ilev_cutoff,2) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) |
|---|
| 247 | ENDIF |
|---|
| 248 | ENDIF |
|---|
| 249 | call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0as,m3as,iw,dtauaer_s,ssa_s(nw),asf_s(nw)) |
|---|
| 250 | |
|---|
| 251 | ! Fractal aerosols : |
|---|
| 252 | !------------------- |
|---|
| 253 | CTYPE = FTYPE |
|---|
| 254 | m0af = pqmo(ilay,3) / 2.0 |
|---|
| 255 | m3af = pqmo(ilay,4) / 2.0 |
|---|
| 256 | ! If not callclouds : must have a cut-off (here for p = 2.7e3Pa / alt = 70km) |
|---|
| 257 | IF (.NOT. callclouds) THEN |
|---|
| 258 | IF (ilay .lt. ilev_cutoff) THEN |
|---|
| 259 | m0af = pqmo(ilev_cutoff,3) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) |
|---|
| 260 | m3af = pqmo(ilev_cutoff,4) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) |
|---|
| 261 | ENDIF |
|---|
| 262 | ENDIF |
|---|
| 263 | call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0af,m3af,iw,dtauaer_f,ssa_f(nw),asf_f(nw)) |
|---|
| 264 | ENDIF |
|---|
| 265 | |
|---|
| 266 | ! Tuning of optical properties for haze : |
|---|
| 267 | !dtauaer_s = dtauaer_s * (FHVIS * (1-ssa_s(nw)) + ssa_s(nw)) |
|---|
| 268 | !ssa_s(nw) = ssa_s(nw) / (FHVIS * (1-ssa_s(nw)) + ssa_s(nw)) |
|---|
| 269 | dtauaer_f = dtauaer_f * (FHVIS * (1-ssa_f(nw)) + ssa_f(nw)) |
|---|
| 270 | ssa_f(nw) = ssa_f(nw) / (FHVIS * (1-ssa_f(nw)) + ssa_f(nw)) |
|---|
| 271 | |
|---|
| 272 | ! Total of Haze opacity (dtau), SSA (w) and ASF (COS) : |
|---|
| 273 | DHAZE_T(k,nw) = dtauaer_s + dtauaer_f |
|---|
| 274 | IF (dtauaer_s + dtauaer_f .GT. 1.D-30) THEN |
|---|
| 275 | SSA_T(k,nw) = ( dtauaer_s*ssa_s(nw) + dtauaer_f*ssa_f(nw) ) / ( dtauaer_s+dtauaer_f ) |
|---|
| 276 | ASF_T(k,nw) = ( dtauaer_s*ssa_s(nw)*asf_s(nw) + dtauaer_f*ssa_f(nw)*asf_f(nw) ) & |
|---|
| 277 | / ( ssa_s(nw)*dtauaer_s + ssa_f(nw)*dtauaer_f ) |
|---|
| 278 | ELSE |
|---|
| 279 | DHAZE_T(k,nw) = 0.D0 |
|---|
| 280 | SSA_T(k,nw) = 1.0 |
|---|
| 281 | ASF_T(k,nw) = 1.0 |
|---|
| 282 | ENDIF |
|---|
| 283 | |
|---|
| 284 | ! Opacity and albedo adjustement below cutoff |
|---|
| 285 | !---------------------------------------------- |
|---|
| 286 | corr_haze=0.85-0.15*TANH((PMID(k)*100.-2500.)/250.) |
|---|
| 287 | corr_alb =0.85-0.15*TANH((PMID(k)*100.-4500.)/1000.) |
|---|
| 288 | IF (.NOT. callclouds) THEN |
|---|
| 289 | IF (ilay .lt. ilev_cutoff) THEN |
|---|
| 290 | DHAZE_T(k,nw) = DHAZE_T(k,nw) * corr_haze |
|---|
| 291 | SSA_T(k,nw) = SSA_T(k,nw) * corr_alb |
|---|
| 292 | ENDIF |
|---|
| 293 | ENDIF |
|---|
| 294 | |
|---|
| 295 | ! Diagnostics for the haze : |
|---|
| 296 | DIAG_OPTH(k,nw,1) = DHAZE_T(k,nw) ! dtau |
|---|
| 297 | DIAG_OPTH(k,nw,2) = SSA_T(k,nw) ! wbar |
|---|
| 298 | DIAG_OPTH(k,nw,3) = ASF_T(k,nw) ! gbar |
|---|
| 299 | |
|---|
| 300 | !--------------------- |
|---|
| 301 | ! CLOUDS (Spherical) : |
|---|
| 302 | !--------------------- |
|---|
| 303 | IF (callclouds) THEN |
|---|
| 304 | CTYPE = 0 |
|---|
| 305 | m0ccn = pqmo(ilay,5) / 2.0 |
|---|
| 306 | m3ccn = pqmo(ilay,6) / 2.0 |
|---|
| 307 | m3cld = pqmo(ilay,6) / 2.0 |
|---|
| 308 | |
|---|
| 309 | ! Clear / Dark column method : |
|---|
| 310 | !----------------------------- |
|---|
| 311 | |
|---|
| 312 | ! CCN's SSA : |
|---|
| 313 | call get_haze_and_cloud_opacity(FTYPE,FTYPE,m0ccn,m3ccn,iw,dtau_ccn,ssa_ccn(nw),asf_ccn(nw)) |
|---|
| 314 | |
|---|
| 315 | ! Clear column (CCN, C2H2, C2H6, HCN, AC6H6) : |
|---|
| 316 | IF (CDCOLUMN == 0) THEN |
|---|
| 317 | DO iq = 2, nice |
|---|
| 318 | m3cld = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0) |
|---|
| 319 | ENDDO |
|---|
| 320 | call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw)) |
|---|
| 321 | |
|---|
| 322 | ! Dark column (CCN, CH4, C2H2, C2H6, HCN, AC6H6) : |
|---|
| 323 | ELSEIF (CDCOLUMN == 1) THEN |
|---|
| 324 | DO iq = 1, nice |
|---|
| 325 | m3cld = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0) |
|---|
| 326 | ENDDO |
|---|
| 327 | call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw)) |
|---|
| 328 | |
|---|
| 329 | ELSE |
|---|
| 330 | WRITE(*,*) 'WARNING OPTCV.F90 : CDCOLUMN MUST BE 0 OR 1' |
|---|
| 331 | WRITE(*,*) 'WE USE DARK COLUMN ...' |
|---|
| 332 | DO iq = 1, nice |
|---|
| 333 | m3cld = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0) |
|---|
| 334 | ENDDO |
|---|
| 335 | call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw)) |
|---|
| 336 | ENDIF |
|---|
| 337 | |
|---|
| 338 | ! For small dropplets, opacity of nucleus dominates |
|---|
| 339 | IF ((m3ccn + m3cld) .le. tiny(m3ccn)) THEN !no cloud !emoisan tests |
|---|
| 340 | dtau_cld = 0. |
|---|
| 341 | ssa_cld(nw) = 0. |
|---|
| 342 | ELSE |
|---|
| 343 | dtau_cld = (dtau_cld*m3ccn + dtau_cld*m3cld) / (m3ccn + m3cld) |
|---|
| 344 | ssa_cld(nw) = (ssa_ccn(nw)*m3ccn + ssa_cld(nw)*m3cld) / (m3ccn + m3cld) |
|---|
| 345 | ENDIF |
|---|
| 346 | |
|---|
| 347 | ! Tuning of optical properties for clouds : |
|---|
| 348 | dtau_cld = dtau_cld * (FCVIS * (1-ssa_cld(nw)) + ssa_cld(nw)) |
|---|
| 349 | ssa_cld(nw) = ssa_cld(nw) / (FCVIS * (1-ssa_cld(nw)) + ssa_cld(nw)) |
|---|
| 350 | |
|---|
| 351 | ! Total of Haze + Clouds opacity (dtau), SSA (w) and ASF (COS) : |
|---|
| 352 | DHAZE_T(k,nw) = dtauaer_s + dtauaer_f + dtau_cld |
|---|
| 353 | IF (DHAZE_T(k,nw) .GT. 1.D-30) THEN |
|---|
| 354 | SSA_T(k,nw) = ( dtauaer_s*ssa_s(nw) + dtauaer_f*ssa_f(nw) + dtau_cld*ssa_cld(nw) ) / ( dtauaer_s+dtauaer_f+dtau_cld ) |
|---|
| 355 | ASF_T(k,nw) = ( dtauaer_s*ssa_s(nw)*asf_s(nw) + dtauaer_f*ssa_f(nw)*asf_f(nw) + dtau_cld*ssa_cld(nw)*asf_cld(nw) ) & |
|---|
| 356 | / ( ssa_s(nw)*dtauaer_s + ssa_f(nw)*dtauaer_f + ssa_cld(nw)*dtau_cld ) |
|---|
| 357 | ELSE |
|---|
| 358 | DHAZE_T(k,nw) = 0.D0 |
|---|
| 359 | SSA_T(k,nw) = 1.0 |
|---|
| 360 | ASF_T(k,nw) = 1.0 |
|---|
| 361 | ENDIF |
|---|
| 362 | |
|---|
| 363 | ! Diagnostics for clouds : |
|---|
| 364 | DIAG_OPTT(k,nw,1) = DHAZE_T(k,nw) ! dtau |
|---|
| 365 | DIAG_OPTT(k,nw,2) = SSA_T(k,nw) ! wbar |
|---|
| 366 | DIAG_OPTT(k,nw,3) = ASF_T(k,nw) ! gbar |
|---|
| 367 | |
|---|
| 368 | ELSE |
|---|
| 369 | ! Diagnostics for clouds : |
|---|
| 370 | DIAG_OPTT(k,nw,1) = 0.D0 ! dtau |
|---|
| 371 | DIAG_OPTT(k,nw,2) = 1.0 ! wbar |
|---|
| 372 | DIAG_OPTT(k,nw,3) = 1.0 ! gbar |
|---|
| 373 | ENDIF |
|---|
| 374 | |
|---|
| 375 | ! Optics and microphysics no coupled : |
|---|
| 376 | ELSE |
|---|
| 377 | ! Call fixed vertical haze profile of extinction - same for all columns |
|---|
| 378 | call disr_haze(dz(k),plev(k),wnov(nw),DHAZE_T(k,nw),SSA_T(k,nw),ASF_T(k,nw)) |
|---|
| 379 | if (seashaze) DHAZE_T(k,nw) = DHAZE_T(k,nw)*seashazefact(k) |
|---|
| 380 | ! Diagnostics for the haze : |
|---|
| 381 | DIAG_OPTH(k,nw,1) = DHAZE_T(k,nw) ! dtau |
|---|
| 382 | DIAG_OPTH(k,nw,2) = SSA_T(k,nw) ! wbar |
|---|
| 383 | DIAG_OPTH(k,nw,3) = ASF_T(k,nw) ! gbar |
|---|
| 384 | ! Diagnostics for clouds : |
|---|
| 385 | DIAG_OPTT(k,nw,1) = 0.D0 ! dtau |
|---|
| 386 | DIAG_OPTT(k,nw,2) = 1.0 ! wbar |
|---|
| 387 | DIAG_OPTT(k,nw,3) = 1.0 ! gbar |
|---|
| 388 | ENDIF ! ENDIF callmufi |
|---|
| 389 | |
|---|
| 390 | !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR |
|---|
| 391 | ! but visible does not handle very well diffusion in first layer. |
|---|
| 392 | ! The tauaero and tauray are thus set to 0 (a small value for rayleigh because the code crashes otherwise) |
|---|
| 393 | ! in the 4 first semilayers in optcv, but not optci. |
|---|
| 394 | ! This solves random variations of the sw heating at the model top. |
|---|
| 395 | if (k<5) DHAZE_T(K,:) = 0.0 |
|---|
| 396 | |
|---|
| 397 | DRAYAER = TRAY(K,NW) |
|---|
| 398 | ! DRAYAER is Tau RAYleigh scattering, plus AERosol opacity |
|---|
| 399 | DRAYAER = DRAYAER + DHAZE_T(K,NW) ! Titan's aerosol |
|---|
| 400 | |
|---|
| 401 | DCONT = 0.0 ! continuum absorption |
|---|
| 402 | |
|---|
| 403 | if(continuum.and.(.not.graybody).and.callgasvis)then |
|---|
| 404 | ! include continua if necessary |
|---|
| 405 | wn_cont = dble(wnov(nw)) |
|---|
| 406 | T_cont = dble(TMID(k)) |
|---|
| 407 | do igas=1,ngasmx |
|---|
| 408 | |
|---|
| 409 | p_cont = dble(PMID(k)*scalep*gfrac(igas,ilay)) |
|---|
| 410 | |
|---|
| 411 | dtemp=0.0 |
|---|
| 412 | if(igas.eq.igas_N2)then |
|---|
| 413 | |
|---|
| 414 | interm = indv(nw,igas,igas) |
|---|
| 415 | ! call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
|---|
| 416 | indv(nw,igas,igas) = interm |
|---|
| 417 | ! only goes to 500 cm^-1, so unless we're around a cold brown dwarf, this is irrelevant in the visible |
|---|
| 418 | |
|---|
| 419 | elseif(igas.eq.igas_H2)then |
|---|
| 420 | |
|---|
| 421 | ! first do self-induced absorption |
|---|
| 422 | interm = indv(nw,igas,igas) |
|---|
| 423 | call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
|---|
| 424 | indv(nw,igas,igas) = interm |
|---|
| 425 | |
|---|
| 426 | ! then cross-interactions with other gases |
|---|
| 427 | do jgas=1,ngasmx |
|---|
| 428 | p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay)) |
|---|
| 429 | dtempc = 0.0 |
|---|
| 430 | if(jgas.eq.igas_N2)then |
|---|
| 431 | interm = indv(nw,igas,jgas) |
|---|
| 432 | call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) |
|---|
| 433 | indv(nw,igas,jgas) = interm |
|---|
| 434 | ! should be irrelevant in the visible |
|---|
| 435 | endif |
|---|
| 436 | dtemp = dtemp + dtempc |
|---|
| 437 | enddo |
|---|
| 438 | |
|---|
| 439 | elseif(igas.eq.igas_CH4)then |
|---|
| 440 | |
|---|
| 441 | ! first do self-induced absorption |
|---|
| 442 | interm = indv(nw,igas,igas) |
|---|
| 443 | call interpolateCH4CH4(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
|---|
| 444 | indv(nw,igas,igas) = interm |
|---|
| 445 | |
|---|
| 446 | ! then cross-interactions with other gases |
|---|
| 447 | do jgas=1,ngasmx |
|---|
| 448 | p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay)) |
|---|
| 449 | dtempc = 0.0 |
|---|
| 450 | if(jgas.eq.igas_N2)then |
|---|
| 451 | interm = indv(nw,igas,jgas) |
|---|
| 452 | call interpolateN2CH4(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) |
|---|
| 453 | indv(nw,igas,jgas) = interm |
|---|
| 454 | endif |
|---|
| 455 | dtemp = dtemp + dtempc |
|---|
| 456 | enddo |
|---|
| 457 | |
|---|
| 458 | endif |
|---|
| 459 | |
|---|
| 460 | DCONT = DCONT + dtemp |
|---|
| 461 | |
|---|
| 462 | enddo |
|---|
| 463 | |
|---|
| 464 | DCONT = DCONT*dz(k) |
|---|
| 465 | |
|---|
| 466 | endif |
|---|
| 467 | |
|---|
| 468 | do ng=1,L_NGAUSS-1 |
|---|
| 469 | |
|---|
| 470 | ! Now compute TAUGAS |
|---|
| 471 | |
|---|
| 472 | ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically |
|---|
| 473 | ! the execution time of optci/v -> ~ factor 2 on the whole radiative |
|---|
| 474 | ! transfer on the tested simulations ! |
|---|
| 475 | |
|---|
| 476 | if (corrk_recombin) then |
|---|
| 477 | tmpk = GASV_RECOMB(MT(K):MT(K)+1,MP(K):MP(K)+1,NW,NG) |
|---|
| 478 | else |
|---|
| 479 | tmpk = GASV(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG) |
|---|
| 480 | endif |
|---|
| 481 | |
|---|
| 482 | KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG) |
|---|
| 483 | KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG) |
|---|
| 484 | KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG) |
|---|
| 485 | KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG) |
|---|
| 486 | |
|---|
| 487 | ! Interpolate the gaseous k-coefficients to the requested T,P values |
|---|
| 488 | |
|---|
| 489 | ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) + & |
|---|
| 490 | LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4) |
|---|
| 491 | |
|---|
| 492 | |
|---|
| 493 | TAUGAS = U(k)*ANS |
|---|
| 494 | |
|---|
| 495 | TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT |
|---|
| 496 | DTAUKV(K,nw,ng) = TAUGAS & |
|---|
| 497 | + DRAYAER & ! DRAYAER includes all scattering contributions |
|---|
| 498 | + DCONT ! For parameterized continuum aborption |
|---|
| 499 | end do |
|---|
| 500 | |
|---|
| 501 | ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS), |
|---|
| 502 | ! which holds continuum opacity only |
|---|
| 503 | |
|---|
| 504 | NG = L_NGAUSS |
|---|
| 505 | DTAUKV(K,nw,ng) = DRAYAER + DCONT ! Scattering + parameterized continuum absorption, including Titan's haze |
|---|
| 506 | |
|---|
| 507 | DIAG_OPTH(K,nw,4) = DRAYAER |
|---|
| 508 | DIAG_OPTH(K,nw,5) = TAUGAS |
|---|
| 509 | DIAG_OPTH(K,nw,6) = DCONT |
|---|
| 510 | DIAG_OPTT(K,nw,4) = DRAYAER |
|---|
| 511 | DIAG_OPTT(K,nw,5) = TAUGAS |
|---|
| 512 | DIAG_OPTT(K,nw,6) = DCONT |
|---|
| 513 | |
|---|
| 514 | end do ! K = L_LEVELS |
|---|
| 515 | end do ! nw = L_NSPECTV |
|---|
| 516 | |
|---|
| 517 | !======================================================================= |
|---|
| 518 | ! Now the full treatment for the layers, where besides the opacity |
|---|
| 519 | ! we need to calculate the scattering albedo and asymmetry factors |
|---|
| 520 | ! ====================================================================== |
|---|
| 521 | |
|---|
| 522 | ! Haze scattering |
|---|
| 523 | !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR |
|---|
| 524 | ! but not in the visible |
|---|
| 525 | ! The dhaze_s is thus set to 0 in the 4 first semilayers in optcv, but not optci. |
|---|
| 526 | ! This solves random variations of the sw heating at the model top. |
|---|
| 527 | DO NW=1,L_NSPECTV |
|---|
| 528 | DHAZES_T(1:4,NW) = 0.d0 |
|---|
| 529 | DO K=5,L_LEVELS |
|---|
| 530 | DHAZES_T(K,NW) = DHAZE_T(K,NW) * SSA_T(K,NW) ! effect of scattering albedo on haze |
|---|
| 531 | ENDDO |
|---|
| 532 | ENDDO |
|---|
| 533 | |
|---|
| 534 | ! NW spectral loop |
|---|
| 535 | DO NW=1,L_NSPECTV |
|---|
| 536 | ! L vertical loop |
|---|
| 537 | DO L=1,L_NLAYRAD-1 |
|---|
| 538 | K = 2*L+1 |
|---|
| 539 | atemp(L,NW) = ASF_T(K,NW)*DHAZES_T(K,NW) + ASF_T(K+1,NW)*DHAZES_T(K+1,NW) |
|---|
| 540 | btemp(L,NW) = DHAZES_T(K,NW) + DHAZES_T(K+1,NW) |
|---|
| 541 | ctemp(L,NW) = btemp(L,NW) + 0.9999*(TRAY(K,NW) + TRAY(K+1,NW)) ! JVO 2017 : does this 0.999 is really meaningful ? |
|---|
| 542 | btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) + TRAY(K+1,NW) |
|---|
| 543 | COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW) |
|---|
| 544 | END DO |
|---|
| 545 | |
|---|
| 546 | ! Last level |
|---|
| 547 | L = L_NLAYRAD |
|---|
| 548 | K = 2*L+1 |
|---|
| 549 | atemp(L,NW) = ASF_T(K,NW)*DHAZES_T(K,NW) |
|---|
| 550 | btemp(L,NW) = DHAZES_T(K,NW) |
|---|
| 551 | ctemp(L,NW) = btemp(L,NW) + 0.9999*TRAY(K,NW) ! JVO 2017 : does this 0.999 is really meaningful ? |
|---|
| 552 | btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) |
|---|
| 553 | COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW) |
|---|
| 554 | END DO |
|---|
| 555 | |
|---|
| 556 | ! NG Gauss loop |
|---|
| 557 | DO NG=1,L_NGAUSS |
|---|
| 558 | ! NW spectral loop |
|---|
| 559 | DO NW=1,L_NSPECTV |
|---|
| 560 | ! L vertical loop |
|---|
| 561 | DO L=1,L_NLAYRAD-1 |
|---|
| 562 | K = 2*L+1 |
|---|
| 563 | DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG) |
|---|
| 564 | WBARV(L,nw,ng) = ctemp(L,NW) / DTAUV(L,nw,ng) |
|---|
| 565 | END DO |
|---|
| 566 | |
|---|
| 567 | ! Last level |
|---|
| 568 | L = L_NLAYRAD |
|---|
| 569 | K = 2*L+1 |
|---|
| 570 | DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) |
|---|
| 571 | WBARV(L,NW,NG) = ctemp(L,NW) / DTAUV(L,NW,NG) |
|---|
| 572 | END DO |
|---|
| 573 | END DO |
|---|
| 574 | |
|---|
| 575 | ! Total extinction optical depths |
|---|
| 576 | !DO NG=1,L_NGAUSS ! full gauss loop |
|---|
| 577 | ! DO NW=1,L_NSPECTV |
|---|
| 578 | ! TAUCUMV(1,NW,NG)=0.0D0 |
|---|
| 579 | ! DO K=2,L_LEVELS |
|---|
| 580 | ! TAUCUMV(K,NW,NG)=TAUCUMV(K-1,NW,NG)+DTAUKV(K,NW,NG) |
|---|
| 581 | ! END DO |
|---|
| 582 | ! DO L=1,L_NLAYRAD |
|---|
| 583 | ! TAUV(L,NW,NG)=TAUCUMV(2*L,NW,NG) |
|---|
| 584 | ! END DO |
|---|
| 585 | ! TAUV(L,NW,NG)=TAUCUMV(2*L_NLAYRAD+1,NW,NG) |
|---|
| 586 | ! END DO |
|---|
| 587 | !END DO ! end full gauss loop |
|---|
| 588 | |
|---|
| 589 | TAUCUMV(:,:,:) = DTAUKV(:,:,:) |
|---|
| 590 | DO L=1,L_NLAYRAD |
|---|
| 591 | TAUV(L,:,:)=TAUCUMV(2*L,:,:) |
|---|
| 592 | END DO |
|---|
| 593 | TAUV(L,:,:)=TAUCUMV(2*L_NLAYRAD+1,:,:) |
|---|
| 594 | |
|---|
| 595 | if(firstcall) firstcall = .false. |
|---|
| 596 | |
|---|
| 597 | return |
|---|
| 598 | |
|---|
| 599 | |
|---|
| 600 | end subroutine optcv |
|---|