Line data Source code
1 : ! path: $Source: /storm/rc1/cvsroot/rc/rrtmg_sw/src/rrtmg_sw_init.f90,v $
2 : ! author: $Author: mike $
3 : ! revision: $Revision: 1.2 $
4 : ! created: $Date: 2007/08/23 20:40:13 $
5 :
6 : module rrtmg_sw_init
7 :
8 : ! --------------------------------------------------------------------------
9 : ! | |
10 : ! | Copyright 2002-2007, Atmospheric & Environmental Research, Inc. (AER). |
11 : ! | This software may be used, copied, or redistributed as long as it is |
12 : ! | not sold and this copyright notice is reproduced on each copy made. |
13 : ! | This model is provided as is without any express or implied warranties. |
14 : ! | (http://www.rtweb.aer.com/) |
15 : ! | |
16 : ! --------------------------------------------------------------------------
17 :
18 : ! ------- Modules -------
19 :
20 : use shr_kind_mod, only: r8 => shr_kind_r8
21 :
22 : ! use parkind, only : jpim, jprb
23 : use rrsw_wvn
24 : use rrtmg_sw_setcoef, only: swatmref
25 :
26 : implicit none
27 :
28 : contains
29 :
30 : ! **************************************************************************
31 1024 : subroutine rrtmg_sw_ini
32 : ! **************************************************************************
33 : !
34 : ! Original version: Michael J. Iacono; February, 2004
35 : ! Revision for F90 formatting: M. J. Iacono, July, 2006
36 : !
37 : ! This subroutine performs calculations necessary for the initialization
38 : ! of the shortwave model. Lookup tables are computed for use in the SW
39 : ! radiative transfer, and input absorption coefficient data for each
40 : ! spectral band are reduced from 224 g-point intervals to 112.
41 : ! **************************************************************************
42 :
43 : use parrrsw, only : mg, nbndsw, ngptsw
44 : use rrsw_tbl, only: ntbl, tblint, pade, bpade, tau_tbl, exp_tbl
45 :
46 : ! ------- Local -------
47 :
48 : integer :: ibnd, igc, ig, ind, ipr
49 : integer :: igcsm, iprsm
50 : integer :: itr
51 :
52 : real(kind=r8) :: wtsum, wtsm(mg)
53 : real(kind=r8) :: tfn
54 :
55 : ! ------- Definitions -------
56 : ! Arrays for 10000-point look-up tables:
57 : ! TAU_TBL Clear-sky optical depth
58 : ! EXP_TBL Exponential lookup table for transmittance
59 : ! PADE Pade approximation constant (= 0.278)
60 : ! BPADE Inverse of the Pade approximation constant
61 : !
62 :
63 : ! Initialize model data
64 1024 : call swdatinit
65 1024 : call swcmbdat ! g-point interval reduction data
66 1024 : call swatmref ! reference MLS profile
67 1024 : call sw_kgb16 ! molecular absorption coefficients
68 1024 : call sw_kgb17
69 1024 : call sw_kgb18
70 1024 : call sw_kgb19
71 1024 : call sw_kgb20
72 1024 : call sw_kgb21
73 1024 : call sw_kgb22
74 1024 : call sw_kgb23
75 1024 : call sw_kgb24
76 1024 : call sw_kgb25
77 1024 : call sw_kgb26
78 1024 : call sw_kgb27
79 1024 : call sw_kgb28
80 1024 : call sw_kgb29
81 :
82 : ! Define exponential lookup tables for transmittance. Tau is
83 : ! computed as a function of the tau transition function, and transmittance
84 : ! is calculated as a function of tau. All tables are computed at intervals
85 : ! of 0.0001. The inverse of the constant used in the Pade approximation to
86 : ! the tau transition function is set to bpade.
87 :
88 1024 : exp_tbl(0) = 1.0_r8
89 1024 : exp_tbl(ntbl) = 0.0_r8
90 1024 : bpade = 1.0_r8 / pade
91 10240000 : do itr = 1, ntbl-1
92 10238976 : tfn = float(itr) / float(ntbl)
93 10238976 : tau_tbl = bpade * tfn / (1._r8 - tfn)
94 10240000 : exp_tbl(itr) = exp(-tau_tbl)
95 : enddo
96 :
97 : ! Perform g-point reduction from 16 per band (224 total points) to
98 : ! a band dependent number (112 total points) for all absorption
99 : ! coefficient input data and Planck fraction input data.
100 : ! Compute relative weighting for new g-point combinations.
101 :
102 1024 : igcsm = 0
103 15360 : do ibnd = 1,nbndsw
104 14336 : iprsm = 0
105 15360 : if (ngc(ibnd).lt.mg) then
106 129024 : do igc = 1,ngc(ibnd)
107 114688 : igcsm = igcsm + 1
108 114688 : wtsum = 0.
109 344064 : do ipr = 1, ngn(igcsm)
110 229376 : iprsm = iprsm + 1
111 344064 : wtsum = wtsum + wt(iprsm)
112 : enddo
113 129024 : wtsm(igc) = wtsum
114 : enddo
115 243712 : do ig = 1, ng(ibnd+15)
116 229376 : ind = (ibnd-1)*mg + ig
117 243712 : rwgt(ind) = wt(ig)/wtsm(ngm(ind))
118 : enddo
119 : else
120 0 : do ig = 1, ng(ibnd+15)
121 0 : igcsm = igcsm + 1
122 0 : ind = (ibnd-1)*mg + ig
123 0 : rwgt(ind) = 1.0_r8
124 : enddo
125 : endif
126 : enddo
127 :
128 : ! Reduce g-points for absorption coefficient data in each LW spectral band.
129 :
130 1024 : call cmbgb16s
131 1024 : call cmbgb17
132 1024 : call cmbgb18
133 1024 : call cmbgb19
134 1024 : call cmbgb20
135 1024 : call cmbgb21
136 1024 : call cmbgb22
137 1024 : call cmbgb23
138 1024 : call cmbgb24
139 1024 : call cmbgb25
140 1024 : call cmbgb26
141 1024 : call cmbgb27
142 1024 : call cmbgb28
143 1024 : call cmbgb29
144 :
145 1024 : end subroutine rrtmg_sw_ini
146 :
147 : !***************************************************************************
148 1024 : subroutine swdatinit
149 : !***************************************************************************
150 :
151 : ! --------- Modules ----------
152 :
153 : use rrsw_con, only: heatfac, grav, planck, boltz, &
154 : clight, avogad, alosmt, gascon, radcn1, radcn2
155 : use rrsw_wvn, only: ng, nspa, nspb, wavenum1, wavenum2, delwave
156 : use shr_const_mod, only: shr_const_avogad
157 : use physconst, only: cday, gravit, cpair
158 :
159 : save
160 :
161 : ! Shortwave spectral band limits (wavenumbers)
162 : wavenum1(:) = (/2600._r8, 3250._r8, 4000._r8, 4650._r8, 5150._r8, 6150._r8, 7700._r8, &
163 1024 : 8050._r8,12850._r8,16000._r8,22650._r8,29000._r8,38000._r8, 820._r8/)
164 : wavenum2(:) = (/3250._r8, 4000._r8, 4650._r8, 5150._r8, 6150._r8, 7700._r8, 8050._r8, &
165 1024 : 12850._r8,16000._r8,22650._r8,29000._r8,38000._r8,50000._r8, 2600._r8/)
166 : delwave(:) = (/ 650._r8, 750._r8, 650._r8, 500._r8, 1000._r8, 1550._r8, 350._r8, &
167 1024 : 4800._r8, 3150._r8, 6650._r8, 6350._r8, 9000._r8,12000._r8, 1780._r8/)
168 :
169 : ! Spectral band information
170 1024 : ng(:) = (/16,16,16,16,16,16,16,16,16,16,16,16,16,16/)
171 1024 : nspa(:) = (/9,9,9,9,1,9,9,1,9,1,0,1,9,1/)
172 1024 : nspb(:) = (/1,5,1,1,1,5,1,0,1,0,0,1,5,1/)
173 :
174 : ! Use constants set in CAM for consistency
175 1024 : grav = gravit
176 1024 : avogad = shr_const_avogad * 1.e-3_r8
177 :
178 : ! Heatfac is the factor by which one must multiply delta-flux/
179 : ! delta-pressure, with flux in w/m-2 and pressure in mbar, to get
180 : ! the heating rate in units of degrees/day. It is equal to
181 : ! (g)x(#sec/day)x(1e-5)/(specific heat of air at const. p)
182 : ! = (9.8066)(86400)(1e-5)/(1.004)
183 : ! heatfac = 8.4391_r8
184 :
185 : ! Calculate heatfac directly from CAM constants:
186 1024 : heatfac = grav * cday * 1.e-5_r8 / (cpair * 1.e-3_r8)
187 :
188 : ! Constants from NIST 01/11/2002
189 :
190 : ! grav = 9.8066_r8
191 1024 : planck = 6.62606876e-27_r8
192 1024 : boltz = 1.3806503e-16_r8
193 1024 : clight = 2.99792458e+10_r8
194 : ! avogad = 6.02214199e+23_r8
195 1024 : alosmt = 2.6867775e+19_r8
196 1024 : gascon = 8.31447200e+07_r8
197 1024 : radcn1 = 1.191042722e-12_r8
198 1024 : radcn2 = 1.4387752_r8
199 :
200 : !
201 : ! units are generally cgs
202 : !
203 : ! The first and second radiation constants are taken from NIST.
204 : ! They were previously obtained from the relations:
205 : ! radcn1 = 2.*planck*clight*clight*1.e-07
206 : ! radcn2 = planck*clight/boltz
207 :
208 1024 : end subroutine swdatinit
209 :
210 : !***************************************************************************
211 1024 : subroutine swcmbdat
212 : !***************************************************************************
213 :
214 : use rrsw_wvn, only: ngc, ngs, ngn, ngb, ngm, wt
215 :
216 : save
217 :
218 : ! ------- Definitions -------
219 : ! Arrays for the g-point reduction from 224 to 112 for the 16 LW bands:
220 : ! This mapping from 224 to 112 points has been carefully selected to
221 : ! minimize the effect on the resulting fluxes and cooling rates, and
222 : ! caution should be used if the mapping is modified. The full 224
223 : ! g-point set can be restored with ngpt=224, ngc=16*16, ngn=224*1., etc.
224 : ! ngpt The total number of new g-points
225 : ! ngc The number of new g-points in each band
226 : ! ngs The cumulative sum of new g-points for each band
227 : ! ngm The index of each new g-point relative to the original
228 : ! 16 g-points for each band.
229 : ! ngn The number of original g-points that are combined to make
230 : ! each new g-point in each band.
231 : ! ngb The band index for each new g-point.
232 : ! wt RRTM weights for 16 g-points.
233 :
234 : ! Use this set for 112 quadrature point (g-point) model
235 : ! ------- Data statements -------
236 1024 : ngc(:) = (/ 6,12, 8, 8,10,10, 2,10, 8, 6, 6, 8, 6,12 /)
237 1024 : ngs(:) = (/ 6,18,26,34,44,54,56,66,74,80,86,94,100,112 /)
238 : ngm(:) = (/ 1,1,2,2,3,3,4,4,5,5,5,5,6,6,6,6, & ! band 16
239 : 1,2,3,4,5,6,6,7,8,8,9,10,10,11,12,12, & ! band 17
240 : 1,2,3,4,5,5,6,6,7,7,7,7,8,8,8,8, & ! band 18
241 : 1,2,3,4,5,5,6,6,7,7,7,7,8,8,8,8, & ! band 19
242 : 1,2,3,4,5,6,7,8,9,9,10,10,10,10,10,10, & ! band 20
243 : 1,2,3,4,5,6,7,8,9,9,10,10,10,10,10,10, & ! band 21
244 : 1,1,1,1,1,1,1,1,2,2,2,2,2,2,2,2, & ! band 22
245 : 1,1,2,2,3,4,5,6,7,8,9,9,10,10,10,10, & ! band 23
246 : 1,1,2,2,3,3,4,4,5,5,6,6,7,7,8,8, & ! band 24
247 : 1,2,3,3,4,4,5,5,5,5,6,6,6,6,6,6, & ! band 25
248 : 1,2,3,3,4,4,5,5,5,5,6,6,6,6,6,6, & ! band 26
249 : 1,2,3,4,5,6,7,7,7,7,8,8,8,8,8,8, & ! band 27
250 : 1,2,3,3,4,4,5,5,5,5,6,6,6,6,6,6, & ! band 28
251 1024 : 1,2,3,4,5,5,6,6,7,7,8,8,9,10,11,12 /) ! band 29
252 : ngn(:) = (/ 2,2,2,2,4,4, & ! band 16
253 : 1,1,1,1,1,2,1,2,1,2,1,2, & ! band 17
254 : 1,1,1,1,2,2,4,4, & ! band 18
255 : 1,1,1,1,2,2,4,4, & ! band 19
256 : 1,1,1,1,1,1,1,1,2,6, & ! band 20
257 : 1,1,1,1,1,1,1,1,2,6, & ! band 21
258 : 8,8, & ! band 22
259 : 2,2,1,1,1,1,1,1,2,4, & ! band 23
260 : 2,2,2,2,2,2,2,2, & ! band 24
261 : 1,1,2,2,4,6, & ! band 25
262 : 1,1,2,2,4,6, & ! band 26
263 : 1,1,1,1,1,1,4,6, & ! band 27
264 : 1,1,2,2,4,6, & ! band 28
265 1024 : 1,1,1,1,2,2,2,2,1,1,1,1 /) ! band 29
266 : ngb(:) = (/ 16,16,16,16,16,16, & ! band 16
267 : 17,17,17,17,17,17,17,17,17,17,17,17, & ! band 17
268 : 18,18,18,18,18,18,18,18, & ! band 18
269 : 19,19,19,19,19,19,19,19, & ! band 19
270 : 20,20,20,20,20,20,20,20,20,20, & ! band 20
271 : 21,21,21,21,21,21,21,21,21,21, & ! band 21
272 : 22,22, & ! band 22
273 : 23,23,23,23,23,23,23,23,23,23, & ! band 23
274 : 24,24,24,24,24,24,24,24, & ! band 24
275 : 25,25,25,25,25,25, & ! band 25
276 : 26,26,26,26,26,26, & ! band 26
277 : 27,27,27,27,27,27,27,27, & ! band 27
278 : 28,28,28,28,28,28, & ! band 28
279 1024 : 29,29,29,29,29,29,29,29,29,29,29,29 /) ! band 29
280 :
281 : ! Use this set for full 224 quadrature point (g-point) model
282 : ! ------- Data statements -------
283 : ! ngc(:) = (/ 16,16,16,16,16,16,16,16,16,16,16,16,16,16 /)
284 : ! ngs(:) = (/ 16,32,48,64,80,96,112,128,144,160,176,192,208,224 /)
285 : ! ngm(:) = (/ 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 16
286 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 17
287 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 18
288 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 19
289 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 20
290 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 21
291 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 22
292 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 23
293 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 24
294 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 25
295 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 26
296 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 27
297 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 28
298 : ! 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16 /) ! band 29
299 : ! ngn(:) = (/ 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 16
300 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 17
301 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 18
302 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 19
303 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 20
304 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 21
305 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 22
306 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 23
307 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 24
308 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 25
309 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 26
310 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 27
311 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 28
312 : ! 1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1 /) ! band 29
313 : ! ngb(:) = (/ 16,16,16,16,16,16,16,16,16,16,16,16,16,16,16,16, & ! band 16
314 : ! 17,17,17,17,17,17,17,17,17,17,17,17,17,17,17,17, & ! band 17
315 : ! 18,18,18,18,18,18,18,18,18,18,18,18,18,18,18,18, & ! band 18
316 : ! 19,19,19,19,19,19,19,19,19,19,19,19,19,19,19,19, & ! band 19
317 : ! 20,20,20,20,20,20,20,20,20,20,20,20,20,20,20,20, & ! band 20
318 : ! 21,21,21,21,21,21,21,21,21,21,21,21,21,21,21,21, & ! band 21
319 : ! 22,22,22,22,22,22,22,22,22,22,22,22,22,22,22,22, & ! band 22
320 : ! 23,23,23,23,23,23,23,23,23,23,23,23,23,23,23,23, & ! band 23
321 : ! 24,24,24,24,24,24,24,24,24,24,24,24,24,24,24,24, & ! band 24
322 : ! 25,25,25,25,25,25,25,25,25,25,25,25,25,25,25,25, & ! band 25
323 : ! 26,26,26,26,26,26,26,26,26,26,26,26,26,26,26,26, & ! band 26
324 : ! 27,27,27,27,27,27,27,27,27,27,27,27,27,27,27,27, & ! band 27
325 : ! 28,28,28,28,28,28,28,28,28,28,28,28,28,28,28,28, & ! band 28
326 : ! 29,29,29,29,29,29,29,29,29,29,29,29,29,29,29,29 /) ! band 29
327 :
328 :
329 : wt(:) = (/ 0.1527534276_r8, 0.1491729617_r8, 0.1420961469_r8, &
330 : 0.1316886544_r8, 0.1181945205_r8, 0.1019300893_r8, &
331 : 0.0832767040_r8, 0.0626720116_r8, 0.0424925000_r8, &
332 : 0.0046269894_r8, 0.0038279891_r8, 0.0030260086_r8, &
333 : 0.0022199750_r8, 0.0014140010_r8, 0.0005330000_r8, &
334 1024 : 0.0000750000_r8 /)
335 :
336 1024 : end subroutine swcmbdat
337 :
338 : !***************************************************************************
339 1024 : subroutine cmbgb16s
340 : !***************************************************************************
341 : !
342 : ! Original version: MJIacono; July 1998
343 : ! Revision for RRTM_SW: MJIacono; November 2002
344 : ! Revision for RRTMG_SW: MJIacono; December 2003
345 : ! Revision for F90 reformatting: MJIacono; July 2006
346 : !
347 : ! The subroutines CMBGB16->CMBGB29 input the absorption coefficient
348 : ! data for each band, which are defined for 16 g-points and 14 spectral
349 : ! bands. The data are combined with appropriate weighting following the
350 : ! g-point mapping arrays specified in RRTMG_SW_INIT. Solar source
351 : ! function data in array SFLUXREF are combined without weighting. All
352 : ! g-point reduced data are put into new arrays for use in RRTMG_SW.
353 : !
354 : ! band 16: 2600-3250 cm-1 (low key- h2o,ch4; high key - ch4)
355 : !
356 : !-----------------------------------------------------------------------
357 :
358 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
359 : use rrsw_kg16, only : kao, kbo, selfrefo, forrefo, sfluxrefo, &
360 : ka, kb, selfref, forref, sfluxref
361 :
362 : ! ------- Local -------
363 : integer :: jn, jt, jp, igc, ipr, iprsm
364 : real(kind=r8) :: sumk, sumf
365 :
366 :
367 10240 : do jn = 1,9
368 56320 : do jt = 1,5
369 654336 : do jp = 1,13
370 599040 : iprsm = 0
371 4239360 : do igc = 1,ngc(1)
372 3594240 : sumk = 0.
373 13178880 : do ipr = 1, ngn(igc)
374 9584640 : iprsm = iprsm + 1
375 13178880 : sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm)
376 : enddo
377 4193280 : ka(jn,jt,jp,igc) = sumk
378 : enddo
379 : enddo
380 : enddo
381 : enddo
382 :
383 6144 : do jt = 1,5
384 246784 : do jp = 13,59
385 240640 : iprsm = 0
386 1689600 : do igc = 1,ngc(1)
387 1443840 : sumk = 0.
388 5294080 : do ipr = 1, ngn(igc)
389 3850240 : iprsm = iprsm + 1
390 5294080 : sumk = sumk + kbo(jt,jp,iprsm)*rwgt(iprsm)
391 : enddo
392 1684480 : kb(jt,jp,igc) = sumk
393 : enddo
394 : enddo
395 : enddo
396 :
397 11264 : do jt = 1,10
398 10240 : iprsm = 0
399 72704 : do igc = 1,ngc(1)
400 61440 : sumk = 0.
401 225280 : do ipr = 1, ngn(igc)
402 163840 : iprsm = iprsm + 1
403 225280 : sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm)
404 : enddo
405 71680 : selfref(jt,igc) = sumk
406 : enddo
407 : enddo
408 :
409 4096 : do jt = 1,3
410 3072 : iprsm = 0
411 22528 : do igc = 1,ngc(1)
412 18432 : sumk = 0.
413 67584 : do ipr = 1, ngn(igc)
414 49152 : iprsm = iprsm + 1
415 67584 : sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm)
416 : enddo
417 21504 : forref(jt,igc) = sumk
418 : enddo
419 : enddo
420 :
421 1024 : iprsm = 0
422 7168 : do igc = 1,ngc(1)
423 6144 : sumf = 0.
424 22528 : do ipr = 1, ngn(igc)
425 16384 : iprsm = iprsm + 1
426 22528 : sumf = sumf + sfluxrefo(iprsm)
427 : enddo
428 7168 : sfluxref(igc) = sumf
429 : enddo
430 :
431 1024 : end subroutine cmbgb16s
432 :
433 : !***************************************************************************
434 1024 : subroutine cmbgb17
435 : !***************************************************************************
436 : !
437 : ! band 17: 3250-4000 cm-1 (low - h2o,co2; high - h2o,co2)
438 : !-----------------------------------------------------------------------
439 :
440 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
441 : use rrsw_kg17, only : kao, kbo, selfrefo, forrefo, sfluxrefo, &
442 : ka, kb, selfref, forref, sfluxref
443 :
444 : ! ------- Local -------
445 : integer :: jn, jt, jp, igc, ipr, iprsm
446 : real(kind=r8) :: sumk, sumf
447 :
448 :
449 10240 : do jn = 1,9
450 56320 : do jt = 1,5
451 654336 : do jp = 1,13
452 599040 : iprsm = 0
453 7833600 : do igc = 1,ngc(2)
454 7188480 : sumk = 0.
455 16773120 : do ipr = 1, ngn(ngs(1)+igc)
456 9584640 : iprsm = iprsm + 1
457 16773120 : sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm+16)
458 : enddo
459 7787520 : ka(jn,jt,jp,igc) = sumk
460 : enddo
461 : enddo
462 : enddo
463 : enddo
464 :
465 6144 : do jn = 1,5
466 31744 : do jt = 1,5
467 1233920 : do jp = 13,59
468 1203200 : iprsm = 0
469 15667200 : do igc = 1,ngc(2)
470 14438400 : sumk = 0.
471 33689600 : do ipr = 1, ngn(ngs(1)+igc)
472 19251200 : iprsm = iprsm + 1
473 33689600 : sumk = sumk + kbo(jn,jt,jp,iprsm)*rwgt(iprsm+16)
474 : enddo
475 15641600 : kb(jn,jt,jp,igc) = sumk
476 : enddo
477 : enddo
478 : enddo
479 : enddo
480 :
481 11264 : do jt = 1,10
482 10240 : iprsm = 0
483 134144 : do igc = 1,ngc(2)
484 122880 : sumk = 0.
485 286720 : do ipr = 1, ngn(ngs(1)+igc)
486 163840 : iprsm = iprsm + 1
487 286720 : sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+16)
488 : enddo
489 133120 : selfref(jt,igc) = sumk
490 : enddo
491 : enddo
492 :
493 5120 : do jt = 1,4
494 4096 : iprsm = 0
495 54272 : do igc = 1,ngc(2)
496 49152 : sumk = 0.
497 114688 : do ipr = 1, ngn(ngs(1)+igc)
498 65536 : iprsm = iprsm + 1
499 114688 : sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+16)
500 : enddo
501 53248 : forref(jt,igc) = sumk
502 : enddo
503 : enddo
504 :
505 6144 : do jp = 1,5
506 5120 : iprsm = 0
507 67584 : do igc = 1,ngc(2)
508 61440 : sumf = 0.
509 143360 : do ipr = 1, ngn(ngs(1)+igc)
510 81920 : iprsm = iprsm + 1
511 143360 : sumf = sumf + sfluxrefo(iprsm,jp)
512 : enddo
513 66560 : sfluxref(igc,jp) = sumf
514 : enddo
515 : enddo
516 :
517 1024 : end subroutine cmbgb17
518 :
519 : !***************************************************************************
520 1024 : subroutine cmbgb18
521 : !***************************************************************************
522 : !
523 : ! band 18: 4000-4650 cm-1 (low - h2o,ch4; high - ch4)
524 : !-----------------------------------------------------------------------
525 :
526 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
527 : use rrsw_kg18, only : kao, kbo, selfrefo, forrefo, sfluxrefo, &
528 : ka, kb, selfref, forref, sfluxref
529 :
530 : ! ------- Local -------
531 : integer :: jn, jt, jp, igc, ipr, iprsm
532 : real(kind=r8) :: sumk, sumf
533 :
534 :
535 10240 : do jn = 1,9
536 56320 : do jt = 1,5
537 654336 : do jp = 1,13
538 599040 : iprsm = 0
539 5437440 : do igc = 1,ngc(3)
540 4792320 : sumk = 0.
541 14376960 : do ipr = 1, ngn(ngs(2)+igc)
542 9584640 : iprsm = iprsm + 1
543 14376960 : sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm+32)
544 : enddo
545 5391360 : ka(jn,jt,jp,igc) = sumk
546 : enddo
547 : enddo
548 : enddo
549 : enddo
550 :
551 6144 : do jt = 1,5
552 246784 : do jp = 13,59
553 240640 : iprsm = 0
554 2170880 : do igc = 1,ngc(3)
555 1925120 : sumk = 0.
556 5775360 : do ipr = 1, ngn(ngs(2)+igc)
557 3850240 : iprsm = iprsm + 1
558 5775360 : sumk = sumk + kbo(jt,jp,iprsm)*rwgt(iprsm+32)
559 : enddo
560 2165760 : kb(jt,jp,igc) = sumk
561 : enddo
562 : enddo
563 : enddo
564 :
565 11264 : do jt = 1,10
566 10240 : iprsm = 0
567 93184 : do igc = 1,ngc(3)
568 81920 : sumk = 0.
569 245760 : do ipr = 1, ngn(ngs(2)+igc)
570 163840 : iprsm = iprsm + 1
571 245760 : sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+32)
572 : enddo
573 92160 : selfref(jt,igc) = sumk
574 : enddo
575 : enddo
576 :
577 4096 : do jt = 1,3
578 3072 : iprsm = 0
579 28672 : do igc = 1,ngc(3)
580 24576 : sumk = 0.
581 73728 : do ipr = 1, ngn(ngs(2)+igc)
582 49152 : iprsm = iprsm + 1
583 73728 : sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+32)
584 : enddo
585 27648 : forref(jt,igc) = sumk
586 : enddo
587 : enddo
588 :
589 10240 : do jp = 1,9
590 9216 : iprsm = 0
591 83968 : do igc = 1,ngc(3)
592 73728 : sumf = 0.
593 221184 : do ipr = 1, ngn(ngs(2)+igc)
594 147456 : iprsm = iprsm + 1
595 221184 : sumf = sumf + sfluxrefo(iprsm,jp)
596 : enddo
597 82944 : sfluxref(igc,jp) = sumf
598 : enddo
599 : enddo
600 :
601 1024 : end subroutine cmbgb18
602 :
603 : !***************************************************************************
604 1024 : subroutine cmbgb19
605 : !***************************************************************************
606 : !
607 : ! band 19: 4650-5150 cm-1 (low - h2o,co2; high - co2)
608 : !-----------------------------------------------------------------------
609 :
610 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
611 : use rrsw_kg19, only : kao, kbo, selfrefo, forrefo, sfluxrefo, &
612 : ka, kb, selfref, forref, sfluxref
613 :
614 : ! ------- Local -------
615 : integer :: jn, jt, jp, igc, ipr, iprsm
616 : real(kind=r8) :: sumk, sumf
617 :
618 :
619 10240 : do jn = 1,9
620 56320 : do jt = 1,5
621 654336 : do jp = 1,13
622 599040 : iprsm = 0
623 5437440 : do igc = 1,ngc(4)
624 4792320 : sumk = 0.
625 14376960 : do ipr = 1, ngn(ngs(3)+igc)
626 9584640 : iprsm = iprsm + 1
627 14376960 : sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm+48)
628 : enddo
629 5391360 : ka(jn,jt,jp,igc) = sumk
630 : enddo
631 : enddo
632 : enddo
633 : enddo
634 :
635 6144 : do jt = 1,5
636 246784 : do jp = 13,59
637 240640 : iprsm = 0
638 2170880 : do igc = 1,ngc(4)
639 1925120 : sumk = 0.
640 5775360 : do ipr = 1, ngn(ngs(3)+igc)
641 3850240 : iprsm = iprsm + 1
642 5775360 : sumk = sumk + kbo(jt,jp,iprsm)*rwgt(iprsm+48)
643 : enddo
644 2165760 : kb(jt,jp,igc) = sumk
645 : enddo
646 : enddo
647 : enddo
648 :
649 11264 : do jt = 1,10
650 10240 : iprsm = 0
651 93184 : do igc = 1,ngc(4)
652 81920 : sumk = 0.
653 245760 : do ipr = 1, ngn(ngs(3)+igc)
654 163840 : iprsm = iprsm + 1
655 245760 : sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+48)
656 : enddo
657 92160 : selfref(jt,igc) = sumk
658 : enddo
659 : enddo
660 :
661 4096 : do jt = 1,3
662 3072 : iprsm = 0
663 28672 : do igc = 1,ngc(4)
664 24576 : sumk = 0.
665 73728 : do ipr = 1, ngn(ngs(3)+igc)
666 49152 : iprsm = iprsm + 1
667 73728 : sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+48)
668 : enddo
669 27648 : forref(jt,igc) = sumk
670 : enddo
671 : enddo
672 :
673 10240 : do jp = 1,9
674 9216 : iprsm = 0
675 83968 : do igc = 1,ngc(4)
676 73728 : sumf = 0.
677 221184 : do ipr = 1, ngn(ngs(3)+igc)
678 147456 : iprsm = iprsm + 1
679 221184 : sumf = sumf + sfluxrefo(iprsm,jp)
680 : enddo
681 82944 : sfluxref(igc,jp) = sumf
682 : enddo
683 : enddo
684 :
685 1024 : end subroutine cmbgb19
686 :
687 : !***************************************************************************
688 1024 : subroutine cmbgb20
689 : !***************************************************************************
690 : !
691 : ! band 20: 5150-6150 cm-1 (low - h2o; high - h2o)
692 : !-----------------------------------------------------------------------
693 :
694 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
695 : use rrsw_kg20, only : kao, kbo, selfrefo, forrefo, sfluxrefo, absch4o, &
696 : ka, kb, selfref, forref, sfluxref, absch4
697 :
698 : ! ------- Local -------
699 : integer :: jt, jp, igc, ipr, iprsm
700 : real(kind=r8) :: sumk, sumf1, sumf2
701 :
702 :
703 6144 : do jt = 1,5
704 71680 : do jp = 1,13
705 66560 : iprsm = 0
706 737280 : do igc = 1,ngc(5)
707 665600 : sumk = 0.
708 1730560 : do ipr = 1, ngn(ngs(4)+igc)
709 1064960 : iprsm = iprsm + 1
710 1730560 : sumk = sumk + kao(jt,jp,iprsm)*rwgt(iprsm+64)
711 : enddo
712 732160 : ka(jt,jp,igc) = sumk
713 : enddo
714 : enddo
715 246784 : do jp = 13,59
716 240640 : iprsm = 0
717 2652160 : do igc = 1,ngc(5)
718 2406400 : sumk = 0.
719 6256640 : do ipr = 1, ngn(ngs(4)+igc)
720 3850240 : iprsm = iprsm + 1
721 6256640 : sumk = sumk + kbo(jt,jp,iprsm)*rwgt(iprsm+64)
722 : enddo
723 2647040 : kb(jt,jp,igc) = sumk
724 : enddo
725 : enddo
726 : enddo
727 :
728 11264 : do jt = 1,10
729 10240 : iprsm = 0
730 113664 : do igc = 1,ngc(5)
731 102400 : sumk = 0.
732 266240 : do ipr = 1, ngn(ngs(4)+igc)
733 163840 : iprsm = iprsm + 1
734 266240 : sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+64)
735 : enddo
736 112640 : selfref(jt,igc) = sumk
737 : enddo
738 : enddo
739 :
740 5120 : do jt = 1,4
741 4096 : iprsm = 0
742 46080 : do igc = 1,ngc(5)
743 40960 : sumk = 0.
744 106496 : do ipr = 1, ngn(ngs(4)+igc)
745 65536 : iprsm = iprsm + 1
746 106496 : sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+64)
747 : enddo
748 45056 : forref(jt,igc) = sumk
749 : enddo
750 : enddo
751 :
752 1024 : iprsm = 0
753 11264 : do igc = 1,ngc(5)
754 10240 : sumf1 = 0.
755 10240 : sumf2 = 0.
756 26624 : do ipr = 1, ngn(ngs(4)+igc)
757 16384 : iprsm = iprsm + 1
758 16384 : sumf1 = sumf1 + sfluxrefo(iprsm)
759 26624 : sumf2 = sumf2 + absch4o(iprsm)*rwgt(iprsm+64)
760 : enddo
761 10240 : sfluxref(igc) = sumf1
762 11264 : absch4(igc) = sumf2
763 : enddo
764 :
765 1024 : end subroutine cmbgb20
766 :
767 : !***************************************************************************
768 1024 : subroutine cmbgb21
769 : !***************************************************************************
770 : !
771 : ! band 21: 6150-7700 cm-1 (low - h2o,co2; high - h2o,co2)
772 : !-----------------------------------------------------------------------
773 :
774 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
775 : use rrsw_kg21, only : kao, kbo, selfrefo, forrefo, sfluxrefo, &
776 : ka, kb, selfref, forref, sfluxref
777 :
778 : ! ------- Local -------
779 : integer :: jn, jt, jp, igc, ipr, iprsm
780 : real(kind=r8) :: sumk, sumf
781 :
782 :
783 10240 : do jn = 1,9
784 56320 : do jt = 1,5
785 654336 : do jp = 1,13
786 599040 : iprsm = 0
787 6635520 : do igc = 1,ngc(6)
788 5990400 : sumk = 0.
789 15575040 : do ipr = 1, ngn(ngs(5)+igc)
790 9584640 : iprsm = iprsm + 1
791 15575040 : sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm+80)
792 : enddo
793 6589440 : ka(jn,jt,jp,igc) = sumk
794 : enddo
795 : enddo
796 : enddo
797 : enddo
798 :
799 6144 : do jn = 1,5
800 31744 : do jt = 1,5
801 1233920 : do jp = 13,59
802 1203200 : iprsm = 0
803 13260800 : do igc = 1,ngc(6)
804 12032000 : sumk = 0.
805 31283200 : do ipr = 1, ngn(ngs(5)+igc)
806 19251200 : iprsm = iprsm + 1
807 31283200 : sumk = sumk + kbo(jn,jt,jp,iprsm)*rwgt(iprsm+80)
808 : enddo
809 13235200 : kb(jn,jt,jp,igc) = sumk
810 : enddo
811 : enddo
812 : enddo
813 : enddo
814 :
815 11264 : do jt = 1,10
816 10240 : iprsm = 0
817 113664 : do igc = 1,ngc(6)
818 102400 : sumk = 0.
819 266240 : do ipr = 1, ngn(ngs(5)+igc)
820 163840 : iprsm = iprsm + 1
821 266240 : sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+80)
822 : enddo
823 112640 : selfref(jt,igc) = sumk
824 : enddo
825 : enddo
826 :
827 5120 : do jt = 1,4
828 4096 : iprsm = 0
829 46080 : do igc = 1,ngc(6)
830 40960 : sumk = 0.
831 106496 : do ipr = 1, ngn(ngs(5)+igc)
832 65536 : iprsm = iprsm + 1
833 106496 : sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+80)
834 : enddo
835 45056 : forref(jt,igc) = sumk
836 : enddo
837 : enddo
838 :
839 10240 : do jp = 1,9
840 9216 : iprsm = 0
841 102400 : do igc = 1,ngc(6)
842 92160 : sumf = 0.
843 239616 : do ipr = 1, ngn(ngs(5)+igc)
844 147456 : iprsm = iprsm + 1
845 239616 : sumf = sumf + sfluxrefo(iprsm,jp)
846 : enddo
847 101376 : sfluxref(igc,jp) = sumf
848 : enddo
849 : enddo
850 :
851 1024 : end subroutine cmbgb21
852 :
853 : !***************************************************************************
854 1024 : subroutine cmbgb22
855 : !***************************************************************************
856 : !
857 : ! band 22: 7700-8050 cm-1 (low - h2o,o2; high - o2)
858 : !-----------------------------------------------------------------------
859 :
860 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
861 : use rrsw_kg22, only : kao, kbo, selfrefo, forrefo, sfluxrefo, &
862 : ka, kb, selfref, forref, sfluxref
863 :
864 : ! ------- Local -------
865 : integer :: jn, jt, jp, igc, ipr, iprsm
866 : real(kind=r8) :: sumk, sumf
867 :
868 :
869 10240 : do jn = 1,9
870 56320 : do jt = 1,5
871 654336 : do jp = 1,13
872 599040 : iprsm = 0
873 1843200 : do igc = 1,ngc(7)
874 1198080 : sumk = 0.
875 10782720 : do ipr = 1, ngn(ngs(6)+igc)
876 9584640 : iprsm = iprsm + 1
877 10782720 : sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm+96)
878 : enddo
879 1797120 : ka(jn,jt,jp,igc) = sumk
880 : enddo
881 : enddo
882 : enddo
883 : enddo
884 :
885 6144 : do jt = 1,5
886 246784 : do jp = 13,59
887 240640 : iprsm = 0
888 727040 : do igc = 1,ngc(7)
889 481280 : sumk = 0.
890 4331520 : do ipr = 1, ngn(ngs(6)+igc)
891 3850240 : iprsm = iprsm + 1
892 4331520 : sumk = sumk + kbo(jt,jp,iprsm)*rwgt(iprsm+96)
893 : enddo
894 721920 : kb(jt,jp,igc) = sumk
895 : enddo
896 : enddo
897 : enddo
898 :
899 11264 : do jt = 1,10
900 10240 : iprsm = 0
901 31744 : do igc = 1,ngc(7)
902 20480 : sumk = 0.
903 184320 : do ipr = 1, ngn(ngs(6)+igc)
904 163840 : iprsm = iprsm + 1
905 184320 : sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+96)
906 : enddo
907 30720 : selfref(jt,igc) = sumk
908 : enddo
909 : enddo
910 :
911 4096 : do jt = 1,3
912 3072 : iprsm = 0
913 10240 : do igc = 1,ngc(7)
914 6144 : sumk = 0.
915 55296 : do ipr = 1, ngn(ngs(6)+igc)
916 49152 : iprsm = iprsm + 1
917 55296 : sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+96)
918 : enddo
919 9216 : forref(jt,igc) = sumk
920 : enddo
921 : enddo
922 :
923 10240 : do jp = 1,9
924 9216 : iprsm = 0
925 28672 : do igc = 1,ngc(7)
926 18432 : sumf = 0.
927 165888 : do ipr = 1, ngn(ngs(6)+igc)
928 147456 : iprsm = iprsm + 1
929 165888 : sumf = sumf + sfluxrefo(iprsm,jp)
930 : enddo
931 27648 : sfluxref(igc,jp) = sumf
932 : enddo
933 : enddo
934 :
935 1024 : end subroutine cmbgb22
936 :
937 : !***************************************************************************
938 1024 : subroutine cmbgb23
939 : !***************************************************************************
940 : !
941 : ! band 23: 8050-12850 cm-1 (low - h2o; high - nothing)
942 : !-----------------------------------------------------------------------
943 :
944 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
945 : use rrsw_kg23, only : kao, selfrefo, forrefo, sfluxrefo, raylo, &
946 : ka, selfref, forref, sfluxref, rayl
947 :
948 : ! ------- Local -------
949 : integer :: jt, jp, igc, ipr, iprsm
950 : real(kind=r8) :: sumk, sumf1, sumf2
951 :
952 :
953 6144 : do jt = 1,5
954 72704 : do jp = 1,13
955 66560 : iprsm = 0
956 737280 : do igc = 1,ngc(8)
957 665600 : sumk = 0.
958 1730560 : do ipr = 1, ngn(ngs(7)+igc)
959 1064960 : iprsm = iprsm + 1
960 1730560 : sumk = sumk + kao(jt,jp,iprsm)*rwgt(iprsm+112)
961 : enddo
962 732160 : ka(jt,jp,igc) = sumk
963 : enddo
964 : enddo
965 : enddo
966 :
967 11264 : do jt = 1,10
968 10240 : iprsm = 0
969 113664 : do igc = 1,ngc(8)
970 102400 : sumk = 0.
971 266240 : do ipr = 1, ngn(ngs(7)+igc)
972 163840 : iprsm = iprsm + 1
973 266240 : sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+112)
974 : enddo
975 112640 : selfref(jt,igc) = sumk
976 : enddo
977 : enddo
978 :
979 4096 : do jt = 1,3
980 3072 : iprsm = 0
981 34816 : do igc = 1,ngc(8)
982 30720 : sumk = 0.
983 79872 : do ipr = 1, ngn(ngs(7)+igc)
984 49152 : iprsm = iprsm + 1
985 79872 : sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+112)
986 : enddo
987 33792 : forref(jt,igc) = sumk
988 : enddo
989 : enddo
990 :
991 1024 : iprsm = 0
992 11264 : do igc = 1,ngc(8)
993 10240 : sumf1 = 0.
994 10240 : sumf2 = 0.
995 26624 : do ipr = 1, ngn(ngs(7)+igc)
996 16384 : iprsm = iprsm + 1
997 16384 : sumf1 = sumf1 + sfluxrefo(iprsm)
998 26624 : sumf2 = sumf2 + raylo(iprsm)*rwgt(iprsm+112)
999 : enddo
1000 10240 : sfluxref(igc) = sumf1
1001 11264 : rayl(igc) = sumf2
1002 : enddo
1003 :
1004 1024 : end subroutine cmbgb23
1005 :
1006 : !***************************************************************************
1007 1024 : subroutine cmbgb24
1008 : !***************************************************************************
1009 : !
1010 : ! band 24: 12850-16000 cm-1 (low - h2o,o2; high - o2)
1011 : !-----------------------------------------------------------------------
1012 :
1013 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
1014 : use rrsw_kg24, only : kao, kbo, selfrefo, forrefo, sfluxrefo, &
1015 : abso3ao, abso3bo, raylao, raylbo, &
1016 : ka, kb, selfref, forref, sfluxref, &
1017 : abso3a, abso3b, rayla, raylb
1018 :
1019 : ! ------- Local -------
1020 : integer :: jn, jt, jp, igc, ipr, iprsm
1021 : real(kind=r8) :: sumk, sumf1, sumf2, sumf3
1022 :
1023 :
1024 10240 : do jn = 1,9
1025 56320 : do jt = 1,5
1026 654336 : do jp = 1,13
1027 599040 : iprsm = 0
1028 5437440 : do igc = 1,ngc(9)
1029 4792320 : sumk = 0.
1030 14376960 : do ipr = 1, ngn(ngs(8)+igc)
1031 9584640 : iprsm = iprsm + 1
1032 14376960 : sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm+128)
1033 : enddo
1034 5391360 : ka(jn,jt,jp,igc) = sumk
1035 : enddo
1036 : enddo
1037 : enddo
1038 : enddo
1039 :
1040 6144 : do jt = 1,5
1041 246784 : do jp = 13,59
1042 240640 : iprsm = 0
1043 2170880 : do igc = 1,ngc(9)
1044 1925120 : sumk = 0.
1045 5775360 : do ipr = 1, ngn(ngs(8)+igc)
1046 3850240 : iprsm = iprsm + 1
1047 5775360 : sumk = sumk + kbo(jt,jp,iprsm)*rwgt(iprsm+128)
1048 : enddo
1049 2165760 : kb(jt,jp,igc) = sumk
1050 : enddo
1051 : enddo
1052 : enddo
1053 :
1054 11264 : do jt = 1,10
1055 10240 : iprsm = 0
1056 93184 : do igc = 1,ngc(9)
1057 81920 : sumk = 0.
1058 245760 : do ipr = 1, ngn(ngs(8)+igc)
1059 163840 : iprsm = iprsm + 1
1060 245760 : sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+128)
1061 : enddo
1062 92160 : selfref(jt,igc) = sumk
1063 : enddo
1064 : enddo
1065 :
1066 4096 : do jt = 1,3
1067 3072 : iprsm = 0
1068 28672 : do igc = 1,ngc(9)
1069 24576 : sumk = 0.
1070 73728 : do ipr = 1, ngn(ngs(8)+igc)
1071 49152 : iprsm = iprsm + 1
1072 73728 : sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+128)
1073 : enddo
1074 27648 : forref(jt,igc) = sumk
1075 : enddo
1076 : enddo
1077 :
1078 1024 : iprsm = 0
1079 9216 : do igc = 1,ngc(9)
1080 8192 : sumf1 = 0.
1081 8192 : sumf2 = 0.
1082 8192 : sumf3 = 0.
1083 24576 : do ipr = 1, ngn(ngs(8)+igc)
1084 16384 : iprsm = iprsm + 1
1085 16384 : sumf1 = sumf1 + raylbo(iprsm)*rwgt(iprsm+128)
1086 16384 : sumf2 = sumf2 + abso3ao(iprsm)*rwgt(iprsm+128)
1087 24576 : sumf3 = sumf3 + abso3bo(iprsm)*rwgt(iprsm+128)
1088 : enddo
1089 8192 : raylb(igc) = sumf1
1090 8192 : abso3a(igc) = sumf2
1091 9216 : abso3b(igc) = sumf3
1092 : enddo
1093 :
1094 10240 : do jp = 1,9
1095 9216 : iprsm = 0
1096 83968 : do igc = 1,ngc(9)
1097 73728 : sumf1 = 0.
1098 73728 : sumf2 = 0.
1099 221184 : do ipr = 1, ngn(ngs(8)+igc)
1100 147456 : iprsm = iprsm + 1
1101 147456 : sumf1 = sumf1 + sfluxrefo(iprsm,jp)
1102 221184 : sumf2 = sumf2 + raylao(iprsm,jp)*rwgt(iprsm+128)
1103 : enddo
1104 73728 : sfluxref(igc,jp) = sumf1
1105 82944 : rayla(igc,jp) = sumf2
1106 : enddo
1107 : enddo
1108 :
1109 1024 : end subroutine cmbgb24
1110 :
1111 : !***************************************************************************
1112 1024 : subroutine cmbgb25
1113 : !***************************************************************************
1114 : !
1115 : ! band 25: 16000-22650 cm-1 (low - h2o; high - nothing)
1116 : !-----------------------------------------------------------------------
1117 :
1118 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
1119 : use rrsw_kg25, only : kao, sfluxrefo, &
1120 : abso3ao, abso3bo, raylo, &
1121 : ka, sfluxref, &
1122 : abso3a, abso3b, rayl
1123 :
1124 : ! ------- Local -------
1125 : integer :: jt, jp, igc, ipr, iprsm
1126 : real(kind=r8) :: sumk, sumf1, sumf2, sumf3, sumf4
1127 :
1128 :
1129 6144 : do jt = 1,5
1130 72704 : do jp = 1,13
1131 66560 : iprsm = 0
1132 471040 : do igc = 1,ngc(10)
1133 399360 : sumk = 0.
1134 1464320 : do ipr = 1, ngn(ngs(9)+igc)
1135 1064960 : iprsm = iprsm + 1
1136 1464320 : sumk = sumk + kao(jt,jp,iprsm)*rwgt(iprsm+144)
1137 : enddo
1138 465920 : ka(jt,jp,igc) = sumk
1139 : enddo
1140 : enddo
1141 : enddo
1142 :
1143 1024 : iprsm = 0
1144 7168 : do igc = 1,ngc(10)
1145 6144 : sumf1 = 0.
1146 6144 : sumf2 = 0.
1147 6144 : sumf3 = 0.
1148 6144 : sumf4 = 0.
1149 22528 : do ipr = 1, ngn(ngs(9)+igc)
1150 16384 : iprsm = iprsm + 1
1151 16384 : sumf1 = sumf1 + sfluxrefo(iprsm)
1152 16384 : sumf2 = sumf2 + abso3ao(iprsm)*rwgt(iprsm+144)
1153 16384 : sumf3 = sumf3 + abso3bo(iprsm)*rwgt(iprsm+144)
1154 22528 : sumf4 = sumf4 + raylo(iprsm)*rwgt(iprsm+144)
1155 : enddo
1156 6144 : sfluxref(igc) = sumf1
1157 6144 : abso3a(igc) = sumf2
1158 6144 : abso3b(igc) = sumf3
1159 7168 : rayl(igc) = sumf4
1160 : enddo
1161 :
1162 1024 : end subroutine cmbgb25
1163 :
1164 : !***************************************************************************
1165 1024 : subroutine cmbgb26
1166 : !***************************************************************************
1167 : !
1168 : ! band 26: 22650-29000 cm-1 (low - nothing; high - nothing)
1169 : !-----------------------------------------------------------------------
1170 :
1171 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
1172 : use rrsw_kg26, only : sfluxrefo, raylo, &
1173 : sfluxref, rayl
1174 :
1175 : ! ------- Local -------
1176 : integer :: igc, ipr, iprsm
1177 : real(kind=r8) :: sumf1, sumf2
1178 :
1179 :
1180 1024 : iprsm = 0
1181 7168 : do igc = 1,ngc(11)
1182 6144 : sumf1 = 0.
1183 6144 : sumf2 = 0.
1184 22528 : do ipr = 1, ngn(ngs(10)+igc)
1185 16384 : iprsm = iprsm + 1
1186 16384 : sumf1 = sumf1 + raylo(iprsm)*rwgt(iprsm+160)
1187 22528 : sumf2 = sumf2 + sfluxrefo(iprsm)
1188 : enddo
1189 6144 : rayl(igc) = sumf1
1190 7168 : sfluxref(igc) = sumf2
1191 : enddo
1192 :
1193 1024 : end subroutine cmbgb26
1194 :
1195 : !***************************************************************************
1196 1024 : subroutine cmbgb27
1197 : !***************************************************************************
1198 : !
1199 : ! band 27: 29000-38000 cm-1 (low - o3; high - o3)
1200 : !-----------------------------------------------------------------------
1201 :
1202 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
1203 : use rrsw_kg27, only : kao, kbo, sfluxrefo, raylo, &
1204 : ka, kb, sfluxref, rayl
1205 :
1206 : ! ------- Local -------
1207 : integer :: jt, jp, igc, ipr, iprsm
1208 : real(kind=r8) :: sumk, sumf1, sumf2
1209 :
1210 :
1211 6144 : do jt = 1,5
1212 71680 : do jp = 1,13
1213 66560 : iprsm = 0
1214 604160 : do igc = 1,ngc(12)
1215 532480 : sumk = 0.
1216 1597440 : do ipr = 1, ngn(ngs(11)+igc)
1217 1064960 : iprsm = iprsm + 1
1218 1597440 : sumk = sumk + kao(jt,jp,iprsm)*rwgt(iprsm+176)
1219 : enddo
1220 599040 : ka(jt,jp,igc) = sumk
1221 : enddo
1222 : enddo
1223 246784 : do jp = 13,59
1224 240640 : iprsm = 0
1225 2170880 : do igc = 1,ngc(12)
1226 1925120 : sumk = 0.
1227 5775360 : do ipr = 1, ngn(ngs(11)+igc)
1228 3850240 : iprsm = iprsm + 1
1229 5775360 : sumk = sumk + kbo(jt,jp,iprsm)*rwgt(iprsm+176)
1230 : enddo
1231 2165760 : kb(jt,jp,igc) = sumk
1232 : enddo
1233 : enddo
1234 : enddo
1235 :
1236 1024 : iprsm = 0
1237 9216 : do igc = 1,ngc(12)
1238 8192 : sumf1 = 0.
1239 8192 : sumf2 = 0.
1240 24576 : do ipr = 1, ngn(ngs(11)+igc)
1241 16384 : iprsm = iprsm + 1
1242 16384 : sumf1 = sumf1 + sfluxrefo(iprsm)
1243 24576 : sumf2 = sumf2 + raylo(iprsm)*rwgt(iprsm+176)
1244 : enddo
1245 8192 : sfluxref(igc) = sumf1
1246 9216 : rayl(igc) = sumf2
1247 : enddo
1248 :
1249 1024 : end subroutine cmbgb27
1250 :
1251 : !***************************************************************************
1252 1024 : subroutine cmbgb28
1253 : !***************************************************************************
1254 : !
1255 : ! band 28: 38000-50000 cm-1 (low - o3,o2; high - o3,o2)
1256 : !-----------------------------------------------------------------------
1257 :
1258 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
1259 : use rrsw_kg28, only : kao, kbo, sfluxrefo, &
1260 : ka, kb, sfluxref
1261 :
1262 : ! ------- Local -------
1263 : integer :: jn, jt, jp, igc, ipr, iprsm
1264 : real(kind=r8) :: sumk, sumf
1265 :
1266 :
1267 10240 : do jn = 1,9
1268 56320 : do jt = 1,5
1269 654336 : do jp = 1,13
1270 599040 : iprsm = 0
1271 4239360 : do igc = 1,ngc(13)
1272 3594240 : sumk = 0.
1273 13178880 : do ipr = 1, ngn(ngs(12)+igc)
1274 9584640 : iprsm = iprsm + 1
1275 13178880 : sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm+192)
1276 : enddo
1277 4193280 : ka(jn,jt,jp,igc) = sumk
1278 : enddo
1279 : enddo
1280 : enddo
1281 : enddo
1282 :
1283 6144 : do jn = 1,5
1284 31744 : do jt = 1,5
1285 1233920 : do jp = 13,59
1286 1203200 : iprsm = 0
1287 8448000 : do igc = 1,ngc(13)
1288 7219200 : sumk = 0.
1289 26470400 : do ipr = 1, ngn(ngs(12)+igc)
1290 19251200 : iprsm = iprsm + 1
1291 26470400 : sumk = sumk + kbo(jn,jt,jp,iprsm)*rwgt(iprsm+192)
1292 : enddo
1293 8422400 : kb(jn,jt,jp,igc) = sumk
1294 : enddo
1295 : enddo
1296 : enddo
1297 : enddo
1298 :
1299 6144 : do jp = 1,5
1300 5120 : iprsm = 0
1301 36864 : do igc = 1,ngc(13)
1302 30720 : sumf = 0.
1303 112640 : do ipr = 1, ngn(ngs(12)+igc)
1304 81920 : iprsm = iprsm + 1
1305 112640 : sumf = sumf + sfluxrefo(iprsm,jp)
1306 : enddo
1307 35840 : sfluxref(igc,jp) = sumf
1308 : enddo
1309 : enddo
1310 :
1311 1024 : end subroutine cmbgb28
1312 :
1313 : !***************************************************************************
1314 1024 : subroutine cmbgb29
1315 : !***************************************************************************
1316 : !
1317 : ! band 29: 820-2600 cm-1 (low - h2o; high - co2)
1318 : !-----------------------------------------------------------------------
1319 :
1320 : use rrsw_wvn, only : ngc, ngs, ngn, wt, rwgt
1321 : use rrsw_kg29, only : kao, kbo, selfrefo, forrefo, sfluxrefo, &
1322 : absh2oo, absco2o, &
1323 : ka, kb, selfref, forref, sfluxref, &
1324 : absh2o, absco2
1325 :
1326 : ! ------- Local -------
1327 : integer :: jt, jp, igc, ipr, iprsm
1328 : real(kind=r8) :: sumk, sumf1, sumf2, sumf3
1329 :
1330 :
1331 6144 : do jt = 1,5
1332 71680 : do jp = 1,13
1333 66560 : iprsm = 0
1334 870400 : do igc = 1,ngc(14)
1335 798720 : sumk = 0.
1336 1863680 : do ipr = 1, ngn(ngs(13)+igc)
1337 1064960 : iprsm = iprsm + 1
1338 1863680 : sumk = sumk + kao(jt,jp,iprsm)*rwgt(iprsm+208)
1339 : enddo
1340 865280 : ka(jt,jp,igc) = sumk
1341 : enddo
1342 : enddo
1343 246784 : do jp = 13,59
1344 240640 : iprsm = 0
1345 3133440 : do igc = 1,ngc(14)
1346 2887680 : sumk = 0.
1347 6737920 : do ipr = 1, ngn(ngs(13)+igc)
1348 3850240 : iprsm = iprsm + 1
1349 6737920 : sumk = sumk + kbo(jt,jp,iprsm)*rwgt(iprsm+208)
1350 : enddo
1351 3128320 : kb(jt,jp,igc) = sumk
1352 : enddo
1353 : enddo
1354 : enddo
1355 :
1356 11264 : do jt = 1,10
1357 10240 : iprsm = 0
1358 134144 : do igc = 1,ngc(14)
1359 122880 : sumk = 0.
1360 286720 : do ipr = 1, ngn(ngs(13)+igc)
1361 163840 : iprsm = iprsm + 1
1362 286720 : sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+208)
1363 : enddo
1364 133120 : selfref(jt,igc) = sumk
1365 : enddo
1366 : enddo
1367 :
1368 5120 : do jt = 1,4
1369 4096 : iprsm = 0
1370 54272 : do igc = 1,ngc(14)
1371 49152 : sumk = 0.
1372 114688 : do ipr = 1, ngn(ngs(13)+igc)
1373 65536 : iprsm = iprsm + 1
1374 114688 : sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+208)
1375 : enddo
1376 53248 : forref(jt,igc) = sumk
1377 : enddo
1378 : enddo
1379 :
1380 1024 : iprsm = 0
1381 13312 : do igc = 1,ngc(14)
1382 12288 : sumf1 = 0.
1383 12288 : sumf2 = 0.
1384 12288 : sumf3 = 0.
1385 28672 : do ipr = 1, ngn(ngs(13)+igc)
1386 16384 : iprsm = iprsm + 1
1387 16384 : sumf1 = sumf1 + sfluxrefo(iprsm)
1388 16384 : sumf2 = sumf2 + absco2o(iprsm)*rwgt(iprsm+208)
1389 28672 : sumf3 = sumf3 + absh2oo(iprsm)*rwgt(iprsm+208)
1390 : enddo
1391 12288 : sfluxref(igc) = sumf1
1392 12288 : absco2(igc) = sumf2
1393 13312 : absh2o(igc) = sumf3
1394 : enddo
1395 :
1396 1024 : end subroutine cmbgb29
1397 :
1398 : !***************************************************************************
1399 :
1400 :
1401 : end module rrtmg_sw_init
1402 :
1403 :
|