1 | SUBROUTINE OPTCV(DTAUV,TAUV,TAUCUMV,PLEV, & |
---|
2 | QXVAER,QSVAER,GVAER,WBARV,COSBV, & |
---|
3 | TAURAY,TAUAERO,TMID,PMID,TAUGSURF,QVAR,MUVAR) |
---|
4 | |
---|
5 | use radinc_h |
---|
6 | use radcommon_h, only: gasv, tlimit, wrefVAR, Cmk, tgasref, pfgasref,wnov,scalep,indv,glat_ig |
---|
7 | use gases_h |
---|
8 | use comcstfi_mod, only: g, r, mugaz |
---|
9 | use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple,callgasvis |
---|
10 | |
---|
11 | implicit none |
---|
12 | |
---|
13 | !================================================================== |
---|
14 | ! |
---|
15 | ! Purpose |
---|
16 | ! ------- |
---|
17 | ! Calculates shortwave optical constants at each level. |
---|
18 | ! |
---|
19 | ! Authors |
---|
20 | ! ------- |
---|
21 | ! Adapted from the NASA Ames code by R. Wordsworth (2009) |
---|
22 | ! |
---|
23 | !================================================================== |
---|
24 | ! |
---|
25 | ! THIS SUBROUTINE SETS THE OPTICAL CONSTANTS IN THE VISUAL |
---|
26 | ! IT CALCULATES FOR EACH LAYER, FOR EACH SPECTRAL INTERVAL IN THE VISUAL |
---|
27 | ! LAYER: WBAR, DTAU, COSBAR |
---|
28 | ! LEVEL: TAU |
---|
29 | ! |
---|
30 | ! TAUV(L,NW,NG) is the cumulative optical depth at the top of radiation code |
---|
31 | ! layer L. NW is spectral wavelength interval, ng the Gauss point index. |
---|
32 | ! |
---|
33 | ! TLEV(L) - Temperature at the layer boundary |
---|
34 | ! PLEV(L) - Pressure at the layer boundary (i.e. level) |
---|
35 | ! GASV(NT,NPS,NW,NG) - Visible k-coefficients |
---|
36 | ! |
---|
37 | !------------------------------------------------------------------- |
---|
38 | |
---|
39 | |
---|
40 | real*8 DTAUV(L_NLAYRAD,L_NSPECTV,L_NGAUSS) |
---|
41 | real*8 DTAUKV(L_LEVELS,L_NSPECTV,L_NGAUSS) |
---|
42 | real*8 TAUV(L_NLEVRAD,L_NSPECTV,L_NGAUSS) |
---|
43 | real*8 TAUCUMV(L_LEVELS,L_NSPECTV,L_NGAUSS) |
---|
44 | real*8 PLEV(L_LEVELS) |
---|
45 | real*8 TMID(L_LEVELS), PMID(L_LEVELS) |
---|
46 | real*8 COSBV(L_NLAYRAD,L_NSPECTV,L_NGAUSS) |
---|
47 | real*8 WBARV(L_NLAYRAD,L_NSPECTV,L_NGAUSS) |
---|
48 | |
---|
49 | ! for aerosols |
---|
50 | real*8 QXVAER(L_LEVELS,L_NSPECTV,NAERKIND) |
---|
51 | real*8 QSVAER(L_LEVELS,L_NSPECTV,NAERKIND) |
---|
52 | real*8 GVAER(L_LEVELS,L_NSPECTV,NAERKIND) |
---|
53 | real*8 TAUAERO(L_LEVELS,NAERKIND) |
---|
54 | real*8 TAUAEROLK(L_LEVELS,L_NSPECTV,NAERKIND) |
---|
55 | real*8 TAEROS(L_LEVELS,L_NSPECTV,NAERKIND) |
---|
56 | |
---|
57 | integer L, NW, NG, K, LK, IAER |
---|
58 | integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS) |
---|
59 | real*8 ANS, TAUGAS |
---|
60 | real*8 TAURAY(L_NSPECTV) |
---|
61 | real*8 TRAY(L_LEVELS,L_NSPECTV) |
---|
62 | real*8 DPR(L_LEVELS), U(L_LEVELS) |
---|
63 | real*8 LCOEF(4), LKCOEF(L_LEVELS,4) |
---|
64 | |
---|
65 | real*8 taugsurf(L_NSPECTV,L_NGAUSS-1) |
---|
66 | real*8 DCONT,DAERO |
---|
67 | real*8 DRAYAER |
---|
68 | double precision wn_cont, p_cont, p_air, T_cont, dtemp, dtempc |
---|
69 | double precision p_cross |
---|
70 | |
---|
71 | ! variable species mixing ratio variables |
---|
72 | real*8 QVAR(L_LEVELS), WRATIO(L_LEVELS), MUVAR(L_LEVELS) |
---|
73 | real*8 KCOEF(4) |
---|
74 | integer NVAR(L_LEVELS) |
---|
75 | |
---|
76 | ! temporary variables to reduce memory access time to gasv |
---|
77 | real*8 tmpk(2,2) |
---|
78 | real*8 tmpkvar(2,2,2) |
---|
79 | |
---|
80 | ! temporary variables for multiple aerosol calculation |
---|
81 | real*8 atemp(L_NLAYRAD,L_NSPECTV) |
---|
82 | real*8 btemp(L_NLAYRAD,L_NSPECTV) |
---|
83 | real*8 ctemp(L_NLAYRAD,L_NSPECTV) |
---|
84 | |
---|
85 | ! variables for k in units m^-1 |
---|
86 | real*8 dz(L_LEVELS) |
---|
87 | |
---|
88 | integer igas, jgas |
---|
89 | |
---|
90 | integer interm |
---|
91 | |
---|
92 | !! AS: to save time in computing continuum (see bilinearbig) |
---|
93 | IF (.not.ALLOCATED(indv)) THEN |
---|
94 | ALLOCATE(indv(L_NSPECTV,ngasmx,ngasmx)) |
---|
95 | indv = -9999 ! this initial value means "to be calculated" |
---|
96 | ENDIF |
---|
97 | |
---|
98 | !======================================================================= |
---|
99 | ! Determine the total gas opacity throughout the column, for each |
---|
100 | ! spectral interval, NW, and each Gauss point, NG. |
---|
101 | ! Calculate the continuum opacities, i.e., those that do not depend on |
---|
102 | ! NG, the Gauss index. |
---|
103 | |
---|
104 | taugsurf(:,:) = 0.0 |
---|
105 | dpr(:) = 0.0 |
---|
106 | lkcoef(:,:) = 0.0 |
---|
107 | |
---|
108 | do K=2,L_LEVELS |
---|
109 | DPR(k) = PLEV(K)-PLEV(K-1) |
---|
110 | |
---|
111 | ! if we have continuum opacities, we need dz |
---|
112 | if(kastprof)then |
---|
113 | dz(k) = dpr(k)*(1000.0d0*8.3145d0/muvar(k))*TMID(K)/(g*PMID(K)) |
---|
114 | U(k) = Cmk*DPR(k)*mugaz/muvar(k) |
---|
115 | else |
---|
116 | dz(k) = dpr(k)*R*TMID(K)/(glat_ig*PMID(K))*mugaz/muvar(k) |
---|
117 | U(k) = Cmk*DPR(k)*mugaz/muvar(k) ! only Cmk line in optci.F |
---|
118 | !JL13 the mugaz/muvar factor takes into account water meanmolecular weight if water is present |
---|
119 | endif |
---|
120 | |
---|
121 | call tpindex(PMID(K),TMID(K),QVAR(K),pfgasref,tgasref,WREFVAR, & |
---|
122 | LCOEF,MT(K),MP(K),NVAR(K),WRATIO(K)) |
---|
123 | |
---|
124 | do LK=1,4 |
---|
125 | LKCOEF(K,LK) = LCOEF(LK) |
---|
126 | end do |
---|
127 | end do ! levels |
---|
128 | |
---|
129 | ! Spectral dependance of aerosol absorption |
---|
130 | do iaer=1,naerkind |
---|
131 | do NW=1,L_NSPECTV |
---|
132 | do K=2,L_LEVELS |
---|
133 | TAEROS(K,NW,IAER) = TAUAERO(K,IAER) * QXVAER(K,NW,IAER) |
---|
134 | end do ! levels |
---|
135 | end do |
---|
136 | end do |
---|
137 | |
---|
138 | ! Rayleigh scattering |
---|
139 | do NW=1,L_NSPECTV |
---|
140 | do K=2,L_LEVELS |
---|
141 | TRAY(K,NW) = TAURAY(NW) * DPR(K) |
---|
142 | end do ! levels |
---|
143 | end do |
---|
144 | |
---|
145 | ! we ignore K=1... |
---|
146 | do K=2,L_LEVELS |
---|
147 | |
---|
148 | do NW=1,L_NSPECTV |
---|
149 | |
---|
150 | DRAYAER = TRAY(K,NW) |
---|
151 | ! DRAYAER is Tau RAYleigh scattering, plus AERosol opacity |
---|
152 | do iaer=1,naerkind |
---|
153 | DRAYAER = DRAYAER + TAEROS(K,NW,IAER) |
---|
154 | end do |
---|
155 | |
---|
156 | DCONT = 0.0 ! continuum absorption |
---|
157 | |
---|
158 | if(continuum.and.(.not.graybody).and.callgasvis)then |
---|
159 | ! include continua if necessary |
---|
160 | wn_cont = dble(wnov(nw)) |
---|
161 | T_cont = dble(TMID(k)) |
---|
162 | do igas=1,ngasmx |
---|
163 | |
---|
164 | if(gfrac(igas).eq.-1)then ! variable |
---|
165 | p_cont = dble(PMID(k)*scalep*QVAR(k)) ! qvar = mol/mol |
---|
166 | else |
---|
167 | p_cont = dble(PMID(k)*scalep*gfrac(igas)*(1.-QVAR(k))) |
---|
168 | endif |
---|
169 | |
---|
170 | dtemp=0.0 |
---|
171 | if(igas.eq.igas_N2)then |
---|
172 | |
---|
173 | interm = indv(nw,igas,igas) |
---|
174 | ! call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
175 | indv(nw,igas,igas) = interm |
---|
176 | ! only goes to 500 cm^-1, so unless we're around a cold brown dwarf, this is irrelevant in the visible |
---|
177 | |
---|
178 | elseif(igas.eq.igas_H2)then |
---|
179 | |
---|
180 | ! first do self-induced absorption |
---|
181 | interm = indv(nw,igas,igas) |
---|
182 | call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) |
---|
183 | indv(nw,igas,igas) = interm |
---|
184 | |
---|
185 | ! then cross-interactions with other gases |
---|
186 | do jgas=1,ngasmx |
---|
187 | p_cross = dble(PMID(k)*scalep*gfrac(jgas)*(1.-QVAR(k))) |
---|
188 | dtempc = 0.0 |
---|
189 | if(jgas.eq.igas_N2)then |
---|
190 | interm = indv(nw,igas,jgas) |
---|
191 | call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) |
---|
192 | indv(nw,igas,jgas) = interm |
---|
193 | ! should be irrelevant in the visible |
---|
194 | elseif(jgas.eq.igas_He)then |
---|
195 | interm = indv(nw,igas,jgas) |
---|
196 | call interpolateH2He(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) |
---|
197 | indv(nw,igas,jgas) = interm |
---|
198 | endif |
---|
199 | dtemp = dtemp + dtempc |
---|
200 | enddo |
---|
201 | |
---|
202 | elseif(igas.eq.igas_H2O.and.T_cont.gt.200.0)then |
---|
203 | |
---|
204 | p_air = dble(PMID(k)*scalep) - p_cont ! note assumes background is air! |
---|
205 | if(H2Ocont_simple)then |
---|
206 | call interpolateH2Ocont_PPC(wn_cont,T_cont,p_cont,p_air,dtemp,.false.) |
---|
207 | else |
---|
208 | interm = indv(nw,igas,igas) |
---|
209 | call interpolateH2Ocont_CKD(wn_cont,T_cont,p_cont,p_air,dtemp,.false.,interm) |
---|
210 | indv(nw,igas,igas) = interm |
---|
211 | endif |
---|
212 | |
---|
213 | endif |
---|
214 | |
---|
215 | DCONT = DCONT + dtemp |
---|
216 | |
---|
217 | enddo |
---|
218 | |
---|
219 | DCONT = DCONT*dz(k) |
---|
220 | |
---|
221 | endif |
---|
222 | |
---|
223 | do ng=1,L_NGAUSS-1 |
---|
224 | |
---|
225 | ! Now compute TAUGAS |
---|
226 | |
---|
227 | ! Interpolate between water mixing ratios |
---|
228 | ! WRATIO = 0.0 if the requested water amount is equal to, or outside the |
---|
229 | ! the water data range |
---|
230 | |
---|
231 | if(L_REFVAR.eq.1)then ! added by RW for special no variable case |
---|
232 | |
---|
233 | ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically |
---|
234 | ! the execution time of optci/v -> ~ factor 2 on the whole radiative |
---|
235 | ! transfer on the tested simulations ! |
---|
236 | |
---|
237 | tmpk = GASV(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG) |
---|
238 | |
---|
239 | KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG) |
---|
240 | KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG) |
---|
241 | KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG) |
---|
242 | KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG) |
---|
243 | |
---|
244 | else |
---|
245 | |
---|
246 | tmpkvar = GASV(MT(K):MT(K)+1,MP(K):MP(K)+1,NVAR(K):NVAR(K)+1,NW,NG) |
---|
247 | |
---|
248 | KCOEF(1) = tmpkvar(1,1,1) + WRATIO(K) * & |
---|
249 | ( tmpkvar(1,1,2)-tmpkvar(1,1,1) ) |
---|
250 | |
---|
251 | KCOEF(2) = tmpkvar(1,2,1) + WRATIO(K) * & |
---|
252 | ( tmpkvar(1,2,2)-tmpkvar(1,2,1) ) |
---|
253 | |
---|
254 | KCOEF(3) = tmpkvar(2,2,1) + WRATIO(K) * & |
---|
255 | ( tmpkvar(2,2,2)-tmpkvar(2,2,1) ) |
---|
256 | |
---|
257 | KCOEF(4) = tmpkvar(2,1,1) + WRATIO(K) * & |
---|
258 | ( tmpkvar(2,1,2)-tmpkvar(2,1,1) ) |
---|
259 | |
---|
260 | |
---|
261 | endif |
---|
262 | |
---|
263 | ! Interpolate the gaseous k-coefficients to the requested T,P values |
---|
264 | |
---|
265 | ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) + & |
---|
266 | LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4) |
---|
267 | |
---|
268 | TAUGAS = U(k)*ANS |
---|
269 | |
---|
270 | TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT |
---|
271 | DTAUKV(K,nw,ng) = TAUGAS & |
---|
272 | + DRAYAER & ! DRAYAER includes all scattering contributions |
---|
273 | + DCONT ! For parameterized continuum aborption |
---|
274 | |
---|
275 | end do |
---|
276 | |
---|
277 | ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS), |
---|
278 | ! which holds continuum opacity only |
---|
279 | |
---|
280 | NG = L_NGAUSS |
---|
281 | DTAUKV(K,nw,ng) = DRAYAER + DCONT ! Scattering + parameterized continuum absorption |
---|
282 | |
---|
283 | end do |
---|
284 | end do |
---|
285 | |
---|
286 | |
---|
287 | !======================================================================= |
---|
288 | ! Now the full treatment for the layers, where besides the opacity |
---|
289 | ! we need to calculate the scattering albedo and asymmetry factors |
---|
290 | |
---|
291 | do iaer=1,naerkind |
---|
292 | DO NW=1,L_NSPECTV |
---|
293 | DO K=2,L_LEVELS |
---|
294 | TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER) * QSVAER(K,NW,IAER) ! effect of scattering albedo |
---|
295 | ENDDO |
---|
296 | ENDDO |
---|
297 | end do |
---|
298 | |
---|
299 | DO NW=1,L_NSPECTV |
---|
300 | DO L=1,L_NLAYRAD-1 |
---|
301 | K = 2*L+1 |
---|
302 | atemp(L,NW) = SUM(GVAER(K,NW,1:naerkind) * TAUAEROLK(K,NW,1:naerkind))+SUM(GVAER(K+1,NW,1:naerkind) * TAUAEROLK(K+1,NW,1:naerkind)) |
---|
303 | btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) + SUM(TAUAEROLK(K+1,NW,1:naerkind)) |
---|
304 | ctemp(L,NW) = btemp(L,NW) + 0.9999*(TRAY(K,NW) + TRAY(K+1,NW)) ! JVO 2017 : does this 0.999 is really meaningful ? |
---|
305 | btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) + TRAY(K+1,NW) |
---|
306 | COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW) |
---|
307 | END DO ! L vertical loop |
---|
308 | |
---|
309 | ! Last level |
---|
310 | L = L_NLAYRAD |
---|
311 | K = 2*L+1 |
---|
312 | atemp(L,NW) = SUM(GVAER(K,NW,1:naerkind) * TAUAEROLK(K,NW,1:naerkind)) |
---|
313 | btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) |
---|
314 | ctemp(L,NW) = btemp(L,NW) + 0.9999*TRAY(K,NW) ! JVO 2017 : does this 0.999 is really meaningful ? |
---|
315 | btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) |
---|
316 | COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW) |
---|
317 | |
---|
318 | |
---|
319 | END DO ! NW spectral loop |
---|
320 | |
---|
321 | DO NG=1,L_NGAUSS |
---|
322 | DO NW=1,L_NSPECTV |
---|
323 | DO L=1,L_NLAYRAD-1 |
---|
324 | |
---|
325 | K = 2*L+1 |
---|
326 | DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG) |
---|
327 | WBARV(L,nw,ng) = ctemp(L,NW) / DTAUV(L,nw,ng) |
---|
328 | |
---|
329 | END DO ! L vertical loop |
---|
330 | |
---|
331 | ! Last level |
---|
332 | |
---|
333 | L = L_NLAYRAD |
---|
334 | K = 2*L+1 |
---|
335 | DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) |
---|
336 | |
---|
337 | WBARV(L,NW,NG) = ctemp(L,NW) / DTAUV(L,NW,NG) |
---|
338 | |
---|
339 | END DO ! NW spectral loop |
---|
340 | END DO ! NG Gauss loop |
---|
341 | |
---|
342 | ! Total extinction optical depths |
---|
343 | |
---|
344 | DO NG=1,L_NGAUSS ! full gauss loop |
---|
345 | DO NW=1,L_NSPECTV |
---|
346 | TAUV(1,NW,NG)=0.0D0 |
---|
347 | DO L=1,L_NLAYRAD |
---|
348 | TAUV(L+1,NW,NG)=TAUV(L,NW,NG)+DTAUV(L,NW,NG) |
---|
349 | END DO |
---|
350 | |
---|
351 | TAUCUMV(1,NW,NG)=0.0D0 |
---|
352 | DO K=2,L_LEVELS |
---|
353 | TAUCUMV(K,NW,NG)=TAUCUMV(K-1,NW,NG)+DTAUKV(K,NW,NG) |
---|
354 | END DO |
---|
355 | END DO |
---|
356 | END DO ! end full gauss loop |
---|
357 | |
---|
358 | |
---|
359 | return |
---|
360 | |
---|
361 | |
---|
362 | end subroutine optcv |
---|