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 to Titan |
---|
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 | |
---|
120 | ! Variables for new optics |
---|
121 | integer iq, iw, FTYPE, CTYPE |
---|
122 | real*8 m0as,m0af,m0ccn,m3as,m3af,m3ccn,m3cld |
---|
123 | real*8 dtauaer_s,dtauaer_f,dtau_ccn,dtau_cld |
---|
124 | real*8,save :: rhoaer_s(L_NSPECTV),ssa_s(L_NSPECTV),asf_s(L_NSPECTV) |
---|
125 | real*8,save :: rhoaer_f(L_NSPECTV),ssa_f(L_NSPECTV),asf_f(L_NSPECTV) |
---|
126 | real*8,save :: ssa_ccn(L_NSPECTV),asf_ccn(L_NSPECTV) |
---|
127 | real*8,save :: ssa_cld(L_NSPECTV),asf_cld(L_NSPECTV) |
---|
128 | !$OMP THREADPRIVATE(rhoaer_s,rhoaer_f,ssa_s,ssa_f,ssa_cld,asf_s,asf_f,asf_cld) |
---|
129 | |
---|
130 | logical,save :: firstcall=.true. |
---|
131 | !$OMP THREADPRIVATE(firstcall) |
---|
132 | |
---|
133 | |
---|
134 | !! AS: to save time in computing continuum (see bilinearbig) |
---|
135 | IF (.not.ALLOCATED(indv)) THEN |
---|
136 | ALLOCATE(indv(L_NSPECTV,ngasmx,ngasmx)) |
---|
137 | indv = -9999 ! this initial value means "to be calculated" |
---|
138 | ENDIF |
---|
139 | |
---|
140 | ! Some initialisation beacause there's a pb with disr_haze at the limits (nw=1) |
---|
141 | ! I should check this - For now we set vars to zero : better than nans - JVO 2017 |
---|
142 | DHAZE_T(:,:) = 0.0 |
---|
143 | SSA_T(:,:) = 0.0 |
---|
144 | ASF_T(:,:) = 0.0 |
---|
145 | |
---|
146 | ! Load tabulated haze optical properties if needed. |
---|
147 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
148 | IF (firstcall .AND. callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN |
---|
149 | OPEN(12,file=TRIM(haze_opt_file),form='formatted') ! The file has been inquired in physiq_mod firstcall |
---|
150 | READ(12,*) ! dummy header |
---|
151 | DO NW=1,L_NSPECTI |
---|
152 | READ(12,*) ! there's IR 1st |
---|
153 | ENDDO |
---|
154 | DO NW=1,L_NSPECTV |
---|
155 | READ(12,*) dumch, dumwvl, rhoaer_f(nw), ssa_f(nw), asf_f(nw), rhoaer_s(nw), ssa_s(nw), asf_s(nw) |
---|
156 | ENDDO |
---|
157 | CLOSE(12) |
---|
158 | ENDIF |
---|
159 | |
---|
160 | !======================================================================= |
---|
161 | ! Determine the total gas opacity throughout the column, for each |
---|
162 | ! spectral interval, NW, and each Gauss point, NG. |
---|
163 | ! Calculate the continuum opacities, i.e., those that do not depend on |
---|
164 | ! NG, the Gauss index. |
---|
165 | |
---|
166 | taugsurf(:,:) = 0.0 |
---|
167 | dpr(:) = 0.0 |
---|
168 | lkcoef(:,:) = 0.0 |
---|
169 | |
---|
170 | do K=2,L_LEVELS |
---|
171 | ilay = L_NLAYRAD+1 - k/2 ! int. arithmetic => gives the gcm layer index (reversed) |
---|
172 | DPR(k) = PLEV(K)-PLEV(K-1) |
---|
173 | |
---|
174 | ! if we have continuum opacities, we need dz |
---|
175 | dz(k) = dpr(k)*R*TMID(K)/(gzlat_ig(ilay)*PMID(K)) |
---|
176 | U(k) = Cmk(ilay)*DPR(k) ! only Cmk line in optcv.F |
---|
177 | |
---|
178 | call tpindex(PMID(K),TMID(K),pfgasref,tgasref,LCOEF,MT(K),MP(K)) |
---|
179 | |
---|
180 | do LK=1,4 |
---|
181 | LKCOEF(K,LK) = LCOEF(LK) |
---|
182 | end do |
---|
183 | end do ! L_LEVELS |
---|
184 | |
---|
185 | ! Rayleigh scattering |
---|
186 | do NW=1,L_NSPECTV |
---|
187 | TRAY(1:4,NW) = 1.d-30 |
---|
188 | do K=5,L_LEVELS |
---|
189 | TRAY(K,NW) = TAURAY(NW) * DPR(K) |
---|
190 | end do ! L_LEVELS |
---|
191 | end do |
---|
192 | |
---|
193 | DIAG_OPTH(:,:,:) = 0.D0 |
---|
194 | DIAG_OPTT(:,:,:) = 0.D0 |
---|
195 | |
---|
196 | do NW=1,L_NSPECTV |
---|
197 | ! We ignore K=1... |
---|
198 | do K=2,L_LEVELS |
---|
199 | ! int. arithmetic => gives the gcm layer index (reversed) |
---|
200 | ilay = L_NLAYRAD+1 - k/2 |
---|
201 | |
---|
202 | ! Optics coupled with the microphysics : |
---|
203 | IF (callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN |
---|
204 | |
---|
205 | !========================================================================================== |
---|
206 | ! Old optics (must have callclouds = .false.): |
---|
207 | !========================================================================================== |
---|
208 | IF (.NOT. opt4clouds) THEN |
---|
209 | m3as = pqmo(ilay,2) / 2.0 |
---|
210 | m3af = pqmo(ilay,4) / 2.0 |
---|
211 | ! Cut-off (here for p = 2.7e3Pa / alt = 70km) |
---|
212 | IF (ilay .lt. 23) THEN |
---|
213 | m3as = pqmo(23,2) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23)) |
---|
214 | m3af = pqmo(23,4) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23)) |
---|
215 | ENDIF |
---|
216 | |
---|
217 | dtauaer_s = m3as*rhoaer_s(nw) |
---|
218 | dtauaer_f = m3af*rhoaer_f(nw) |
---|
219 | |
---|
220 | !========================================================================================== |
---|
221 | ! New optics : |
---|
222 | !========================================================================================== |
---|
223 | ELSE |
---|
224 | iw = (L_NSPECTV + 1) - NW ! Visible first and return |
---|
225 | !----------------------------- |
---|
226 | ! HAZE (Spherical + Fractal) : |
---|
227 | !----------------------------- |
---|
228 | FTYPE = 1 |
---|
229 | |
---|
230 | ! Spherical aerosols : |
---|
231 | !--------------------- |
---|
232 | CTYPE = 5 |
---|
233 | m0as = pqmo(ilay,1) / 2.0 |
---|
234 | m3as = pqmo(ilay,2) / 2.0 |
---|
235 | ! If not callclouds : must have a cut-off (here for p = 2.7e3Pa / alt = 70km) |
---|
236 | IF (.NOT. callclouds) THEN |
---|
237 | IF (ilay .lt. 23) THEN |
---|
238 | m0as = pqmo(23,1) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23)) |
---|
239 | m3as = pqmo(23,2) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23)) |
---|
240 | ENDIF |
---|
241 | ENDIF |
---|
242 | call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0as,m3as,iw,dtauaer_s,ssa_s(nw),asf_s(nw)) |
---|
243 | |
---|
244 | ! Fractal aerosols : |
---|
245 | !------------------- |
---|
246 | CTYPE = FTYPE |
---|
247 | m0af = pqmo(ilay,3) / 2.0 |
---|
248 | m3af = pqmo(ilay,4) / 2.0 |
---|
249 | ! If not callclouds : must have a cut-off (here for p = 2.7e3Pa / alt = 70km) |
---|
250 | IF (.NOT. callclouds) THEN |
---|
251 | IF (ilay .lt. 23) THEN |
---|
252 | m0af = pqmo(23,3) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23)) |
---|
253 | m3af = pqmo(23,4) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23)) |
---|
254 | ENDIF |
---|
255 | ENDIF |
---|
256 | call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0af,m3af,iw,dtauaer_f,ssa_f(nw),asf_f(nw)) |
---|
257 | ENDIF |
---|
258 | |
---|
259 | ! Tuning of optical properties for haze : |
---|
260 | !dtauaer_s = dtauaer_s * (FHVIS * (1-ssa_s(nw)) + ssa_s(nw)) |
---|
261 | !ssa_s(nw) = ssa_s(nw) / (FHVIS * (1-ssa_s(nw)) + ssa_s(nw)) |
---|
262 | dtauaer_f = dtauaer_f * (FHVIS * (1-ssa_f(nw)) + ssa_f(nw)) |
---|
263 | ssa_f(nw) = ssa_f(nw) / (FHVIS * (1-ssa_f(nw)) + ssa_f(nw)) |
---|
264 | |
---|
265 | ! Total of Haze opacity (dtau), SSA (w) and ASF (COS) : |
---|
266 | DHAZE_T(k,nw) = dtauaer_s + dtauaer_f |
---|
267 | IF (dtauaer_s + dtauaer_f .GT. 1.D-30) THEN |
---|
268 | SSA_T(k,nw) = ( dtauaer_s*ssa_s(nw) + dtauaer_f*ssa_f(nw) ) / ( dtauaer_s+dtauaer_f ) |
---|
269 | ASF_T(k,nw) = ( dtauaer_s*ssa_s(nw)*asf_s(nw) + dtauaer_f*ssa_f(nw)*asf_f(nw) ) & |
---|
270 | / ( ssa_s(nw)*dtauaer_s + ssa_f(nw)*dtauaer_f ) |
---|
271 | ELSE |
---|
272 | DHAZE_T(k,nw) = 0.D0 |
---|
273 | SSA_T(k,nw) = 1.0 |
---|
274 | ASF_T(k,nw) = 1.0 |
---|
275 | ENDIF |
---|
276 | |
---|
277 | ! Diagnostics for the haze : |
---|
278 | DIAG_OPTH(k,nw,1) = DHAZE_T(k,nw) ! dtau |
---|
279 | DIAG_OPTH(k,nw,2) = SSA_T(k,nw) ! wbar |
---|
280 | DIAG_OPTH(k,nw,3) = ASF_T(k,nw) ! gbar |
---|
281 | |
---|
282 | !--------------------- |
---|
283 | ! CLOUDS (Spherical) : |
---|
284 | !--------------------- |
---|
285 | IF (callclouds) THEN |
---|
286 | CTYPE = 0 |
---|
287 | m0ccn = pqmo(ilay,5) / 2.0 |
---|
288 | m3ccn = pqmo(ilay,6) / 2.0 |
---|
289 | m3cld = pqmo(ilay,6) / 2.0 |
---|
290 | |
---|
291 | ! Clear / Dark column method : |
---|
292 | !----------------------------- |
---|
293 | |
---|
294 | ! CCN's SSA : |
---|
295 | call get_haze_and_cloud_opacity(FTYPE,FTYPE,m0ccn,m3ccn,iw,dtau_ccn,ssa_ccn(nw),asf_ccn(nw)) |
---|
296 | |
---|
297 | ! Clear column (CCN, C2H2, C2H6, HCN) : |
---|
298 | IF (CDCOLUMN == 0) THEN |
---|
299 | DO iq = 2, nice |
---|
300 | m3cld = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0) |
---|
301 | ENDDO |
---|
302 | call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw)) |
---|
303 | |
---|
304 | ! Dark column (CCN, CH4, C2H2, C2H6, HCN) : |
---|
305 | ELSEIF (CDCOLUMN == 1) THEN |
---|
306 | DO iq = 1, nice |
---|
307 | m3cld = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0) |
---|
308 | ENDDO |
---|
309 | call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw)) |
---|
310 | |
---|
311 | ELSE |
---|
312 | WRITE(*,*) 'WARNING OPTCV.F90 : CDCOLUMN MUST BE 0 OR 1' |
---|
313 | WRITE(*,*) 'WE USE DARK COLUMN ...' |
---|
314 | DO iq = 1, nice |
---|
315 | m3cld = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0) |
---|
316 | ENDDO |
---|
317 | call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw)) |
---|
318 | ENDIF |
---|
319 | |
---|
320 | ! For small dropplets, opacity of nucleus dominates |
---|
321 | dtau_cld = (dtau_cld*m3ccn + dtau_cld*m3cld) / (m3ccn + m3cld) |
---|
322 | ssa_cld(nw) = (ssa_ccn(nw)*m3ccn + ssa_cld(nw)*m3cld) / (m3ccn + m3cld) |
---|
323 | |
---|
324 | ! Tuning of optical properties for clouds : |
---|
325 | dtau_cld = dtau_cld * (FCVIS * (1-ssa_cld(nw)) + ssa_cld(nw)) |
---|
326 | ssa_cld(nw) = ssa_cld(nw) / (FCVIS * (1-ssa_cld(nw)) + ssa_cld(nw)) |
---|
327 | |
---|
328 | ! Total of Haze + Clouds opacity (dtau), SSA (w) and ASF (COS) : |
---|
329 | DHAZE_T(k,nw) = dtauaer_s + dtauaer_f + dtau_cld |
---|
330 | IF (DHAZE_T(k,nw) .GT. 1.D-30) THEN |
---|
331 | 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 ) |
---|
332 | 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) ) & |
---|
333 | / ( ssa_s(nw)*dtauaer_s + ssa_f(nw)*dtauaer_f + ssa_cld(nw)*dtau_cld ) |
---|
334 | ELSE |
---|
335 | DHAZE_T(k,nw) = 0.D0 |
---|
336 | SSA_T(k,nw) = 1.0 |
---|
337 | ASF_T(k,nw) = 1.0 |
---|
338 | ENDIF |
---|
339 | |
---|
340 | ! Diagnostics for clouds : |
---|
341 | DIAG_OPTT(k,nw,1) = DHAZE_T(k,nw) ! dtau |
---|
342 | DIAG_OPTT(k,nw,2) = SSA_T(k,nw) ! wbar |
---|
343 | DIAG_OPTT(k,nw,3) = ASF_T(k,nw) ! gbar |
---|
344 | |
---|
345 | ELSE |
---|
346 | ! Diagnostics for clouds : |
---|
347 | DIAG_OPTT(k,nw,1) = 0.D0 ! dtau |
---|
348 | DIAG_OPTT(k,nw,2) = 1.0 ! wbar |
---|
349 | DIAG_OPTT(k,nw,3) = 1.0 ! gbar |
---|
350 | ENDIF |
---|
351 | |
---|
352 | ! Optics and microphysics no coupled : |
---|
353 | ELSE |
---|
354 | ! Call fixed vertical haze profile of extinction - same for all columns |
---|
355 | call disr_haze(dz(k),plev(k),wnov(nw),DHAZE_T(k,nw),SSA_T(k,nw),ASF_T(k,nw)) |
---|
356 | if (seashaze) DHAZE_T(k,nw) = DHAZE_T(k,nw)*seashazefact(k) |
---|
357 | ! Diagnostics for the haze : |
---|
358 | DIAG_OPTH(k,nw,1) = DHAZE_T(k,nw) ! dtau |
---|
359 | DIAG_OPTH(k,nw,2) = SSA_T(k,nw) ! wbar |
---|
360 | DIAG_OPTH(k,nw,3) = ASF_T(k,nw) ! gbar |
---|
361 | ! Diagnostics for clouds : |
---|
362 | DIAG_OPTT(k,nw,1) = 0.D0 ! dtau |
---|
363 | DIAG_OPTT(k,nw,2) = 1.0 ! wbar |
---|
364 | DIAG_OPTT(k,nw,3) = 1.0 ! gbar |
---|
365 | ENDIF ! ENDIF callmufi |
---|
366 | |
---|
367 | !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR |
---|
368 | ! but visible does not handle very well diffusion in first layer. |
---|
369 | ! The tauaero and tauray are thus set to 0 (a small value for rayleigh because the code crashes otherwise) |
---|
370 | ! in the 4 first semilayers in optcv, but not optci. |
---|
371 | ! This solves random variations of the sw heating at the model top. |
---|
372 | if (k<5) DHAZE_T(K,:) = 0.0 |
---|
373 | |
---|
374 | DRAYAER = TRAY(K,NW) |
---|
375 | ! DRAYAER is Tau RAYleigh scattering, plus AERosol opacity |
---|
376 | DRAYAER = DRAYAER + DHAZE_T(K,NW) ! Titan's aerosol |
---|
377 | |
---|
378 | DCONT = 0.0 ! continuum absorption |
---|
379 | |
---|
380 | if(continuum.and.(.not.graybody).and.callgasvis)then |
---|
381 | ! include continua if necessary |
---|
382 | wn_cont = dble(wnov(nw)) |
---|
383 | T_cont = dble(TMID(k)) |
---|
384 | do igas=1,ngasmx |
---|
385 | |
---|
386 | p_cont = dble(PMID(k)*scalep*gfrac(igas,ilay)) |
---|
387 | |
---|
388 | dtemp=0.0 |
---|
389 | if(igas.eq.igas_N2)then |
---|
390 | |
---|
391 | interm = indv(nw,igas,igas) |
---|
392 | ! call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
393 | indv(nw,igas,igas) = interm |
---|
394 | ! only goes to 500 cm^-1, so unless we're around a cold brown dwarf, this is irrelevant in the visible |
---|
395 | |
---|
396 | elseif(igas.eq.igas_H2)then |
---|
397 | |
---|
398 | ! first do self-induced absorption |
---|
399 | interm = indv(nw,igas,igas) |
---|
400 | call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
401 | indv(nw,igas,igas) = interm |
---|
402 | |
---|
403 | ! then cross-interactions with other gases |
---|
404 | do jgas=1,ngasmx |
---|
405 | p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay)) |
---|
406 | dtempc = 0.0 |
---|
407 | if(jgas.eq.igas_N2)then |
---|
408 | interm = indv(nw,igas,jgas) |
---|
409 | call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) |
---|
410 | indv(nw,igas,jgas) = interm |
---|
411 | ! should be irrelevant in the visible |
---|
412 | endif |
---|
413 | dtemp = dtemp + dtempc |
---|
414 | enddo |
---|
415 | |
---|
416 | elseif(igas.eq.igas_CH4)then |
---|
417 | |
---|
418 | ! first do self-induced absorption |
---|
419 | interm = indv(nw,igas,igas) |
---|
420 | call interpolateCH4CH4(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
421 | indv(nw,igas,igas) = interm |
---|
422 | |
---|
423 | ! then cross-interactions with other gases |
---|
424 | do jgas=1,ngasmx |
---|
425 | p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay)) |
---|
426 | dtempc = 0.0 |
---|
427 | if(jgas.eq.igas_N2)then |
---|
428 | interm = indv(nw,igas,jgas) |
---|
429 | call interpolateN2CH4(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) |
---|
430 | indv(nw,igas,jgas) = interm |
---|
431 | endif |
---|
432 | dtemp = dtemp + dtempc |
---|
433 | enddo |
---|
434 | |
---|
435 | endif |
---|
436 | |
---|
437 | DCONT = DCONT + dtemp |
---|
438 | |
---|
439 | enddo |
---|
440 | |
---|
441 | DCONT = DCONT*dz(k) |
---|
442 | |
---|
443 | endif |
---|
444 | |
---|
445 | do ng=1,L_NGAUSS-1 |
---|
446 | |
---|
447 | ! Now compute TAUGAS |
---|
448 | |
---|
449 | ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically |
---|
450 | ! the execution time of optci/v -> ~ factor 2 on the whole radiative |
---|
451 | ! transfer on the tested simulations ! |
---|
452 | |
---|
453 | if (corrk_recombin) then |
---|
454 | tmpk = GASV_RECOMB(MT(K):MT(K)+1,MP(K):MP(K)+1,NW,NG) |
---|
455 | else |
---|
456 | tmpk = GASV(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG) |
---|
457 | endif |
---|
458 | |
---|
459 | KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG) |
---|
460 | KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG) |
---|
461 | KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG) |
---|
462 | KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG) |
---|
463 | |
---|
464 | ! Interpolate the gaseous k-coefficients to the requested T,P values |
---|
465 | |
---|
466 | ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) + & |
---|
467 | LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4) |
---|
468 | |
---|
469 | |
---|
470 | TAUGAS = U(k)*ANS |
---|
471 | |
---|
472 | TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT |
---|
473 | DTAUKV(K,nw,ng) = TAUGAS & |
---|
474 | + DRAYAER & ! DRAYAER includes all scattering contributions |
---|
475 | + DCONT ! For parameterized continuum aborption |
---|
476 | end do |
---|
477 | |
---|
478 | ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS), |
---|
479 | ! which holds continuum opacity only |
---|
480 | |
---|
481 | NG = L_NGAUSS |
---|
482 | DTAUKV(K,nw,ng) = DRAYAER + DCONT ! Scattering + parameterized continuum absorption, including Titan's haze |
---|
483 | |
---|
484 | DIAG_OPTH(K,nw,4) = DRAYAER |
---|
485 | DIAG_OPTH(K,nw,5) = TAUGAS |
---|
486 | DIAG_OPTH(K,nw,6) = DCONT |
---|
487 | DIAG_OPTT(K,nw,4) = DRAYAER |
---|
488 | DIAG_OPTT(K,nw,5) = TAUGAS |
---|
489 | DIAG_OPTT(K,nw,6) = DCONT |
---|
490 | |
---|
491 | end do ! K = L_LEVELS |
---|
492 | end do ! nw = L_NSPECTV |
---|
493 | |
---|
494 | !======================================================================= |
---|
495 | ! Now the full treatment for the layers, where besides the opacity |
---|
496 | ! we need to calculate the scattering albedo and asymmetry factors |
---|
497 | ! ====================================================================== |
---|
498 | |
---|
499 | ! Haze scattering |
---|
500 | !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR |
---|
501 | ! but not in the visible |
---|
502 | ! The dhaze_s is thus set to 0 in the 4 first semilayers in optcv, but not optci. |
---|
503 | ! This solves random variations of the sw heating at the model top. |
---|
504 | DO NW=1,L_NSPECTV |
---|
505 | DHAZES_T(1:4,NW) = 0.d0 |
---|
506 | DO K=5,L_LEVELS |
---|
507 | DHAZES_T(K,NW) = DHAZE_T(K,NW) * SSA_T(K,NW) ! effect of scattering albedo on haze |
---|
508 | ENDDO |
---|
509 | ENDDO |
---|
510 | |
---|
511 | ! NW spectral loop |
---|
512 | DO NW=1,L_NSPECTV |
---|
513 | ! L vertical loop |
---|
514 | DO L=1,L_NLAYRAD-1 |
---|
515 | K = 2*L+1 |
---|
516 | atemp(L,NW) = ASF_T(K,NW)*DHAZES_T(K,NW) + ASF_T(K+1,NW)*DHAZES_T(K+1,NW) |
---|
517 | btemp(L,NW) = DHAZES_T(K,NW) + DHAZES_T(K+1,NW) |
---|
518 | ctemp(L,NW) = btemp(L,NW) + 0.9999*(TRAY(K,NW) + TRAY(K+1,NW)) ! JVO 2017 : does this 0.999 is really meaningful ? |
---|
519 | btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) + TRAY(K+1,NW) |
---|
520 | COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW) |
---|
521 | END DO |
---|
522 | |
---|
523 | ! Last level |
---|
524 | L = L_NLAYRAD |
---|
525 | K = 2*L+1 |
---|
526 | atemp(L,NW) = ASF_T(K,NW)*DHAZES_T(K,NW) |
---|
527 | btemp(L,NW) = DHAZES_T(K,NW) |
---|
528 | ctemp(L,NW) = btemp(L,NW) + 0.9999*TRAY(K,NW) ! JVO 2017 : does this 0.999 is really meaningful ? |
---|
529 | btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) |
---|
530 | COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW) |
---|
531 | END DO |
---|
532 | |
---|
533 | ! NG Gauss loop |
---|
534 | DO NG=1,L_NGAUSS |
---|
535 | ! NW spectral loop |
---|
536 | DO NW=1,L_NSPECTV |
---|
537 | ! L vertical loop |
---|
538 | DO L=1,L_NLAYRAD-1 |
---|
539 | K = 2*L+1 |
---|
540 | DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG) |
---|
541 | WBARV(L,nw,ng) = ctemp(L,NW) / DTAUV(L,nw,ng) |
---|
542 | END DO |
---|
543 | |
---|
544 | ! Last level |
---|
545 | L = L_NLAYRAD |
---|
546 | K = 2*L+1 |
---|
547 | DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) |
---|
548 | WBARV(L,NW,NG) = ctemp(L,NW) / DTAUV(L,NW,NG) |
---|
549 | END DO |
---|
550 | END DO |
---|
551 | |
---|
552 | ! Total extinction optical depths |
---|
553 | !DO NG=1,L_NGAUSS ! full gauss loop |
---|
554 | ! DO NW=1,L_NSPECTV |
---|
555 | ! TAUCUMV(1,NW,NG)=0.0D0 |
---|
556 | ! DO K=2,L_LEVELS |
---|
557 | ! TAUCUMV(K,NW,NG)=TAUCUMV(K-1,NW,NG)+DTAUKV(K,NW,NG) |
---|
558 | ! END DO |
---|
559 | ! DO L=1,L_NLAYRAD |
---|
560 | ! TAUV(L,NW,NG)=TAUCUMV(2*L,NW,NG) |
---|
561 | ! END DO |
---|
562 | ! TAUV(L,NW,NG)=TAUCUMV(2*L_NLAYRAD+1,NW,NG) |
---|
563 | ! END DO |
---|
564 | !END DO ! end full gauss loop |
---|
565 | |
---|
566 | TAUCUMV(:,:,:) = DTAUKV(:,:,:) |
---|
567 | DO L=1,L_NLAYRAD |
---|
568 | TAUV(L,:,:)=TAUCUMV(2*L,:,:) |
---|
569 | END DO |
---|
570 | TAUV(L,:,:)=TAUCUMV(2*L_NLAYRAD+1,:,:) |
---|
571 | |
---|
572 | if(firstcall) firstcall = .false. |
---|
573 | |
---|
574 | return |
---|
575 | |
---|
576 | |
---|
577 | end subroutine optcv |
---|